source: branches/iLoveclim/SOURCES/New-remplimat/diagno-L2_mod.f90 @ 100

Last change on this file since 100 was 77, checked in by dumas, 8 years ago

Merge branche iLOVECLIM sur rev 76

File size: 12.5 KB
Line 
1!module diagno_L2_mod               ! Nouvelle version, compatible remplimat 2008 Cat
2module diagno_mod                   ! nom pendant les tests
3 !$ USE OMP_LIB
4use module3D_phy
5use module_choix
6     
7implicit none
8
9
10
11real                   :: somint,test,delp,prec
12real, dimension(nx,ny) :: uxb1
13real, dimension(nx,ny) :: uyb1
14
15integer, dimension(nx,ny) :: imx_diag
16integer, dimension(nx,ny) :: imy_diag
17
18integer :: nxd1,nxd2     ! domaine selon x Dans l'appel rempli_L2
19integer :: nyd1,nyd2     ! domaine selon y
20
21integer :: itour_pvi
22
23integer :: ifail_diagno  ! pour recuperation d'erreur
24integer :: iplus1,jplus1
25integer ::  ctvisco,iumax,jumax
26real  :: delumax,errmax
27real  :: phiphi,bt2,d02,discr,ttau
28real :: sf3,sf1,epsxxm,epsyym,epsm,sf01,sf03    ! pour le calcul de la viscosite
29real :: viscm
30real :: sf_shelf                                ! coef mult enhancement factor pour shelves
31
32logical :: stopvisco,viscolin
33logical :: test_visc
34
35contains   
36
37!------------------------------------------------------------------------------------
38subroutine init_diagno
39
40namelist/diagno_rheol/sf01,sf03,pvimin
41
42! attribution des coefficients de viscosite
43
44! formats pour les ecritures dans 42
45428 format(A)
46
47! lecture des parametres du run                      block draghwat
48!--------------------------------------------------------------------
49
50rewind(num_param)        ! pour revenir au debut du fichier param_list.dat
51read(num_param,diagno_rheol)
52
53write(num_rep_42,428)'!___________________________________________________________' 
54write(num_rep_42,428) '&diagno_rheol              ! nom du bloc  diagno_rheol'
55write(num_rep_42,*)
56write(num_rep_42,*) 'sf01           = ',sf01 
57write(num_rep_42,*) 'sf03           = ',sf03
58write(num_rep_42,*) 'pvimin         = ',pvimin
59write(num_rep_42,*)'/'     
60write(num_rep_42,428) '! coefficients par rapport a la loi glace posee '                       
61write(num_rep_42,428) '! sf01 : coefficient viscosite loi lineaire '
62write(num_rep_42,428) '! sf03 : coefficient viscosite loi n=3 '
63write(num_rep_42,428) '! pvimin : valeur de pvi pour les noeuds fictifs ~ 1.e3'
64write(num_rep_42,428) '! tres petit par rapport aux valeurs standards ~ 1.e10'
65
66write(num_rep_42,*)
67
68
69!      Precision utilisee dans de calcul
70prec = 1.e-2
71itour_pvi=1   ! si prend les valeurs analytiques dans le shelf
72
73if (geoplace(1:5).eq.'mism3') then
74   sf_shelf = 1.
75   itour_pvi= 0
76
77else
78   sf_shelf = 0.4
79end if
80
81
82return
83end subroutine init_diagno
84
85!------------------------------------------------------------------------------------
86subroutine diagnoshelf !      Resolution numerique des equations diagnostiques
87
88
89  if (itracebug.eq.1)  call tracebug(' Entree dans diagnoshelf')
90
91
92  itour_pvi=itour_pvi+1       ! boucle sur la viscosite (pour l'instant pas actif)
93
94!  pvi(:,:)=0.
95  Taushelf(:,:)=0.
96
97  ! attention le bloc suivant est pour debug
98!!$gzmx(:,:)=.false.
99!!$gzmy(:,:)=.false.
100!!$ilemx(:,:)=.false.
101!!$ilemy(:,:)=.false.
102!!$flgzmx(:,:)=flotmx(:,:)
103!!$flgzmy(:,:)=flotmy(:,:)
104
105  call dragging                    ! doit etre appele avant imx_imy
106
107
108
109if (itour_pvi.le.1) then
110  call calc_pvi                    ! calcule les viscosites integrees
111
112!$OMP PARALLEL
113!$OMP WORKSHARE
114   where (flot(:,:).and.(H.gt.1))  ! valeur analytique pour les shelfs
115      pvi(:,:) = (4./coef_Sflot/rog)**2/btt(:,:,1,1)/H(:,:)
116   end where
117!$OMP END WORKSHARE
118!$OMP END PARALLEL
119  ! avec couplage thermomecanique
120!  write(166,*) ' apres call calc_pvi',itour_pvi
121
122else
123   call calc_pvi   
124end if
125
126  call imx_imy_nx_ny         ! pour rempli_L2 : calcule les masques imx et imy qui
127  ! donnent les cas de conditions aux limites
128  !
129  ! version pour travailler sur tout le domaine nx ny
130
131  if (geoplace(1:5).eq.'mism3') call mismip_boundary_cond
132
133
134  ! appel a la routine rempl_L2 -------------------domaine nx x ny ------------
135  !
136
137  ! pour tout le domaine
138  nxd1=1
139  nxd2=nx
140  nyd1=1
141  nyd2=ny
142
143  !call rempli_L2(1,nx,1,ny,uxbar,uybar,uxb1,uyb1,imx_diag,imy_diag,ifail_diagno)   
144  !nxd1=15
145  !nxd2=19
146  !nyd1=30
147  !nyd2=34
148
149  !nxd1=35
150  !nxd2=60
151  !nyd1=35
152  !nyd2=60
153
154
155  call rempli_L2(nxd1,nxd2,nyd1,nyd2,uxbar(nxd1:nxd2,nyd1:nyd2),uybar(nxd1:nxd2,nyd1:nyd2), & 
156       uxb1(nxd1:nxd2,nyd1:nyd2),uyb1(nxd1:nxd2,nyd1:nyd2),  &
157       imx_diag(nxd1:nxd2,nyd1:nyd2),imy_diag(nxd1:nxd2,nyd1:nyd2),ifail_diagno)   
158
159  ! Dans rempli_L2
160
161  !      uxprex(n1,n2)                                                 
162  !      uyprec(n1,n2) vitesses de l'iteration precedente     
163  !                                                                   
164  !      uxnew(n1,n2)                                                   
165  !      uynew(n1,n2) uynew resultat de cette iteration
166  !
167  !      imx(n1,n2) masque pour imposer les vitesses ou leur dérivee   
168  !      imy(n1,n2) masque pour imposer les vitesses ou leur dérivée   
169
170  ! eventuellement le domaine n1,n2 peut etre un sous-domaine de nx,ny
171  ! attention il faudra alors appeler avec des sous-tableaux
172
173  ! Dans l'appel uxbar -> uxprec   et Uxb1 -> Uxnew
174
175  !---------------------------------------------------------------------------
176
177  if (ifail_diagno.gt.0) then
178     write(6,*) ' Probleme resolution systeme L2. ifail=',ifail_diagno
179     STOP
180  endif
181
182  !      nouvelles vitesses         
183  !$OMP PARALLEL
184  !$OMP WORKSHARE
185  uxbar(:,:)=uxb1(:,:)
186  uybar(:,:)=uyb1(:,:)
187  !$OMP END WORKSHARE
188
189
190  !    calcul de tobmx et tobmy (frottement basal) apres calcul des vitesses
191  !    ---------------------------------------------------------------------
192  !$OMP DO
193  do j=1,ny
194     do i=1,nx
195        tobmx(i,j)=-betamx(i,j)*uxbar(i,j)
196        tobmy(i,j)=-betamy(i,j)*uybar(i,j)
197     enddo
198  enddo
199  !$OMP END DO
200
201  ! Mise ne réserve des vitesses flgzmx et flgzmy
202  !$OMP WORKSHARE
203  where (flgzmx(:,:))
204     uxflgz(:,:)=uxbar(:,:)
205  elsewhere
206     uxflgz(:,:)=0.
207  endwhere
208
209  where (flgzmy(:,:))
210     uyflgz(:,:)=uybar(:,:)
211  elsewhere
212     uyflgz(:,:)=0.
213  endwhere
214  !$OMP END WORKSHARE
215  !$OMP END PARALLEL
216  !i=92
217  !j=152
218  !write(6,*) 'time',time, uxbar(92,152),gzmx(92,152),ilemx(92,152),flotmx(92,152), flgzmx(92,152)
219
220  return
221end subroutine diagnoshelf
222
223
224!-------------------------------------------------------------------
225subroutine calc_pvi
226
227! calcule les viscosites integrees pvi et pvm
228! loi polynomiale + couplage thermomécanique
229!
230!     Attention ne marche que si la loi est la loi en n=3 + n=1
231!     y compris le pur glen (n=3) ou le pur Newtonien (n=1)
232! --------------------------------------------------------------------
233
234!     les deformations sont supposées indépendantes de la profondeur
235!     et sont calculées dans strain_rate (appelé par main)       
236
237!     eps(i,j) -> eps
238!     ttau -> tau (2eme invariant du deviateur des contraintes)
239!     BT2 loi en n=3, phiphi loi en n=1 calculés dans flowlaw
240
241!
242     
243
244! La viscosité est calculée partout y compris pour la glace posée (ou elle est moins
245! précise qu'un calcul direct avec la loi de déformation.)
246! Pour les noeuds posés mais ayant un voisin stream ou flottant, on calcule
247! la viscosité avec stream/shelves
248! le calcul se fait sur les noeuds majeurs
249
250!$  integer :: rang ,nb_taches
251!$  logical :: paral
252
253  integer :: t1,t2,ir
254  real                           :: temps, t_cpu_0, t_cpu_1, t_cpu, norme 
255
256if (itracebug.eq.1)  call tracebug(' Calc pvi')
257
258!$OMP PARALLEL PRIVATE(rang,iplus1,jplus1,sf3,sf1,BT2,phiphi,ttau,d02,discr)
259!$ paral = OMP_IN_PARALLEL()
260!$ rang=OMP_GET_THREAD_NUM()
261!$ nb_taches=OMP_GET_NUM_THREADS()
262
263!$OMP WORKSHARE
264pvi(:,:)  = pvimin
265Abar(:,:) = 0.
266!$OMP END WORKSHARE
267
268!$OMP DO
269do j=1,ny
270   do i=1,nx
271      iplus1=min(i+1,nx)
272      jplus1=min(j+1,ny)
273
274
275      if (flot(i,j)) then    ! noeuds flottants
276         sf3=sf03*sf_shelf     
277         sf1=sf01*sf_shelf
278     
279      else if (gzmx(i,j).or.gzmx(iplus1,j).or.gzmy(i,j).or.gzmy(i,jplus1)) then
280         sf1=sf01
281         sf3=max(sf03,0.01)   ! pour les fleuves de glace, un peu de Glen
282     
283      else if (ilemx(i,j).or.ilemx(iplus1,j).or.ilemy(i,j).or.ilemy(i,jplus1)) then
284         sf1=sf01
285         sf3=max(sf03,0.01)   ! pour les iles aussi
286
287      else
288!         sf1=1
289!         sf3=1
290         sf1=sf01             ! pour la viscosite anisotrope (ici en longitudinal)
291         sf3=sf03
292      endif
293
294
295      do k=1,nz 
296
297         BT2=BTT(i,j,k,1)*sf3    ! changement du sf
298         phiphi=BTT(i,j,k,2)*sf1 !  changement du sf
299
300
301
302         if (BT2.lt.1.e-25) then      ! pur newtonien
303            visc(i,j,k)=1./phiphi
304            ttau=2.*visc(i,j,k)*eps(i,j)
305         else                         ! polynomial
306
307
308!  en mettant Bt2 en facteur
309            d02=eps(i,j) 
310            discr=((phiphi/3.)**3.)/Bt2+d02**2
311            discr=discr**0.5
312
313            ttau=(d02+discr)**(1./3.)-(discr-d02)**(1./3.)
314
315
316            ttau=ttau*Bt2**(-1./3.)
317
318
319            visc(i,j,k)=Bt2*ttau*ttau+phiphi
320
321            if (visc(i,j,k).gt.1.e-15) then
322               visc(i,j,k)=1./visc(i,j,k)
323            else
324               visc(i,j,k)=1.e15
325            endif
326         endif
327         pvi(i,j)=pvi(i,j)+visc(i,j,k)
328         Abar(i,j) =(Bt2/2.)**(-1./3.) + Abar(i,j) 
329
330         Taushelf(i,j)=Taushelf(i,j)+ttau
331
332
333      end do
334     
335
336   
337      pvi(i,j)  = pvi(i,j)*H(i,j)/nz
338      Abar(i,j) = (Abar(i,j) /nz)**(-3.)
339     
340
341      Taushelf(i,j)=Taushelf(i,j)/nz
342   
343   
344   end do
345end do
346!$OMP END DO
347
348! cas des noeuds fictifs, si l'épaisseur est très petite
349! pvimin est très petit
350!$OMP WORKSHARE
351where (H(:,:).le.1.)
352   pvi(:,:) = pvimin
353end where
354
355where (ramollo(:,:).ge..5)
356   pvi(:,:) = pvimin
357end where
358
359
360debug_3D(:,:,27)=pvi(:,:)
361!$OMP END WORKSHARE
362! attention run 35
363!--------------------
364!!$if (time.gt.10.) then
365!!$   where (flot(:,:))
366!!$      pvi(:,:)=pvimin
367!!$   end where
368!!$end if
369
370!  calcul de la viscosite integree au milieu des mailles (pvm)
371!$OMP DO
372do i=2,nx
373   do j=2,ny
374
375! les lignes suivantes pour un pvm moyenne des pvi
376      pvm(i,j)=0.25*((pvi(i,j)+pvi(i-1,j-1))+    &   
377           (pvi(i,j-1)+pvi(i-1,j)))
378
379   end do
380end do
381!$OMP END DO
382!$OMP END PARALLEL
383
384end subroutine calc_pvi
385!------------------------------------------------------------------
386
387subroutine imx_imy_nx_ny
388
389! definition des masques
390! pour rempli_L2 : calcule les masques imx et imy qui
391! donnent les cas de conditions aux limites
392! version pour travailler sur tout le domaine nx ny
393!----------------------------------------------------
394
395
396
397!   -34                -3   Nord              -23   
398!     !----------------------------------------!
399!     !                                        !
400!     !              1 (prescrite)             !
401!  -4 !                   ou                   ! -2
402! West!              2 (L2)                    !  Est
403!     !                                        !
404!     !                                        !
405!     !----------------------------------------!
406!    -41               -1  Sud                -12
407
408!$OMP PARALLEL
409!$OMP WORKSHARE
410imx_diag(:,:)=0
411imy_diag(:,:)=0
412
413! a l'interieur du domaine
414!-------------------------
415
416where (flgzmx(:,:))
417   imx_diag(:,:)=2             ! shelf ou stream       
418elsewhere
419   imx_diag(:,:)=1             ! vitesse imposee
420end where
421
422where (flgzmy(:,:))     
423   imy_diag(:,:)=2             ! shelf ou stream
424elsewhere
425   imy_diag(:,:)=1             ! vitesse imposee
426end where
427
428! bord sud
429imx_diag(:,1)=-1
430imy_diag(:,2)=-1
431
432! bord nord
433imx_diag(:,ny)=-3
434imy_diag(:,ny)=-3
435
436! bord Est
437imx_diag(1,:)=0    ! hors domaine a cause des mailles alternees
438imx_diag(2,:)=-4
439imy_diag(1,:)=-4
440
441! bord West
442imx_diag(nx,:)=-2
443imy_diag(nx,:)=-2
444
445! Coins
446imx_diag(2,1)=-41       ! SW
447imy_diag(1,2)=-41
448
449imx_diag(nx,1)=-12      ! SE
450imy_diag(nx,2)=-12
451
452imx_diag(nx,ny)=-23     ! NE
453imy_diag(nx,ny)=-23
454
455imx_diag(2,ny)=-34      ! NW
456imy_diag(1,ny)=-34
457
458! hors domaine
459imx_diag(1,:)=0     ! hors domaine a cause des mailles alternees
460imy_diag(:,1)=0     ! hors domaine a cause des mailles alternees
461!$OMP END WORKSHARE
462!$OMP END PARALLEL
463 
464end subroutine imx_imy_nx_ny
465!___________________________________________________________________________
466! pour imposer les conditions de mismip sur les bords du fleuve
467! a appeler apres imx_imy_nx_ny
468
469subroutine mismip_boundary_cond
470if (itracebug.eq.1)  call tracebug(' Subroutine mismip_boundray_cond')
471
472! Condition pas de flux sur les bords nord et sud
473
474 imy_diag(:,2)      = 1
475 imy_diag(:,3)      = 1
476 imy_diag(:,ny-1)   = 1
477 imy_diag(:,ny)     = 1
478
479
480 Uybar(:,2)    = 0.
481 Uybar(:,3)    = 0.
482 Uybar(:,ny-1) = 0.
483 Uybar(:,ny)   = 0.
484
485
486! condition pas de cisaillement sur les bords nord et sud
487imx_diag(:,2)    = -1
488imx_diag(:,ny-1) = -3
489
490! coins
491imx_diag(2,2)     =  -41
492imx_diag(nx,2)    =  -12
493imx_diag(2,ny-1)  =  -34
494imx_diag(nx,ny-1) =  -23
495imx_diag(1,:)     =    0     ! ces points sont hors grille
496
497end subroutine mismip_boundary_cond
498
499!end module diagno_L2_mod
500end module diagno_mod
Note: See TracBrowser for help on using the repository browser.