[4] | 1 | !> \file resol_adv_diff_2D-juin2009.f90 |
---|
| 2 | !!Resoud l equation adv-diffusion par une methode de relaxation |
---|
| 3 | !< |
---|
| 4 | |
---|
| 5 | !> \namespace reso_adv_diff_2D_vect |
---|
| 6 | !! Resoud l equation adv-diffusion par une methode de relaxation |
---|
| 7 | !! @note - passage du vecteur et de l epaisseur en dummy |
---|
| 8 | !! @note - epaisseur prescrite selon un masque |
---|
| 9 | !! @note use module3D_phy |
---|
| 10 | !! @todo essayer dans le futur la methode de Richard dUX/dx = H dU/dx + U dH/dx |
---|
| 11 | |
---|
| 12 | |
---|
| 13 | module reso_adv_diff_2D_vect |
---|
| 14 | |
---|
| 15 | use module3D_phy |
---|
| 16 | |
---|
| 17 | |
---|
| 18 | |
---|
| 19 | implicit none |
---|
| 20 | real :: omega !< parametre schema temporel de la resolution partie diffusion |
---|
| 21 | real :: mu_adv !< parametre schema temporel de la resolution advection |
---|
| 22 | real :: upwind !< schema spatial pour l'advection |
---|
| 23 | |
---|
| 24 | |
---|
| 25 | integer, dimension(2) :: ijmax ! position de maxval |
---|
| 26 | integer :: iFAIL |
---|
| 27 | |
---|
| 28 | contains |
---|
| 29 | |
---|
| 30 | !> initialise init_reso_adv_diff_2D |
---|
| 31 | !! definition des parametres qui gerent la resolution |
---|
| 32 | !! @return omega, mu_adv, upwind |
---|
| 33 | |
---|
| 34 | subroutine init_reso_adv_diff_2D |
---|
| 35 | |
---|
| 36 | |
---|
| 37 | write(num_rep_42,*)'partie diffusive' |
---|
| 38 | write(num_rep_42,*)'le type de schema temporel diffusif depend de omega' |
---|
| 39 | write(num_rep_42,*)'0 -> explicite, 0.5 -> Crank-Nicolson' |
---|
| 40 | write(num_rep_42,*)'1 -> semi implicite, >1 -> over-implicite' |
---|
| 41 | |
---|
| 42 | omega = 2.5 |
---|
| 43 | |
---|
| 44 | write(num_rep_42,*)'omega = ',omega |
---|
| 45 | |
---|
| 46 | write(num_rep_42,*) |
---|
| 47 | write(num_rep_42,*)'partie advective' |
---|
| 48 | write(num_rep_42,*)' le schéma temporel advectif dépend de mu' |
---|
| 49 | write(num_rep_42,*)'0 -> explicite, 0.5 -> Crank-Nicolson' |
---|
| 50 | write(num_rep_42,*)'1 -> semi implicite, >1 -> over-implicite' |
---|
| 51 | |
---|
| 52 | mu_adv=1. |
---|
| 53 | write(num_rep_42,*)'mu_adv = ',mu_adv |
---|
| 54 | |
---|
| 55 | write(num_rep_42,*)' le schéma spatial dépend de upwind' |
---|
| 56 | write(num_rep_42,*)'1 -> schema upwind, 0.5 -> schema centre' |
---|
| 57 | |
---|
| 58 | upwind=1 |
---|
| 59 | write(num_rep_42,*)'upwind = ',upwind |
---|
| 60 | |
---|
| 61 | |
---|
| 62 | write(num_rep_42,*)'------------------------------------------------------' |
---|
| 63 | |
---|
| 64 | return |
---|
| 65 | end subroutine init_reso_adv_diff_2D |
---|
| 66 | !------------------------------------------------------------------ |
---|
| 67 | |
---|
| 68 | !> subroutine resol_adv_diff_2D_Hpresc |
---|
| 69 | !! definition des parametres qui gerent la resolution |
---|
| 70 | !! @param Dfx,Dfy termes diffusifs |
---|
| 71 | !! @param advx,advy termes advectifs |
---|
| 72 | !! @param vieuxH ancienne valeur de H |
---|
| 73 | !! @param H_presc la valeur H prescrit pour certains noeuds |
---|
| 74 | !! @param i_Hpresc le masque ou la valeur de H est prescrite |
---|
| 75 | !! @param bil le bilan pour la colonne |
---|
| 76 | !! @return newH la valeur de H apres le pas de temps |
---|
| 77 | |
---|
| 78 | subroutine resol_adv_diff_2D_vect(Dfx,Dfy,advx,advy,H_presc,i_Hpresc,bil,vieuxH,newH) |
---|
| 79 | |
---|
| 80 | implicit none |
---|
| 81 | |
---|
| 82 | real,dimension(nx,ny), intent(in) :: Dfx !< terme diffusif selon x |
---|
| 83 | real,dimension(nx,ny), intent(in) :: Dfy !< terme diffusif selon y |
---|
| 84 | real,dimension(nx,ny), intent(in) :: Advx !< terme advectif selon x |
---|
| 85 | real,dimension(nx,ny), intent(in) :: Advy !< terme advectif selon y |
---|
| 86 | real,dimension(nx,ny), intent(in) :: vieuxH !< ancienne valeur de H |
---|
| 87 | real,dimension(nx,ny), intent(out):: newH !< nouvelle valeur de H |
---|
| 88 | |
---|
| 89 | real,dimension(nx,ny), intent(in) :: bil !< bilan de masse pour la colonne |
---|
| 90 | |
---|
| 91 | real,dimension(nx,ny), intent(in) :: H_presc !< H value if prescribed |
---|
| 92 | integer,dimension(nx,ny), intent(in) :: i_Hpresc !< 1 if H is prescribed on this node, else 0 |
---|
| 93 | |
---|
| 94 | |
---|
| 95 | ! tableaux de travail. resolution M H = Frelax |
---|
| 96 | real,dimension(nx,ny) :: crelax !< diagnonale de M |
---|
| 97 | real,dimension(nx,ny) :: arelax !< sous diagonale selon x |
---|
| 98 | real,dimension(nx,ny) :: brelax !< sur diagonale selon x |
---|
| 99 | real,dimension(nx,ny) :: drelax !< sous diagonale selon y |
---|
| 100 | real,dimension(nx,ny) :: erelax !< sur diagonale selon y |
---|
| 101 | real,dimension(nx,ny) :: frelax !< vecteur |
---|
| 102 | real,dimension(nx,ny) :: c_west !< sur demi mailles Ux |
---|
| 103 | real,dimension(nx,ny) :: c_east !< sur demi mailles Ux |
---|
| 104 | real,dimension(nx,ny) :: c_north !< sur demi mailles Uy |
---|
| 105 | real,dimension(nx,ny) :: c_south !<sur demi mailles Ux |
---|
| 106 | |
---|
| 107 | real,dimension(nx,ny) :: bdx !< pente socle |
---|
| 108 | real,dimension(nx,ny) :: bdy !< pente socle |
---|
| 109 | |
---|
| 110 | real,dimension(nx,ny) :: hdx !< pente epaisseur |
---|
| 111 | real,dimension(nx,ny) :: hdy !< pente epaisseur |
---|
| 112 | |
---|
| 113 | real :: frdx,frdy !< pour calcul frelax : termes diffusion |
---|
| 114 | real :: fraxw,fraxe,frays,frayn !< termes advection |
---|
| 115 | |
---|
| 116 | real,dimension(nx,ny) :: deltah ! dans calcul relax |
---|
| 117 | real :: delh ! dans calcul relax |
---|
| 118 | real :: testh ! dans calcul relax |
---|
| 119 | |
---|
| 120 | logical :: stopp |
---|
| 121 | integer :: ntour |
---|
| 122 | real :: reste |
---|
| 123 | |
---|
| 124 | |
---|
| 125 | integer :: it1,it2,jt1,jt2 ! pour des tests d'asymétrie |
---|
| 126 | |
---|
| 127 | |
---|
| 128 | if (itracebug.eq.1) call tracebug(' Entree dans routine resolution_diffusion') |
---|
| 129 | |
---|
| 130 | ! calcul de bdx et hdx |
---|
| 131 | hdx(:,:)=dx1*(vieuxH(:,:)-eoshift(vieuxH(:,:),shift=-1,boundary=0.,dim=1)) |
---|
| 132 | hdy(:,:)=dx1*(vieuxH(:,:)-eoshift(vieuxH(:,:),shift=-1,boundary=0.,dim=2)) |
---|
| 133 | bdx(:,:)=dx1*(B(:,:)-eoshift(B(:,:),shift=-1,boundary=0.,dim=1)) |
---|
| 134 | bdy(:,:)=dx1*(B(:,:)-eoshift(B(:,:),shift=-1,boundary=0.,dim=2)) |
---|
| 135 | |
---|
| 136 | |
---|
| 137 | ! initialisations (qui feront aussi les conditions aux limites) |
---|
| 138 | Arelax(:,:)=0. |
---|
| 139 | Brelax(:,:)=0. |
---|
| 140 | Drelax(:,:)=0. |
---|
| 141 | Erelax(:,:)=0. |
---|
| 142 | Crelax(:,:)=1. |
---|
| 143 | Frelax(:,:)=0. |
---|
| 144 | DeltaH(:,:)=0. |
---|
| 145 | |
---|
| 146 | |
---|
| 147 | ! schema spatial |
---|
| 148 | |
---|
| 149 | if (upwind.eq.1.) then !schema amont |
---|
| 150 | |
---|
| 151 | where (Advx(:,:).ge.0.) |
---|
| 152 | c_west(:,:)=1. |
---|
| 153 | c_east(:,:)=0. |
---|
| 154 | elsewhere |
---|
| 155 | c_west(:,:)=0. |
---|
| 156 | c_east(:,:)=1. |
---|
| 157 | end where |
---|
| 158 | |
---|
| 159 | where (Advy(:,:).ge.0.) |
---|
| 160 | c_south(:,:)=1. |
---|
| 161 | c_north(:,:)=0. |
---|
| 162 | elsewhere |
---|
| 163 | c_south(:,:)=0. |
---|
| 164 | c_north(:,:)=1. |
---|
| 165 | end where |
---|
| 166 | |
---|
| 167 | else if (upwind.lt.1.) then ! schema centre |
---|
| 168 | c_west(:,:)=0.5 |
---|
| 169 | c_east(:,:)=0.5 |
---|
| 170 | c_south(:,:)=0.5 |
---|
| 171 | c_north(:,:)=0.5 |
---|
| 172 | end if |
---|
| 173 | |
---|
| 174 | ! attribution des elements des diagonales |
---|
| 175 | do j=2,ny-1 |
---|
| 176 | do i=2,nx-1 |
---|
| 177 | |
---|
| 178 | ! sous diagonale en x |
---|
| 179 | arelax(i,j)=-omega*Dtdx2*Dfx(i,j) & ! partie diffusive en x |
---|
| 180 | -mu_adv*dtdx*c_west(i,j)*Advx(i,j) ! partie advective en x |
---|
| 181 | |
---|
| 182 | ! sur diagonale en x |
---|
| 183 | brelax(i,j)=-omega*Dtdx2*Dfx(i+1,j) & ! partie diffusive |
---|
| 184 | +mu_adv*dtdx*c_east(i+1,j)*Advx(i+1,j) ! partie advective |
---|
| 185 | |
---|
| 186 | ! sous diagonale en y |
---|
| 187 | drelax(i,j)=-omega*Dtdx2*Dfy(i,j) & ! partie diffusive en y |
---|
| 188 | -mu_adv*dtdx*c_south(i,j)*Advy(i,j) ! partie advective en y |
---|
| 189 | |
---|
| 190 | ! sur diagonale en y |
---|
| 191 | erelax(i,j)=-omega*Dtdx2*Dfy(i,j+1) & ! partie diffusive |
---|
| 192 | +mu_adv*dtdx*c_north(i,j+1)*Advy(i,j+1) ! partie advective |
---|
| 193 | |
---|
| 194 | |
---|
| 195 | |
---|
| 196 | ! diagonale |
---|
| 197 | crelax(i,j)=omega*Dtdx2*((Dfx(i+1,j)+Dfx(i,j)) & ! partie diffusive en x |
---|
| 198 | +(Dfy(i,j+1)+Dfy(i,j))) ! partie diffusive en y |
---|
| 199 | crelax(i,j)=crelax(i,j)+mu_adv*dtdx* & |
---|
| 200 | ( (c_west(i+1,j)*Advx(i+1,j)-c_east(i,j)*Advx(i,j)) & !partie advective en x |
---|
| 201 | +(c_south(i,j+1)*Advy(i,j+1)-c_north(i,j)*Advy(i,j))) !partie advective en y |
---|
| 202 | crelax(i,j)=1.+crelax(i,j) ! partie temporelle |
---|
| 203 | |
---|
| 204 | |
---|
| 205 | ! terme du vecteur |
---|
| 206 | |
---|
| 207 | frdx= -Dfx(i,j) * (Bdx(i,j) +(1.-omega)*Hdx(i,j)) & ! partie diffusive en x |
---|
| 208 | +Dfx(i+1,j)*(Bdx(i+1,j)+(1.-omega)*Hdx(i+1,j)) |
---|
| 209 | |
---|
| 210 | frdy= -Dfy(i,j) * (Bdy(i,j) +(1.-omega)*Hdy(i,j)) & ! partie diffusive en y |
---|
| 211 | +Dfy(i,j+1)*(Bdy(i,j+1)+(1.-omega)*Hdy(i,j+1)) |
---|
| 212 | |
---|
| 213 | fraxw= -c_west(i,j)* Advx(i,j) * vieuxH(i-1,j) & ! partie advective en x |
---|
| 214 | +c_west(i+1,j)*Advx(i+1,j)*vieuxH(i,j) ! venant de l'west |
---|
| 215 | |
---|
| 216 | fraxe= -c_east(i,j) * Advx(i,j) * vieuxH(i,j) & ! partie advective en x |
---|
| 217 | +c_east(i+1,j)*Advx(i+1,j)*vieuxH(i+1,j) ! venant de l'est |
---|
| 218 | |
---|
| 219 | frays= -c_south(i,j) * Advy(i,j) * vieuxH(i,j-1) & ! partie advective en y |
---|
| 220 | +c_south(i,j+1)*Advy(i,j+1)*vieuxH(i,j) ! venant du sud |
---|
| 221 | |
---|
| 222 | frayn= -c_north(i,j) * Advy(i,j) * vieuxH(i,j) & ! partie advective en y |
---|
| 223 | +c_north(i,j+1)*Advy(i,j+1)*vieuxH(i,j+1) ! venant du nord |
---|
| 224 | |
---|
| 225 | |
---|
| 226 | |
---|
| 227 | |
---|
| 228 | |
---|
| 229 | frelax(i,j)=vieuxH(i,j)+Dt*bil(i,j)+dtdx*(frdx+frdy) & |
---|
| 230 | + (1.-mu_adv)*dtdx*((fraxw+fraxe)+(frays+frayn)) |
---|
| 231 | |
---|
| 232 | end do |
---|
| 233 | end do |
---|
| 234 | |
---|
| 235 | |
---|
| 236 | where (i_hpresc(:,:) .eq.1) ! thickness prescribed |
---|
| 237 | frelax(:,:) = H_presc(:,:) |
---|
| 238 | arelax(:,:) = 0. |
---|
| 239 | brelax(:,:) = 0. |
---|
| 240 | crelax(:,:) = 1. |
---|
| 241 | drelax(:,:) = 0. |
---|
| 242 | erelax(:,:) = 0. |
---|
| 243 | end where |
---|
| 244 | |
---|
| 245 | stopp = .false. |
---|
| 246 | ntour=0 |
---|
| 247 | |
---|
| 248 | |
---|
| 249 | |
---|
| 250 | relax_loop: do while(.not.stopp) |
---|
| 251 | ntour=ntour+1 |
---|
| 252 | |
---|
| 253 | do j=2,ny-1 |
---|
| 254 | do i=2,nx-1 |
---|
| 255 | |
---|
| 256 | reste = (((arelax(i,j)*newH(i-1,j) +drelax(i,j)*newH(i,j-1)) & |
---|
| 257 | + (brelax(i,j)*newH(i+1,j) + erelax(i,j)*newH(i,j+1))) & |
---|
| 258 | + crelax(i,j)*newH(i,j))- frelax(i,j) |
---|
| 259 | |
---|
| 260 | if (ntour.eq.1) debug_3D(i,j,49)=(((arelax(i,j)*newH(i-1,j) +drelax(i,j)*newH(i,j-1)) & |
---|
| 261 | + (brelax(i,j)*newH(i+1,j) + erelax(i,j)*newH(i,j+1))) & |
---|
| 262 | + crelax(i,j)*newH(i,j)) |
---|
| 263 | |
---|
| 264 | |
---|
| 265 | deltaH(i,j) = reste/crelax(i,j) |
---|
| 266 | |
---|
| 267 | end do |
---|
| 268 | end do |
---|
| 269 | |
---|
| 270 | debug_3D(:,:,50)=arelax(:,:) |
---|
| 271 | debug_3D(:,:,51)=brelax(:,:) |
---|
| 272 | debug_3D(:,:,52)=crelax(:,:) |
---|
| 273 | debug_3D(:,:,53)=drelax(:,:) |
---|
| 274 | debug_3D(:,:,54)=erelax(:,:) |
---|
| 275 | debug_3D(:,:,55)=frelax(:,:) |
---|
| 276 | |
---|
| 277 | |
---|
| 278 | |
---|
| 279 | |
---|
| 280 | newH(:,:)=newH(:,:)-deltaH(:,:) |
---|
| 281 | |
---|
| 282 | |
---|
| 283 | ! critere d'arret: |
---|
| 284 | ! ---------------- |
---|
| 285 | |
---|
| 286 | delh=0 |
---|
| 287 | |
---|
| 288 | |
---|
| 289 | do j=2,ny-1 |
---|
| 290 | do i=2,nx-1 |
---|
| 291 | delh=delh+deltaH(i,j)**2 |
---|
| 292 | end do |
---|
| 293 | end do |
---|
| 294 | |
---|
| 295 | if (delh.gt.0.) then |
---|
| 296 | testh=sqrt(delh)/((nx-2)*(ny-2)) |
---|
| 297 | else |
---|
| 298 | testh=0. |
---|
| 299 | endif |
---|
| 300 | stopp = (testh.lt.1.e-4).or.(ntour.gt.100) |
---|
| 301 | |
---|
| 302 | |
---|
| 303 | end do relax_loop |
---|
| 304 | |
---|
| 305 | |
---|
| 306 | |
---|
| 307 | return |
---|
| 308 | |
---|
| 309 | end subroutine resol_adv_diff_2D_vect |
---|
| 310 | |
---|
| 311 | end module reso_adv_diff_2D_vect |
---|
| 312 | |
---|
| 313 | |
---|
| 314 | |
---|
| 315 | |
---|