1 | !> \file resol_adv_diff_2D-paral.f90 |
---|
2 | !!Resoud l equation adv-diffusion par une methode de relaxation |
---|
3 | !< |
---|
4 | |
---|
5 | !> \namespace reso_adv_diff_2D_paral |
---|
6 | !! Resoud l equation adv-diffusion par une methode de relaxation |
---|
7 | !! @note - adapté a un domaine parallele pour les experiences MISMIP |
---|
8 | !! @note - passage du vecteur et de l epaisseur en dummy |
---|
9 | !! @note - epaisseur prescrite selon un masque |
---|
10 | !! @note use module3D_phy |
---|
11 | |
---|
12 | !< |
---|
13 | |
---|
14 | module reso_adv_diff_2D_paral |
---|
15 | |
---|
16 | use module3D_phy |
---|
17 | implicit none |
---|
18 | |
---|
19 | |
---|
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 | ! tableaux de travail. resolution M H = Frelax |
---|
25 | real,dimension(nx,ny) :: crelax !< diagnonale de M |
---|
26 | real,dimension(nx,ny) :: arelax !< sous diagonale selon x |
---|
27 | real,dimension(nx,ny) :: brelax !< sur diagonale selon x |
---|
28 | real,dimension(nx,ny) :: drelax !< sous diagonale selon y |
---|
29 | real,dimension(nx,ny) :: erelax !< sur diagonale selon y |
---|
30 | real,dimension(nx,ny) :: frelax !< vecteur |
---|
31 | real,dimension(nx,ny) :: c_west !< sur demi mailles Ux |
---|
32 | real,dimension(nx,ny) :: c_east !< sur demi mailles Ux |
---|
33 | real,dimension(nx,ny) :: c_north !< sur demi mailles Uy |
---|
34 | real,dimension(nx,ny) :: c_south !< sur demi mailles Uy |
---|
35 | |
---|
36 | real,dimension(nx,ny) :: bdx !< pente socle |
---|
37 | real,dimension(nx,ny) :: bdy !< pente socle |
---|
38 | |
---|
39 | real,dimension(nx,ny) :: hdx !< pente epaisseur |
---|
40 | real,dimension(nx,ny) :: hdy !< pente epaisseur |
---|
41 | |
---|
42 | |
---|
43 | integer, dimension(2) :: ijmax ! position de maxval |
---|
44 | integer :: iFAIL |
---|
45 | integer :: paral_x ! 1-> bords paralleles a x |
---|
46 | integer :: paral_y ! 1-> bords paralleles a y |
---|
47 | |
---|
48 | contains |
---|
49 | |
---|
50 | !> initialise init_reso_adv_diff_2D |
---|
51 | !! definition des parametres qui gerent la resolution |
---|
52 | !! @return omega, mu_adv, upwind |
---|
53 | |
---|
54 | subroutine init_reso_adv_diff_2D |
---|
55 | |
---|
56 | |
---|
57 | write(num_rep_42,*)'partie diffusive' |
---|
58 | write(num_rep_42,*)'le type de schema temporel diffusif depend de omega' |
---|
59 | write(num_rep_42,*)'0 -> explicite, 0.5 -> Crank-Nicolson' |
---|
60 | write(num_rep_42,*)'1 -> semi implicite, >1 -> over-implicite' |
---|
61 | |
---|
62 | ! over implicite propose par Hindmasrsh |
---|
63 | |
---|
64 | omega = 2.5 |
---|
65 | |
---|
66 | write(num_rep_42,*)'omega = ',omega |
---|
67 | |
---|
68 | write(num_rep_42,*) |
---|
69 | write(num_rep_42,*)'partie advective' |
---|
70 | write(num_rep_42,*)' le schéma temporel advectif dépend de mu' |
---|
71 | write(num_rep_42,*)'0 -> explicite, 0.5 -> Crank-Nicolson' |
---|
72 | write(num_rep_42,*)'1 -> semi implicite, >1 -> over-implicite' |
---|
73 | |
---|
74 | mu_adv=1. ! l'over implicite ne s'applique pas a l'equation advective |
---|
75 | write(num_rep_42,*)'mu_adv = ',mu_adv |
---|
76 | |
---|
77 | write(num_rep_42,*)' le schéma spatial dépend de upwind' |
---|
78 | write(num_rep_42,*)'1 -> schema upwind, 0.5 -> schema centre' |
---|
79 | |
---|
80 | upwind=1 |
---|
81 | |
---|
82 | write(num_rep_42,*)'upwind = ',upwind |
---|
83 | |
---|
84 | |
---|
85 | write(num_rep_42,*)'------------------------------------------------------' |
---|
86 | |
---|
87 | |
---|
88 | ! recherche la plus grande dimension horizontale |
---|
89 | ! c'est le long de la plus grande dimension que c'est parallele |
---|
90 | |
---|
91 | if (nx.eq.max(nx,ny)) then ! bords paralleles a x |
---|
92 | paral_x = 1 |
---|
93 | paral_y = 0 |
---|
94 | else ! bords paralleles a y |
---|
95 | paral_x = 0 |
---|
96 | paral_y = 1 |
---|
97 | end if |
---|
98 | |
---|
99 | |
---|
100 | |
---|
101 | return |
---|
102 | end subroutine init_reso_adv_diff_2D |
---|
103 | !------------------------------------------------------------------ |
---|
104 | |
---|
105 | !> subroutine resol_adv_diff_2D_Hpresc |
---|
106 | !! definition des parametres qui gerent la resolution |
---|
107 | !! @param Dfx,Dfy termes diffusifs |
---|
108 | !! @param advx,advy termes advectifs |
---|
109 | !! @param vieuxH ancienne valeur de H |
---|
110 | !! @param H_presc la valeur H prescrit pour certains noeuds |
---|
111 | !! @param i_Hpresc le masque ou la valeur de H est prescrite |
---|
112 | !! @param bil le bilan pour la colonne |
---|
113 | !! @return newH la valeur de H apres le pas de temps |
---|
114 | |
---|
115 | subroutine resol_adv_diff_2D_vect(Dfx,Dfy,advx,advy,H_presc,i_Hpresc,bil,vieuxH,newH) |
---|
116 | |
---|
117 | implicit none |
---|
118 | |
---|
119 | real,dimension(nx,ny), intent(in) :: Dfx !< terme diffusif selon x |
---|
120 | real,dimension(nx,ny), intent(in) :: Dfy !< terme diffusif selon y |
---|
121 | real,dimension(nx,ny), intent(in) :: Advx !< terme advectif selon x |
---|
122 | real,dimension(nx,ny), intent(in) :: Advy !< terme advectif selon y |
---|
123 | real,dimension(nx,ny), intent(in) :: vieuxH !< ancienne valeur de H |
---|
124 | real,dimension(nx,ny), intent(out):: newH !< nouvelle valeur de H |
---|
125 | |
---|
126 | real,dimension(nx,ny), intent(in) :: bil !< bilan de masse pour la colonne |
---|
127 | |
---|
128 | real,dimension(nx,ny), intent(in) :: H_presc !< H value if prescribed |
---|
129 | integer,dimension(nx,ny), intent(in) :: i_Hpresc !< 1 if H is prescribed on this node, else 0 |
---|
130 | |
---|
131 | |
---|
132 | real :: frdx,frdy !< pour calcul frelax : termes diffusion |
---|
133 | real :: fraxw,fraxe,frays,frayn !< termes advection |
---|
134 | |
---|
135 | real,dimension(nx,ny) :: deltah ! dans calcul relax |
---|
136 | real :: delh ! dans calcul relax |
---|
137 | real :: testh ! dans calcul relax |
---|
138 | |
---|
139 | logical :: stopp |
---|
140 | integer :: ntour |
---|
141 | real :: reste |
---|
142 | |
---|
143 | |
---|
144 | integer :: it1,it2,jt1,jt2 ! pour des tests d'asymétrie |
---|
145 | |
---|
146 | |
---|
147 | if (itracebug.eq.1) call tracebug(' Entree dans routine resolution_diffusion') |
---|
148 | |
---|
149 | ! calcul de bdx et hdx |
---|
150 | hdx(:,:)=dx1*(vieuxH(:,:)-eoshift(vieuxH(:,:),shift=-1,boundary=0.,dim=1)) |
---|
151 | hdy(:,:)=dx1*(vieuxH(:,:)-eoshift(vieuxH(:,:),shift=-1,boundary=0.,dim=2)) |
---|
152 | bdx(:,:)=dx1*(B(:,:)-eoshift(B(:,:),shift=-1,boundary=0.,dim=1)) |
---|
153 | bdy(:,:)=dx1*(B(:,:)-eoshift(B(:,:),shift=-1,boundary=0.,dim=2)) |
---|
154 | |
---|
155 | |
---|
156 | ! initialisations (qui feront aussi les conditions aux limites) |
---|
157 | Arelax(:,:)=0. |
---|
158 | Brelax(:,:)=0. |
---|
159 | Drelax(:,:)=0. |
---|
160 | Erelax(:,:)=0. |
---|
161 | Crelax(:,:)=1. |
---|
162 | Frelax(:,:)=0. |
---|
163 | DeltaH(:,:)=0. |
---|
164 | |
---|
165 | |
---|
166 | ! schema spatial |
---|
167 | |
---|
168 | if (upwind.eq.1.) then !schema amont |
---|
169 | |
---|
170 | where (Advx(:,:).ge.0.) |
---|
171 | c_west(:,:)=1. |
---|
172 | c_east(:,:)=0. |
---|
173 | elsewhere |
---|
174 | c_west(:,:)=0. |
---|
175 | c_east(:,:)=1. |
---|
176 | end where |
---|
177 | |
---|
178 | where (Advy(:,:).ge.0.) |
---|
179 | c_south(:,:)=1. |
---|
180 | c_north(:,:)=0. |
---|
181 | elsewhere |
---|
182 | c_south(:,:)=0. |
---|
183 | c_north(:,:)=1. |
---|
184 | end where |
---|
185 | |
---|
186 | else if (upwind.lt.1.) then ! schema centre |
---|
187 | c_west(:,:)=0.5 |
---|
188 | c_east(:,:)=0.5 |
---|
189 | c_south(:,:)=0.5 |
---|
190 | c_north(:,:)=0.5 |
---|
191 | end if |
---|
192 | |
---|
193 | ! attribution des elements des diagonales |
---|
194 | do j=2,ny-1 |
---|
195 | do i=2,nx-1 |
---|
196 | |
---|
197 | ! sous diagonale en x |
---|
198 | arelax(i,j)=-omega*Dtdx2*Dfx(i,j) & ! partie diffusive en x |
---|
199 | -mu_adv*dtdx*c_west(i,j)*Advx(i,j) ! partie advective en x |
---|
200 | |
---|
201 | ! sur diagonale en x |
---|
202 | brelax(i,j)=-omega*Dtdx2*Dfx(i+1,j) & ! partie diffusive |
---|
203 | +mu_adv*dtdx*c_east(i+1,j)*Advx(i+1,j) ! partie advective |
---|
204 | |
---|
205 | ! sous diagonale en y |
---|
206 | drelax(i,j)=-omega*Dtdx2*Dfy(i,j) & ! partie diffusive en y |
---|
207 | -mu_adv*dtdx*c_south(i,j)*Advy(i,j) ! partie advective en y |
---|
208 | |
---|
209 | ! sur diagonale en y |
---|
210 | erelax(i,j)=-omega*Dtdx2*Dfy(i,j+1) & ! partie diffusive |
---|
211 | +mu_adv*dtdx*c_north(i,j+1)*Advy(i,j+1) ! partie advective |
---|
212 | |
---|
213 | |
---|
214 | |
---|
215 | ! diagonale |
---|
216 | crelax(i,j)=omega*Dtdx2*((Dfx(i+1,j)+Dfx(i,j)) & ! partie diffusive en x |
---|
217 | +(Dfy(i,j+1)+Dfy(i,j))) ! partie diffusive en y |
---|
218 | crelax(i,j)=crelax(i,j)+mu_adv*dtdx* & |
---|
219 | ( (c_west(i+1,j)*Advx(i+1,j)-c_east(i,j)*Advx(i,j)) & !partie advective en x |
---|
220 | +(c_south(i,j+1)*Advy(i,j+1)-c_north(i,j)*Advy(i,j))) !partie advective en y |
---|
221 | crelax(i,j)=1.+crelax(i,j) ! partie temporelle |
---|
222 | |
---|
223 | |
---|
224 | ! terme du vecteur |
---|
225 | |
---|
226 | frdx= -Dfx(i,j) * (Bdx(i,j) +(1.-omega)*Hdx(i,j)) & ! partie diffusive en x |
---|
227 | +Dfx(i+1,j)*(Bdx(i+1,j)+(1.-omega)*Hdx(i+1,j)) |
---|
228 | |
---|
229 | frdy= -Dfy(i,j) * (Bdy(i,j) +(1.-omega)*Hdy(i,j)) & ! partie diffusive en y |
---|
230 | +Dfy(i,j+1)*(Bdy(i,j+1)+(1.-omega)*Hdy(i,j+1)) |
---|
231 | |
---|
232 | fraxw= -c_west(i,j)* Advx(i,j) * vieuxH(i-1,j) & ! partie advective en x |
---|
233 | +c_west(i+1,j)*Advx(i+1,j)*vieuxH(i,j) ! venant de l'west |
---|
234 | |
---|
235 | fraxe= -c_east(i,j) * Advx(i,j) * vieuxH(i,j) & ! partie advective en x |
---|
236 | +c_east(i+1,j)*Advx(i+1,j)*vieuxH(i+1,j) ! venant de l'est |
---|
237 | |
---|
238 | frays= -c_south(i,j) * Advy(i,j) * vieuxH(i,j-1) & ! partie advective en y |
---|
239 | +c_south(i,j+1)*Advy(i,j+1)*vieuxH(i,j) ! venant du sud |
---|
240 | |
---|
241 | frayn= -c_north(i,j) * Advy(i,j) * vieuxH(i,j) & ! partie advective en y |
---|
242 | +c_north(i,j+1)*Advy(i,j+1)*vieuxH(i,j+1) ! venant du nord |
---|
243 | |
---|
244 | |
---|
245 | |
---|
246 | frelax(i,j) = vieuxH(i,j) + Dt*bil(i,j) & |
---|
247 | + dtdx*(frdx+frdy) & ! termes lies a la diffusion (socle, omega) |
---|
248 | + (1.-mu_adv)*dtdx*((fraxw+fraxe)+(frays+frayn)) |
---|
249 | ! le dernier terme serait pour |
---|
250 | ! de l'over implicite en diffusion |
---|
251 | ! (en pratique pas utilise) |
---|
252 | end do |
---|
253 | end do |
---|
254 | |
---|
255 | |
---|
256 | where (i_hpresc(:,:) .eq.1) ! thickness prescribed |
---|
257 | frelax(:,:) = H_presc(:,:) |
---|
258 | arelax(:,:) = 0. |
---|
259 | brelax(:,:) = 0. |
---|
260 | crelax(:,:) = 1. |
---|
261 | drelax(:,:) = 0. |
---|
262 | erelax(:,:) = 0. |
---|
263 | end where |
---|
264 | |
---|
265 | |
---|
266 | ! conditions aux limites bords paralleles |
---|
267 | |
---|
268 | if (paral_x .eq.1) then ! bords paralleles a x : pas de termes selon y |
---|
269 | |
---|
270 | do i = 2,nx-1 |
---|
271 | |
---|
272 | ! sous diagonale en y |
---|
273 | drelax(i,2) = 0. |
---|
274 | drelax(i,ny-1) = 0. |
---|
275 | |
---|
276 | ! sur diagonale en y |
---|
277 | erelax(i,2) = 0. |
---|
278 | erelax(i,ny-1) = 0. |
---|
279 | |
---|
280 | ! diagonale |
---|
281 | j = 2 |
---|
282 | crelax(i,j)=omega*Dtdx2*((Dfx(i+1,j)+Dfx(i,j)) & ! partie diffusive en x |
---|
283 | |
---|
284 | crelax(i,j)=crelax(i,j)+mu_adv*dtdx* & |
---|
285 | ( (c_west(i+1,j)*Advx(i+1,j)-c_east(i,j)*Advx(i,j)) & ! partie advective en x |
---|
286 | |
---|
287 | crelax(i,j)=1.+crelax(i,j) ! partie temporelle |
---|
288 | |
---|
289 | |
---|
290 | j = ny-1 |
---|
291 | crelax(i,j)=omega*Dtdx2*((Dfx(i+1,j)+Dfx(i,j)) & ! partie diffusive en x |
---|
292 | crelax(i,j)=crelax(i,j)+mu_adv*dtdx* & |
---|
293 | ( (c_west(i+1,j)*Advx(i+1,j)-c_east(i,j)*Advx(i,j)) & ! partie advective en x |
---|
294 | |
---|
295 | crelax(i,j)=1.+crelax(i,j) ! partie temporelle |
---|
296 | |
---|
297 | end do |
---|
298 | |
---|
299 | |
---|
300 | else if (paral_y .eq.1) then ! bords paralleles a y : pas de termes selon x |
---|
301 | |
---|
302 | do j = 2,ny-1 |
---|
303 | |
---|
304 | ! sous diagonale en x |
---|
305 | arelax(2,j) = 0 |
---|
306 | |
---|
307 | ! sur diagonale en x |
---|
308 | brelax(nx-1,j) = 0. |
---|
309 | |
---|
310 | |
---|
311 | ! diagonale |
---|
312 | i = 2 |
---|
313 | crelax(i,j)=omega*Dtdx2*(Dfy(i,j+1)+Dfy(i,j)) ! partie diffusive en y |
---|
314 | crelax(i,j)=crelax(i,j)+mu_adv*dtdx* & |
---|
315 | (c_south(i,j+1)*Advy(i,j+1)-c_north(i,j)*Advy(i,j)) !partie advective en y |
---|
316 | crelax(i,j)=1.+crelax(i,j) ! partie temporelle |
---|
317 | |
---|
318 | i = nx-1 |
---|
319 | crelax(i,j)=omega*Dtdx2*(Dfy(i,j+1)+Dfy(i,j)) ! partie diffusive en y |
---|
320 | crelax(i,j)=crelax(i,j)+mu_adv*dtdx* & |
---|
321 | (c_south(i,j+1)*Advy(i,j+1)-c_north(i,j)*Advy(i,j)) !partie advective en y |
---|
322 | crelax(i,j)=1.+crelax(i,j) ! partie temporelle |
---|
323 | |
---|
324 | end do |
---|
325 | end if |
---|
326 | |
---|
327 | stopp = .false. |
---|
328 | ntour=0 |
---|
329 | |
---|
330 | |
---|
331 | |
---|
332 | relax_loop: do while(.not.stopp) |
---|
333 | ntour=ntour+1 |
---|
334 | |
---|
335 | do j=2,ny-1 |
---|
336 | do i=2,nx-1 |
---|
337 | |
---|
338 | reste = (((arelax(i,j)*newH(i-1,j) +drelax(i,j)*newH(i,j-1)) & |
---|
339 | + (brelax(i,j)*newH(i+1,j) + erelax(i,j)*newH(i,j+1))) & |
---|
340 | + crelax(i,j)*newH(i,j))- frelax(i,j) |
---|
341 | |
---|
342 | if (ntour.eq.1) debug_3D(i,j,49)=(((arelax(i,j)*newH(i-1,j) +drelax(i,j)*newH(i,j-1)) & |
---|
343 | + (brelax(i,j)*newH(i+1,j) + erelax(i,j)*newH(i,j+1))) & |
---|
344 | + crelax(i,j)*newH(i,j)) |
---|
345 | |
---|
346 | |
---|
347 | deltaH(i,j) = reste/crelax(i,j) |
---|
348 | |
---|
349 | end do |
---|
350 | end do |
---|
351 | |
---|
352 | debug_3D(:,:,50)=arelax(:,:) |
---|
353 | debug_3D(:,:,51)=brelax(:,:) |
---|
354 | debug_3D(:,:,52)=crelax(:,:) |
---|
355 | debug_3D(:,:,53)=drelax(:,:) |
---|
356 | debug_3D(:,:,54)=erelax(:,:) |
---|
357 | debug_3D(:,:,55)=frelax(:,:) |
---|
358 | |
---|
359 | |
---|
360 | |
---|
361 | newH(:,:)=newH(:,:)-deltaH(:,:) |
---|
362 | |
---|
363 | |
---|
364 | ! critere d'arret: |
---|
365 | ! ---------------- |
---|
366 | |
---|
367 | delh=0 |
---|
368 | |
---|
369 | |
---|
370 | do j=2,ny-1 |
---|
371 | do i=2,nx-1 |
---|
372 | delh=delh+deltaH(i,j)**2 |
---|
373 | end do |
---|
374 | end do |
---|
375 | |
---|
376 | if (delh.gt.0.) then |
---|
377 | testh=sqrt(delh)/((nx-2)*(ny-2)) |
---|
378 | else |
---|
379 | testh=0. |
---|
380 | endif |
---|
381 | stopp = (testh.lt.1.e-4).or.(ntour.gt.100) |
---|
382 | |
---|
383 | |
---|
384 | end do relax_loop |
---|
385 | |
---|
386 | |
---|
387 | ! thickness at the upwind node |
---|
388 | do j = 1, ny |
---|
389 | do i = 2, nx |
---|
390 | debug_3D(i,j,92) = c_west(i,j) * newH(i-1,j) + c_east(i,j) * newH(i,j) |
---|
391 | end do |
---|
392 | end do |
---|
393 | do j = 2, ny |
---|
394 | do i = 1, nx |
---|
395 | debug_3D(i,j,93) = c_south(i,j) * newH(i,j-1) + c_north(i,j) * newH(i,j) |
---|
396 | end do |
---|
397 | end do |
---|
398 | |
---|
399 | |
---|
400 | return |
---|
401 | |
---|
402 | end subroutine resol_adv_diff_2D_vect |
---|
403 | |
---|
404 | end module reso_adv_diff_2D_paral |
---|
405 | |
---|
406 | |
---|
407 | |
---|
408 | |
---|