1 | !module diagno_L2_mod ! Nouvelle version, compatible remplimat 2008 Cat |
---|
2 | module diagno_mod ! nom pendant les tests |
---|
3 | !$ USE OMP_LIB |
---|
4 | use module3D_phy |
---|
5 | use module_choix |
---|
6 | |
---|
7 | implicit none |
---|
8 | |
---|
9 | |
---|
10 | |
---|
11 | real :: somint,test,delp,prec |
---|
12 | real, dimension(nx,ny) :: uxb1 |
---|
13 | real, dimension(nx,ny) :: uyb1 |
---|
14 | |
---|
15 | real, dimension(nx,ny) :: uxb1ramollo |
---|
16 | real, dimension(nx,ny) :: uyb1ramollo |
---|
17 | real, dimension(nx,ny) :: pvi_keep |
---|
18 | |
---|
19 | !cdc transfere dans module3d pour compatibilite avec furst_schoof_mod |
---|
20 | !cdc integer, dimension(nx,ny) :: imx_diag |
---|
21 | !cdc integer, dimension(nx,ny) :: imy_diag |
---|
22 | |
---|
23 | integer :: nneigh=max(int(400000./dx),1) !nb. points around the grline kept for back force comp. in case we reduce the extent |
---|
24 | integer, dimension(nx,ny) :: imx_reduce !afq -- to limit the number of points where vel is computed for backforce |
---|
25 | integer, dimension(nx,ny) :: imy_reduce !afq -- to limit the number of points where vel is computed for backforce |
---|
26 | |
---|
27 | integer :: nxd1,nxd2 ! domaine selon x Dans l'appel rempli_L2 |
---|
28 | integer :: nyd1,nyd2 ! domaine selon y |
---|
29 | |
---|
30 | integer :: itour_pvi |
---|
31 | |
---|
32 | integer :: ifail_diagno ! pour recuperation d'erreur |
---|
33 | integer :: ifail_diagno_ramollo ! pour recuperation d'erreur shelf ramollo |
---|
34 | integer :: iplus1,jplus1 |
---|
35 | integer :: ctvisco,iumax,jumax |
---|
36 | real :: delumax,errmax |
---|
37 | real :: phiphi,bt2,d02,discr,ttau |
---|
38 | real :: sf3,sf1,epsxxm,epsyym,epsm,sf01,sf03 ! pour le calcul de la viscosite |
---|
39 | real :: viscm |
---|
40 | real :: sf_shelf ! coef mult enhancement factor pour shelves |
---|
41 | |
---|
42 | logical :: stopvisco,viscolin |
---|
43 | logical :: test_visc |
---|
44 | |
---|
45 | contains |
---|
46 | |
---|
47 | !------------------------------------------------------------------------------------ |
---|
48 | subroutine init_diagno |
---|
49 | |
---|
50 | namelist/diagno_rheol/sf01,sf03,pvimin |
---|
51 | |
---|
52 | ! attribution des coefficients de viscosite |
---|
53 | |
---|
54 | ! formats pour les ecritures dans 42 |
---|
55 | 428 format(A) |
---|
56 | |
---|
57 | ! lecture des parametres du run block draghwat |
---|
58 | !-------------------------------------------------------------------- |
---|
59 | |
---|
60 | rewind(num_param) ! pour revenir au debut du fichier param_list.dat |
---|
61 | read(num_param,diagno_rheol) |
---|
62 | |
---|
63 | write(num_rep_42,428)'!___________________________________________________________' |
---|
64 | write(num_rep_42,428) '&diagno_rheol ! nom du bloc diagno_rheol' |
---|
65 | write(num_rep_42,*) |
---|
66 | write(num_rep_42,*) 'sf01 = ',sf01 |
---|
67 | write(num_rep_42,*) 'sf03 = ',sf03 |
---|
68 | write(num_rep_42,*) 'pvimin = ',pvimin |
---|
69 | write(num_rep_42,*)'/' |
---|
70 | write(num_rep_42,428) '! coefficients par rapport a la loi glace posee ' |
---|
71 | write(num_rep_42,428) '! sf01 : coefficient viscosite loi lineaire ' |
---|
72 | write(num_rep_42,428) '! sf03 : coefficient viscosite loi n=3 ' |
---|
73 | write(num_rep_42,428) '! pvimin : valeur de pvi pour les noeuds fictifs ~ 1.e3' |
---|
74 | write(num_rep_42,428) '! tres petit par rapport aux valeurs standards ~ 1.e10' |
---|
75 | |
---|
76 | write(num_rep_42,*) |
---|
77 | |
---|
78 | |
---|
79 | ! Precision utilisee dans de calcul |
---|
80 | prec = 1.e-2 |
---|
81 | itour_pvi=1 ! si prend les valeurs analytiques dans le shelf |
---|
82 | |
---|
83 | if (geoplace(1:5).eq.'mism3') then |
---|
84 | sf_shelf = 1. |
---|
85 | itour_pvi= 0 |
---|
86 | |
---|
87 | else |
---|
88 | sf_shelf = 0.4 |
---|
89 | end if |
---|
90 | |
---|
91 | return |
---|
92 | end subroutine init_diagno |
---|
93 | |
---|
94 | !------------------------------------------------------------------------------------ |
---|
95 | subroutine diagnoshelf ! Resolution numerique des equations diagnostiques |
---|
96 | |
---|
97 | |
---|
98 | integer :: diagno_grline |
---|
99 | |
---|
100 | |
---|
101 | if (itracebug.eq.1) call tracebug(' Entree dans diagnoshelf') |
---|
102 | |
---|
103 | |
---|
104 | itour_pvi=itour_pvi+1 ! boucle sur la viscosite (pour l'instant pas actif) |
---|
105 | |
---|
106 | ! pvi(:,:)=0. |
---|
107 | Taushelf(:,:)=0. |
---|
108 | |
---|
109 | ! attention le bloc suivant est pour debug |
---|
110 | !!$gzmx(:,:)=.false. |
---|
111 | !!$gzmy(:,:)=.false. |
---|
112 | !!$ilemx(:,:)=.false. |
---|
113 | !!$ilemy(:,:)=.false. |
---|
114 | !!$flgzmx(:,:)=flotmx(:,:) |
---|
115 | !!$flgzmy(:,:)=flotmy(:,:) |
---|
116 | |
---|
117 | call dragging ! doit etre appele avant imx_imy |
---|
118 | |
---|
119 | |
---|
120 | |
---|
121 | if (itour_pvi.le.1) then |
---|
122 | call calc_pvi ! calcule les viscosites integrees |
---|
123 | |
---|
124 | !$OMP PARALLEL |
---|
125 | !$OMP WORKSHARE |
---|
126 | where (flot(:,:).and.(H.gt.1)) ! valeur analytique pour les shelfs |
---|
127 | pvi(:,:) = (4./coef_Sflot/rog)**2/btt(:,:,1,1)/H(:,:) |
---|
128 | end where |
---|
129 | !$OMP END WORKSHARE |
---|
130 | !$OMP END PARALLEL |
---|
131 | ! avec couplage thermomecanique |
---|
132 | ! write(166,*) ' apres call calc_pvi',itour_pvi |
---|
133 | |
---|
134 | else |
---|
135 | call calc_pvi |
---|
136 | end if |
---|
137 | |
---|
138 | call imx_imy_nx_ny ! pour rempli_L2 : calcule les masques imx et imy qui |
---|
139 | |
---|
140 | diagno_grline=sum(gr_line(2:(nx-1),2:(ny-1))) |
---|
141 | if ((Schoof.ge.1).and.(diagno_grline.gt.0)) then |
---|
142 | call imx_imy_nx_ny_reduce(1) !afq: provides imx_reduce and imy_reduce |
---|
143 | endif |
---|
144 | |
---|
145 | !cdc debug Schoof !!!!!!!!!!!! |
---|
146 | !~ do j=1,ny |
---|
147 | !~ do i=1,nx |
---|
148 | !~ write(578,*) uxbar(i,j) |
---|
149 | !~ write(579,*) uybar(i,j) |
---|
150 | !~ enddo |
---|
151 | !~ enddo |
---|
152 | |
---|
153 | !if (Schoof.eq.1.and.nt.GT.15000) then ! flux grounding line Schoof |
---|
154 | ! afq -- below: if (Schoof.eq.1) then ! flux grounding line Schoof |
---|
155 | ! afq -- below: call interpol_glflux ! calcul flux GL + interpolation sur voisins |
---|
156 | ! afq -- below: endif |
---|
157 | |
---|
158 | !~ do j=1,ny |
---|
159 | !~ do i=1,nx |
---|
160 | !~ write(588,*) uxbar(i,j) |
---|
161 | !~ write(589,*) uybar(i,j) |
---|
162 | !~ enddo |
---|
163 | !~ enddo |
---|
164 | !~ print*,'ecriteure termineee !!!!!!' |
---|
165 | !~ read(*,*) |
---|
166 | |
---|
167 | ! donnent les cas de conditions aux limites |
---|
168 | ! |
---|
169 | ! version pour travailler sur tout le domaine nx ny |
---|
170 | |
---|
171 | if (geoplace(1:5).eq.'mism3') call mismip_boundary_cond |
---|
172 | |
---|
173 | |
---|
174 | ! appel a la routine rempl_L2 -------------------domaine nx x ny ------------ |
---|
175 | ! |
---|
176 | |
---|
177 | ! pour tout le domaine |
---|
178 | nxd1=1 |
---|
179 | nxd2=nx |
---|
180 | nyd1=1 |
---|
181 | nyd2=ny |
---|
182 | |
---|
183 | !call rempli_L2(1,nx,1,ny,uxbar,uybar,uxb1,uyb1,imx_diag,imy_diag,ifail_diagno) |
---|
184 | !nxd1=15 |
---|
185 | !nxd2=19 |
---|
186 | !nyd1=30 |
---|
187 | !nyd2=34 |
---|
188 | |
---|
189 | !nxd1=35 |
---|
190 | !nxd2=60 |
---|
191 | !nyd1=35 |
---|
192 | !nyd2=60 |
---|
193 | |
---|
194 | |
---|
195 | !if (Schoof.eq.1.and.nt.GT.15000) then ! flux grounding line Schoof avec calcul de la back force par shelf ramollo |
---|
196 | if ((Schoof.ge.1).and.(diagno_grline.gt.0)) then ! flux grounding line Schoof avec calcul de la back force par shelf ramollo |
---|
197 | |
---|
198 | call rempli_L2(nxd1,nxd2,nyd1,nyd2,uxbar(nxd1:nxd2,nyd1:nyd2),uybar(nxd1:nxd2,nyd1:nyd2), & |
---|
199 | uxb1(nxd1:nxd2,nyd1:nyd2),uyb1(nxd1:nxd2,nyd1:nyd2), & |
---|
200 | imx_reduce(nxd1:nxd2,nyd1:nyd2),imy_reduce(nxd1:nxd2,nyd1:nyd2),ifail_diagno) |
---|
201 | |
---|
202 | pvi_keep(:,:)=pvi(:,:) |
---|
203 | where (flot(:,:).and.H(:,:).GT.2.) |
---|
204 | ! pvi(:,:)=1.e5 |
---|
205 | pvi(:,:)=pvimin |
---|
206 | endwhere |
---|
207 | |
---|
208 | call rempli_L2(nxd1,nxd2,nyd1,nyd2,uxbar(nxd1:nxd2,nyd1:nyd2),uybar(nxd1:nxd2,nyd1:nyd2), & |
---|
209 | uxb1ramollo(nxd1:nxd2,nyd1:nyd2),uyb1ramollo(nxd1:nxd2,nyd1:nyd2), & |
---|
210 | imx_reduce(nxd1:nxd2,nyd1:nyd2),imy_reduce(nxd1:nxd2,nyd1:nyd2),ifail_diagno_ramollo) |
---|
211 | |
---|
212 | pvi(:,:)=pvi_keep(:,:) |
---|
213 | |
---|
214 | where (abs(uxb1ramollo(:,:)) .GT.1.e-5) |
---|
215 | back_force_x(:,:) = 1.0 * abs(uxb1(:,:)) / abs(uxb1ramollo(:,:)) |
---|
216 | elsewhere |
---|
217 | back_force_x(:,:)=1. |
---|
218 | endwhere |
---|
219 | where (abs(uyb1ramollo(:,:)) .GT.1.e-5) |
---|
220 | back_force_y(:,:) = 1.0 * abs(uyb1(:,:)) / abs(uyb1ramollo(:,:)) |
---|
221 | elsewhere |
---|
222 | back_force_y(:,:)=1. |
---|
223 | endwhere |
---|
224 | back_force_x(:,:) = min ( back_force_x(:,:), 1. ) |
---|
225 | back_force_y(:,:) = min ( back_force_y(:,:), 1. ) |
---|
226 | debug_3D(:,:,64) = back_force_x(:,:) |
---|
227 | debug_3D(:,:,65) = back_force_y(:,:) |
---|
228 | |
---|
229 | if (ifail_diagno_ramollo.gt.0) then |
---|
230 | ! write(6,*) ' Probleme resolution systeme L2. ramollo ifail=',ifail_diagno_ramollo |
---|
231 | ! STOP |
---|
232 | write(*,*) ' Probleme resolution systeme L2. ramollo ifail=',ifail_diagno_ramollo |
---|
233 | write(*,*) ' ... we go on anyway!' |
---|
234 | endif |
---|
235 | !~ do j=1,ny |
---|
236 | !~ do i=1,nx |
---|
237 | !~ if (sqrt(uxb1(i,j)**2+ uyb1(i,j)**2).gt.0..and..not.flot(i,j)) then |
---|
238 | !~ write(1034,*) sqrt(uxb1(i,j)**2+ uyb1(i,j)**2) / sqrt(uxb1ramollo(i,j)**2 + uyb1ramollo(i,j)**2) |
---|
239 | !~ else |
---|
240 | !~ write(1034,*) 1. |
---|
241 | !~ endif |
---|
242 | !~ enddo |
---|
243 | !~ enddo |
---|
244 | |
---|
245 | !~ print*,'apres calcul rempli_L2' |
---|
246 | !~ read(*,*) |
---|
247 | |
---|
248 | call interpol_glflux ! calcul flux GL + interpolation sur voisins |
---|
249 | |
---|
250 | endif |
---|
251 | |
---|
252 | call rempli_L2(nxd1,nxd2,nyd1,nyd2,uxbar(nxd1:nxd2,nyd1:nyd2),uybar(nxd1:nxd2,nyd1:nyd2), & |
---|
253 | uxb1(nxd1:nxd2,nyd1:nyd2),uyb1(nxd1:nxd2,nyd1:nyd2), & |
---|
254 | imx_diag(nxd1:nxd2,nyd1:nyd2),imy_diag(nxd1:nxd2,nyd1:nyd2),ifail_diagno) |
---|
255 | |
---|
256 | |
---|
257 | ! Dans rempli_L2 |
---|
258 | |
---|
259 | ! uxprex(n1,n2) |
---|
260 | ! uyprec(n1,n2) vitesses de l'iteration precedente |
---|
261 | ! |
---|
262 | ! uxnew(n1,n2) |
---|
263 | ! uynew(n1,n2) uynew resultat de cette iteration |
---|
264 | ! |
---|
265 | ! imx(n1,n2) masque pour imposer les vitesses ou leur dérivee |
---|
266 | ! imy(n1,n2) masque pour imposer les vitesses ou leur dérivée |
---|
267 | |
---|
268 | ! eventuellement le domaine n1,n2 peut etre un sous-domaine de nx,ny |
---|
269 | ! attention il faudra alors appeler avec des sous-tableaux |
---|
270 | |
---|
271 | ! Dans l'appel uxbar -> uxprec et Uxb1 -> Uxnew |
---|
272 | |
---|
273 | !--------------------------------------------------------------------------- |
---|
274 | |
---|
275 | if (ifail_diagno.gt.0) then |
---|
276 | write(6,*) ' Probleme resolution systeme L2. ifail=',ifail_diagno |
---|
277 | STOP |
---|
278 | endif |
---|
279 | |
---|
280 | ! nouvelles vitesses |
---|
281 | !$OMP PARALLEL |
---|
282 | !$OMP WORKSHARE |
---|
283 | uxbar(:,:)=max(min(uxb1(:,:),V_limit),-V_limit) |
---|
284 | uybar(:,:)=max(min(uyb1(:,:),V_limit),-V_limit) |
---|
285 | !$OMP END WORKSHARE |
---|
286 | |
---|
287 | |
---|
288 | ! calcul de tobmx et tobmy (frottement basal) apres calcul des vitesses |
---|
289 | ! --------------------------------------------------------------------- |
---|
290 | |
---|
291 | if (gr_select.eq.1) then !Tsai |
---|
292 | !$OMP WORKSHARE |
---|
293 | where(gr_line_schoof(:,:).gt.0) |
---|
294 | tobmx(:,:) = sign( neffmx(:,:) * frot_coef , -uxbar(:,:)) |
---|
295 | tobmy(:,:) = sign( neffmy(:,:) * frot_coef , -uybar(:,:)) |
---|
296 | elsewhere |
---|
297 | tobmx(:,:) = -betamx(:,:) * uxbar(:,:) |
---|
298 | tobmy(:,:) = -betamy(:,:) * uybar(:,:) |
---|
299 | end where |
---|
300 | !$OMP END WORKSHARE |
---|
301 | else |
---|
302 | !$OMP WORKSHARE |
---|
303 | tobmx(:,:) = -betamx(:,:) * uxbar(:,:) |
---|
304 | tobmy(:,:) = -betamy(:,:) * uybar(:,:) |
---|
305 | !$OMP END WORKSHARE |
---|
306 | end if |
---|
307 | !$OMP BARRIER |
---|
308 | |
---|
309 | ! afq -- initMIP outputs: |
---|
310 | !$OMP DO |
---|
311 | do j=2,ny-1 |
---|
312 | do i=2,nx-1 |
---|
313 | taub(i,j) = sqrt( ((tobmx(i,j)+tobmx(i+1,j))*0.5)**2 & |
---|
314 | + ((tobmy(i,j)+tobmy(i,j+1))*0.5)**2 ) |
---|
315 | end do |
---|
316 | end do |
---|
317 | !$OMP END DO |
---|
318 | !$OMP WORKSHARE |
---|
319 | debug_3d(:,:,117) = taub(:,:) |
---|
320 | !$OMP END WORKSHARE |
---|
321 | |
---|
322 | |
---|
323 | ! Mise ne réserve des vitesses flgzmx et flgzmy |
---|
324 | !$OMP WORKSHARE |
---|
325 | where (flgzmx(:,:)) |
---|
326 | uxflgz(:,:)=uxbar(:,:) |
---|
327 | elsewhere |
---|
328 | uxflgz(:,:)=0. |
---|
329 | endwhere |
---|
330 | |
---|
331 | where (flgzmy(:,:)) |
---|
332 | uyflgz(:,:)=uybar(:,:) |
---|
333 | elsewhere |
---|
334 | uyflgz(:,:)=0. |
---|
335 | endwhere |
---|
336 | !$OMP END WORKSHARE |
---|
337 | !$OMP END PARALLEL |
---|
338 | !i=92 |
---|
339 | !j=152 |
---|
340 | !write(6,*) 'time',time, uxbar(92,152),gzmx(92,152),ilemx(92,152),flotmx(92,152), flgzmx(92,152) |
---|
341 | |
---|
342 | return |
---|
343 | end subroutine diagnoshelf |
---|
344 | |
---|
345 | |
---|
346 | !------------------------------------------------------------------- |
---|
347 | subroutine calc_pvi |
---|
348 | |
---|
349 | ! calcule les viscosites integrees pvi et pvm |
---|
350 | ! loi polynomiale + couplage thermomécanique |
---|
351 | ! |
---|
352 | ! Attention ne marche que si la loi est la loi en n=3 + n=1 |
---|
353 | ! y compris le pur glen (n=3) ou le pur Newtonien (n=1) |
---|
354 | ! -------------------------------------------------------------------- |
---|
355 | |
---|
356 | ! les deformations sont supposées indépendantes de la profondeur |
---|
357 | ! et sont calculées dans strain_rate (appelé par main) |
---|
358 | |
---|
359 | ! eps(i,j) -> eps |
---|
360 | ! ttau -> tau (2eme invariant du deviateur des contraintes) |
---|
361 | ! BT2 loi en n=3, phiphi loi en n=1 calculés dans flowlaw |
---|
362 | |
---|
363 | ! |
---|
364 | |
---|
365 | |
---|
366 | ! La viscosité est calculée partout y compris pour la glace posée (ou elle est moins |
---|
367 | ! précise qu'un calcul direct avec la loi de déformation.) |
---|
368 | ! Pour les noeuds posés mais ayant un voisin stream ou flottant, on calcule |
---|
369 | ! la viscosité avec stream/shelves |
---|
370 | ! le calcul se fait sur les noeuds majeurs |
---|
371 | |
---|
372 | !$ integer :: rang ,nb_taches |
---|
373 | !$ logical :: paral |
---|
374 | |
---|
375 | if (itracebug.eq.1) call tracebug(' Calc pvi') |
---|
376 | |
---|
377 | !$OMP PARALLEL PRIVATE(rang,iplus1,jplus1,sf3,sf1,BT2,phiphi,ttau,d02,discr) |
---|
378 | !$ paral = OMP_IN_PARALLEL() |
---|
379 | !$ rang=OMP_GET_THREAD_NUM() |
---|
380 | !$ nb_taches=OMP_GET_NUM_THREADS() |
---|
381 | |
---|
382 | !$OMP WORKSHARE |
---|
383 | pvi(:,:) = pvimin |
---|
384 | Abar(:,:) = 0. |
---|
385 | !$OMP END WORKSHARE |
---|
386 | |
---|
387 | !$OMP DO |
---|
388 | do j=1,ny |
---|
389 | do i=1,nx |
---|
390 | iplus1=min(i+1,nx) |
---|
391 | jplus1=min(j+1,ny) |
---|
392 | |
---|
393 | |
---|
394 | if (flot(i,j)) then ! noeuds flottants |
---|
395 | sf3=sf03*sf_shelf |
---|
396 | sf1=sf01*sf_shelf |
---|
397 | |
---|
398 | else if (gzmx(i,j).or.gzmx(iplus1,j).or.gzmy(i,j).or.gzmy(i,jplus1)) then |
---|
399 | sf1=sf01 |
---|
400 | sf3=max(sf03,0.01) ! pour les fleuves de glace, un peu de Glen |
---|
401 | |
---|
402 | else if (ilemx(i,j).or.ilemx(iplus1,j).or.ilemy(i,j).or.ilemy(i,jplus1)) then |
---|
403 | sf1=sf01 |
---|
404 | sf3=max(sf03,0.01) ! pour les iles aussi |
---|
405 | |
---|
406 | else |
---|
407 | ! sf1=1 |
---|
408 | ! sf3=1 |
---|
409 | sf1=sf01 ! pour la viscosite anisotrope (ici en longitudinal) |
---|
410 | sf3=sf03 |
---|
411 | endif |
---|
412 | |
---|
413 | |
---|
414 | do k=1,nz |
---|
415 | |
---|
416 | BT2=BTT(i,j,k,1)*sf3 ! changement du sf |
---|
417 | phiphi=BTT(i,j,k,2)*sf1 ! changement du sf |
---|
418 | |
---|
419 | |
---|
420 | |
---|
421 | if (BT2.lt.1.e-25) then ! pur newtonien |
---|
422 | visc(i,j,k)=1./phiphi |
---|
423 | ttau=2.*visc(i,j,k)*eps(i,j) |
---|
424 | else ! polynomial |
---|
425 | |
---|
426 | |
---|
427 | ! en mettant Bt2 en facteur |
---|
428 | d02=eps(i,j) |
---|
429 | discr=((phiphi/3.)**3.)/Bt2+d02**2 |
---|
430 | discr=discr**0.5 |
---|
431 | |
---|
432 | ttau=(d02+discr)**(1./3.)-(discr-d02)**(1./3.) |
---|
433 | |
---|
434 | |
---|
435 | ttau=ttau*Bt2**(-1./3.) |
---|
436 | |
---|
437 | |
---|
438 | visc(i,j,k)=Bt2*ttau*ttau+phiphi |
---|
439 | |
---|
440 | if (visc(i,j,k).gt.1.e-15) then |
---|
441 | visc(i,j,k)=1./visc(i,j,k) |
---|
442 | else |
---|
443 | visc(i,j,k)=1.e15 |
---|
444 | endif |
---|
445 | endif |
---|
446 | pvi(i,j)=pvi(i,j)+visc(i,j,k) |
---|
447 | Abar(i,j) =(Bt2/2.)**(-1./3.) + Abar(i,j) |
---|
448 | |
---|
449 | Taushelf(i,j)=Taushelf(i,j)+ttau |
---|
450 | |
---|
451 | |
---|
452 | end do |
---|
453 | |
---|
454 | |
---|
455 | |
---|
456 | pvi(i,j) = pvi(i,j)*H(i,j)/nz |
---|
457 | Abar(i,j) = (Abar(i,j) /nz)**(-3.) |
---|
458 | |
---|
459 | |
---|
460 | Taushelf(i,j)=Taushelf(i,j)/nz |
---|
461 | |
---|
462 | |
---|
463 | end do |
---|
464 | end do |
---|
465 | !$OMP END DO |
---|
466 | |
---|
467 | ! cas des noeuds fictifs, si l'épaisseur est très petite |
---|
468 | ! pvimin est très petit |
---|
469 | !$OMP WORKSHARE |
---|
470 | where (H(:,:).le.1.) |
---|
471 | pvi(:,:) = pvimin |
---|
472 | end where |
---|
473 | |
---|
474 | where (ramollo(:,:).ge..5) |
---|
475 | pvi(:,:) = pvimin |
---|
476 | end where |
---|
477 | |
---|
478 | |
---|
479 | debug_3D(:,:,27)=pvi(:,:) |
---|
480 | !$OMP END WORKSHARE |
---|
481 | ! attention run 35 |
---|
482 | !-------------------- |
---|
483 | !!$if (time.gt.10.) then |
---|
484 | !!$ where (flot(:,:)) |
---|
485 | !!$ pvi(:,:)=pvimin |
---|
486 | !!$ end where |
---|
487 | !!$end if |
---|
488 | |
---|
489 | ! calcul de la viscosite integree au milieu des mailles (pvm) |
---|
490 | !$OMP DO |
---|
491 | do i=2,nx |
---|
492 | do j=2,ny |
---|
493 | |
---|
494 | ! les lignes suivantes pour un pvm moyenne des pvi |
---|
495 | pvm(i,j)=0.25*((pvi(i,j)+pvi(i-1,j-1))+ & |
---|
496 | (pvi(i,j-1)+pvi(i-1,j))) |
---|
497 | |
---|
498 | end do |
---|
499 | end do |
---|
500 | !$OMP END DO |
---|
501 | !$OMP END PARALLEL |
---|
502 | |
---|
503 | end subroutine calc_pvi |
---|
504 | !------------------------------------------------------------------ |
---|
505 | |
---|
506 | subroutine imx_imy_nx_ny |
---|
507 | |
---|
508 | ! definition des masques |
---|
509 | ! pour rempli_L2 : calcule les masques imx et imy qui |
---|
510 | ! donnent les cas de conditions aux limites |
---|
511 | ! version pour travailler sur tout le domaine nx ny |
---|
512 | !---------------------------------------------------- |
---|
513 | |
---|
514 | |
---|
515 | |
---|
516 | ! -34 -3 Nord -23 |
---|
517 | ! !----------------------------------------! |
---|
518 | ! ! ! |
---|
519 | ! ! 1 (prescrite) ! |
---|
520 | ! -4 ! ou ! -2 |
---|
521 | ! West! 2 (L2) ! Est |
---|
522 | ! ! ! |
---|
523 | ! ! ! |
---|
524 | ! !----------------------------------------! |
---|
525 | ! -41 -1 Sud -12 |
---|
526 | |
---|
527 | !$OMP PARALLEL |
---|
528 | !$OMP WORKSHARE |
---|
529 | imx_diag(:,:)=0 |
---|
530 | imy_diag(:,:)=0 |
---|
531 | |
---|
532 | ! a l'interieur du domaine |
---|
533 | !------------------------- |
---|
534 | |
---|
535 | where (flgzmx(:,:)) |
---|
536 | imx_diag(:,:)=2 ! shelf ou stream |
---|
537 | elsewhere |
---|
538 | imx_diag(:,:)=1 ! vitesse imposee |
---|
539 | end where |
---|
540 | |
---|
541 | where (flgzmy(:,:)) |
---|
542 | imy_diag(:,:)=2 ! shelf ou stream |
---|
543 | elsewhere |
---|
544 | imy_diag(:,:)=1 ! vitesse imposee |
---|
545 | end where |
---|
546 | |
---|
547 | ! bord sud |
---|
548 | imx_diag(:,1)=-1 |
---|
549 | imy_diag(:,2)=-1 |
---|
550 | |
---|
551 | ! bord nord |
---|
552 | imx_diag(:,ny)=-3 |
---|
553 | imy_diag(:,ny)=-3 |
---|
554 | |
---|
555 | ! bord Est |
---|
556 | imx_diag(1,:)=0 ! hors domaine a cause des mailles alternees |
---|
557 | imx_diag(2,:)=-4 |
---|
558 | imy_diag(1,:)=-4 |
---|
559 | |
---|
560 | ! bord West |
---|
561 | imx_diag(nx,:)=-2 |
---|
562 | imy_diag(nx,:)=-2 |
---|
563 | |
---|
564 | ! Coins |
---|
565 | imx_diag(2,1)=-41 ! SW |
---|
566 | imy_diag(1,2)=-41 |
---|
567 | |
---|
568 | imx_diag(nx,1)=-12 ! SE |
---|
569 | imy_diag(nx,2)=-12 |
---|
570 | |
---|
571 | imx_diag(nx,ny)=-23 ! NE |
---|
572 | imy_diag(nx,ny)=-23 |
---|
573 | |
---|
574 | imx_diag(2,ny)=-34 ! NW |
---|
575 | imy_diag(1,ny)=-34 |
---|
576 | |
---|
577 | ! hors domaine |
---|
578 | imx_diag(1,:)=0 ! hors domaine a cause des mailles alternees |
---|
579 | imy_diag(:,1)=0 ! hors domaine a cause des mailles alternees |
---|
580 | !$OMP END WORKSHARE |
---|
581 | !$OMP END PARALLEL |
---|
582 | |
---|
583 | end subroutine imx_imy_nx_ny |
---|
584 | |
---|
585 | |
---|
586 | subroutine imx_imy_nx_ny_reduce(choix) |
---|
587 | |
---|
588 | !afq -- For the backforce computation we do not need to compute the velocities everywhere |
---|
589 | ! We simply compute the velocities around (nneigh) the grounding line |
---|
590 | |
---|
591 | integer, intent(in) :: choix |
---|
592 | |
---|
593 | integer i,j,nvx,nvy |
---|
594 | integer bordouest,bordest,bordsud,bordnord |
---|
595 | |
---|
596 | if (choix.eq.1) then ! we reduce the extent of velocity computation |
---|
597 | |
---|
598 | imx_reduce(:,:) = 1 |
---|
599 | imy_reduce(:,:) = 1 |
---|
600 | |
---|
601 | where (flot(:,:)) |
---|
602 | imx_reduce(:,:) = imx_diag(:,:) |
---|
603 | imy_reduce(:,:) = imy_diag(:,:) |
---|
604 | endwhere |
---|
605 | |
---|
606 | do j=1,ny |
---|
607 | do i=1,nx |
---|
608 | if (gr_line(i,j).eq.1) then |
---|
609 | bordouest=max(i-nneigh,1) |
---|
610 | bordest=min(i+nneigh,nx) |
---|
611 | bordsud=max(j-nneigh,1) |
---|
612 | bordnord=min(j+nneigh,ny) |
---|
613 | do nvx=bordouest,bordest |
---|
614 | do nvy=bordsud,bordnord |
---|
615 | imx_reduce(nvx,nvy)=imx_diag(nvx,nvy) |
---|
616 | imy_reduce(nvx,nvy)=imy_diag(nvx,nvy) |
---|
617 | enddo |
---|
618 | enddo |
---|
619 | endif |
---|
620 | enddo |
---|
621 | enddo |
---|
622 | |
---|
623 | ! -- we need to keep imx_diag for the domain edges for land terminating geometries |
---|
624 | imx_reduce(:,1) = imx_diag(:,1) |
---|
625 | imx_reduce(:,2) = imx_diag(:,2) |
---|
626 | imx_reduce(:,ny) = imx_diag(:,ny) |
---|
627 | imx_reduce(1,:) = imx_diag(1,:) |
---|
628 | imx_reduce(2,:) = imx_diag(2,:) |
---|
629 | imx_reduce(nx,:) = imx_diag(nx,:) |
---|
630 | |
---|
631 | imy_reduce(:,1) = imy_diag(:,1) |
---|
632 | imy_reduce(:,2) = imy_diag(:,2) |
---|
633 | imy_reduce(:,ny) = imy_diag(:,ny) |
---|
634 | imy_reduce(1,:) = imy_diag(1,:) |
---|
635 | imy_reduce(2,:) = imy_diag(2,:) |
---|
636 | imy_reduce(nx,:) = imy_diag(nx,:) |
---|
637 | |
---|
638 | |
---|
639 | else ! we do not reduce the extent |
---|
640 | imx_reduce(:,:)=imx_diag(:,:) |
---|
641 | imy_reduce(:,:)=imy_diag(:,:) |
---|
642 | endif |
---|
643 | |
---|
644 | end subroutine imx_imy_nx_ny_reduce |
---|
645 | |
---|
646 | |
---|
647 | !___________________________________________________________________________ |
---|
648 | ! pour imposer les conditions de mismip sur les bords du fleuve |
---|
649 | ! a appeler apres imx_imy_nx_ny |
---|
650 | |
---|
651 | subroutine mismip_boundary_cond |
---|
652 | if (itracebug.eq.1) call tracebug(' Subroutine mismip_boundray_cond') |
---|
653 | |
---|
654 | ! Condition pas de flux sur les bords nord et sud |
---|
655 | |
---|
656 | imy_diag(:,2) = 1 |
---|
657 | imy_diag(:,3) = 1 |
---|
658 | imy_diag(:,ny-1) = 1 |
---|
659 | imy_diag(:,ny) = 1 |
---|
660 | |
---|
661 | |
---|
662 | Uybar(:,2) = 0. |
---|
663 | Uybar(:,3) = 0. |
---|
664 | Uybar(:,ny-1) = 0. |
---|
665 | Uybar(:,ny) = 0. |
---|
666 | |
---|
667 | |
---|
668 | ! condition pas de cisaillement sur les bords nord et sud |
---|
669 | imx_diag(:,2) = -1 |
---|
670 | imx_diag(:,ny-1) = -3 |
---|
671 | |
---|
672 | ! coins |
---|
673 | imx_diag(2,2) = -41 |
---|
674 | imx_diag(nx,2) = -12 |
---|
675 | imx_diag(2,ny-1) = -34 |
---|
676 | imx_diag(nx,ny-1) = -23 |
---|
677 | imx_diag(1,:) = 0 ! ces points sont hors grille |
---|
678 | |
---|
679 | end subroutine mismip_boundary_cond |
---|
680 | |
---|
681 | !end module diagno_L2_mod |
---|
682 | end module diagno_mod |
---|