source: trunk/SOURCES/New-remplimat/diagno-L2_mod.f90 @ 65

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

Deleting unused variables and move old sources

File size: 11.8 KB
Line 
1!module diagno_L2_mod               ! Nouvelle version, compatible remplimat 2008 Cat
2module diagno_mod                   ! nom pendant les tests
3
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   where (flot(:,:).and.(H.gt.1))  ! valeur analytique pour les shelfs
113      pvi(:,:) = (4./coef_Sflot/rog)**2/btt(:,:,1,1)/H(:,:)
114   end where
115
116  ! avec couplage thermomecanique
117!  write(166,*) ' apres call calc_pvi',itour_pvi
118
119else
120   call calc_pvi   
121end if
122
123  call imx_imy_nx_ny         ! pour rempli_L2 : calcule les masques imx et imy qui
124  ! donnent les cas de conditions aux limites
125  !
126  ! version pour travailler sur tout le domaine nx ny
127
128  if (geoplace(1:5).eq.'mism3') call mismip_boundary_cond
129
130
131  ! appel a la routine rempl_L2 -------------------domaine nx x ny ------------
132  !
133
134  ! pour tout le domaine
135  nxd1=1
136  nxd2=nx
137  nyd1=1
138  nyd2=ny
139
140  !call rempli_L2(1,nx,1,ny,uxbar,uybar,uxb1,uyb1,imx_diag,imy_diag,ifail_diagno)   
141  !nxd1=15
142  !nxd2=19
143  !nyd1=30
144  !nyd2=34
145
146  !nxd1=35
147  !nxd2=60
148  !nyd1=35
149  !nyd2=60
150
151
152  call rempli_L2(nxd1,nxd2,nyd1,nyd2,uxbar(nxd1:nxd2,nyd1:nyd2),uybar(nxd1:nxd2,nyd1:nyd2), & 
153       uxb1(nxd1:nxd2,nyd1:nyd2),uyb1(nxd1:nxd2,nyd1:nyd2),  &
154       imx_diag(nxd1:nxd2,nyd1:nyd2),imy_diag(nxd1:nxd2,nyd1:nyd2),ifail_diagno)   
155
156  ! Dans rempli_L2
157
158  !      uxprex(n1,n2)                                                 
159  !      uyprec(n1,n2) vitesses de l'iteration precedente     
160  !                                                                   
161  !      uxnew(n1,n2)                                                   
162  !      uynew(n1,n2) uynew resultat de cette iteration
163  !
164  !      imx(n1,n2) masque pour imposer les vitesses ou leur dérivee   
165  !      imy(n1,n2) masque pour imposer les vitesses ou leur dérivée   
166
167  ! eventuellement le domaine n1,n2 peut etre un sous-domaine de nx,ny
168  ! attention il faudra alors appeler avec des sous-tableaux
169
170  ! Dans l'appel uxbar -> uxprec   et Uxb1 -> Uxnew
171
172  !---------------------------------------------------------------------------
173
174  if (ifail_diagno.gt.0) then
175     write(6,*) ' Probleme resolution systeme L2. ifail=',ifail_diagno
176     STOP
177  endif
178
179  !      nouvelles vitesses         
180
181  uxbar(:,:)=uxb1(:,:)
182  uybar(:,:)=uyb1(:,:)
183
184
185
186  !    calcul de tobmx et tobmy (frottement basal) apres calcul des vitesses
187  !    ---------------------------------------------------------------------
188  do j=1,ny
189     do i=1,nx
190        tobmx(i,j)=-betamx(i,j)*uxbar(i,j)
191        tobmy(i,j)=-betamy(i,j)*uybar(i,j)
192     enddo
193  enddo
194
195  ! Mise ne réserve des vitesses flgzmx et flgzmy
196  where (flgzmx(:,:))
197     uxflgz(:,:)=uxbar(:,:)
198  elsewhere
199     uxflgz(:,:)=0.
200  endwhere
201
202  where (flgzmy(:,:))
203     uyflgz(:,:)=uybar(:,:)
204  elsewhere
205     uyflgz(:,:)=0.
206  endwhere
207
208  !i=92
209  !j=152
210  !write(6,*) 'time',time, uxbar(92,152),gzmx(92,152),ilemx(92,152),flotmx(92,152), flgzmx(92,152)
211
212  return
213end subroutine diagnoshelf
214
215
216!-------------------------------------------------------------------
217subroutine calc_pvi
218
219! calcule les viscosites integrees pvi et pvm
220! loi polynomiale + couplage thermomécanique
221!
222!     Attention ne marche que si la loi est la loi en n=3 + n=1
223!     y compris le pur glen (n=3) ou le pur Newtonien (n=1)
224! --------------------------------------------------------------------
225
226!     les deformations sont supposées indépendantes de la profondeur
227!     et sont calculées dans strain_rate (appelé par main)       
228
229!     eps(i,j) -> eps
230!     ttau -> tau (2eme invariant du deviateur des contraintes)
231!     BT2 loi en n=3, phiphi loi en n=1 calculés dans flowlaw
232
233!
234     
235
236! La viscosité est calculée partout y compris pour la glace posée (ou elle est moins
237! précise qu'un calcul direct avec la loi de déformation.)
238! Pour les noeuds posés mais ayant un voisin stream ou flottant, on calcule
239! la viscosité avec stream/shelves
240
241! le calcul se fait sur les noeuds majeurs
242if (itracebug.eq.1)  call tracebug(' Calc pvi')
243pvi(:,:)  = pvimin
244Abar(:,:) = 0.
245
246
247do j=1,ny
248   do i=1,nx
249      iplus1=min(i+1,nx)
250      jplus1=min(j+1,ny)
251
252
253      if (flot(i,j)) then    ! noeuds flottants
254         sf3=sf03*sf_shelf     
255         sf1=sf01*sf_shelf
256     
257      else if (gzmx(i,j).or.gzmx(iplus1,j).or.gzmy(i,j).or.gzmy(i,jplus1)) then
258         sf1=sf01
259         sf3=max(sf03,0.01)   ! pour les fleuves de glace, un peu de Glen
260     
261      else if (ilemx(i,j).or.ilemx(iplus1,j).or.ilemy(i,j).or.ilemy(i,jplus1)) then
262         sf1=sf01
263         sf3=max(sf03,0.01)   ! pour les iles aussi
264
265      else
266!         sf1=1
267!         sf3=1
268         sf1=sf01             ! pour la viscosite anisotrope (ici en longitudinal)
269         sf3=sf03
270      endif
271
272
273      do k=1,nz 
274
275         BT2=BTT(i,j,k,1)*sf3    ! changement du sf
276         phiphi=BTT(i,j,k,2)*sf1 !  changement du sf
277
278
279
280         if (BT2.lt.1.e-25) then      ! pur newtonien
281            visc(i,j,k)=1./phiphi
282            ttau=2.*visc(i,j,k)*eps(i,j)
283         else                         ! polynomial
284
285
286!  en mettant Bt2 en facteur
287            d02=eps(i,j) 
288            discr=((phiphi/3.)**3.)/Bt2+d02**2
289            discr=discr**0.5
290
291            ttau=(d02+discr)**(1./3.)-(discr-d02)**(1./3.)
292
293
294            ttau=ttau*Bt2**(-1./3.)
295
296
297            visc(i,j,k)=Bt2*ttau*ttau+phiphi
298
299            if (visc(i,j,k).gt.1.e-15) then
300               visc(i,j,k)=1./visc(i,j,k)
301            else
302               visc(i,j,k)=1.e15
303            endif
304         endif
305         pvi(i,j)=pvi(i,j)+visc(i,j,k)
306         Abar(i,j) =(Bt2/2.)**(-1./3.) + Abar(i,j) 
307
308         Taushelf(i,j)=Taushelf(i,j)+ttau
309
310
311      end do
312     
313
314   
315      pvi(i,j)  = pvi(i,j)*H(i,j)/nz
316      Abar(i,j) = (Abar(i,j) /nz)**(-3.)
317     
318
319      Taushelf(i,j)=Taushelf(i,j)/nz
320   
321   
322   end do
323end do
324
325! cas des noeuds fictifs, si l'épaisseur est très petite
326! pvimin est très petit
327
328where (H(:,:).le.1.)
329   pvi(:,:) = pvimin
330end where
331
332where (ramollo(:,:).ge..5)
333   pvi(:,:) = pvimin
334end where
335
336debug_3D(:,:,27)=pvi(:,:)
337
338! attention run 35
339!--------------------
340!!$if (time.gt.10.) then
341!!$   where (flot(:,:))
342!!$      pvi(:,:)=pvimin
343!!$   end where
344!!$end if
345
346!  calcul de la viscosite integree au milieu des mailles (pvm)
347
348do i=2,nx
349   do j=2,ny
350
351! les lignes suivantes pour un pvm moyenne des pvi
352      pvm(i,j)=0.25*((pvi(i,j)+pvi(i-1,j-1))+    &   
353           (pvi(i,j-1)+pvi(i-1,j)))
354
355   end do
356end do
357end subroutine calc_pvi
358!------------------------------------------------------------------
359
360subroutine imx_imy_nx_ny
361
362! definition des masques
363! pour rempli_L2 : calcule les masques imx et imy qui
364! donnent les cas de conditions aux limites
365! version pour travailler sur tout le domaine nx ny
366!----------------------------------------------------
367
368
369
370!   -34                -3   Nord              -23   
371!     !----------------------------------------!
372!     !                                        !
373!     !              1 (prescrite)             !
374!  -4 !                   ou                   ! -2
375! West!              2 (L2)                    !  Est
376!     !                                        !
377!     !                                        !
378!     !----------------------------------------!
379!    -41               -1  Sud                -12
380
381
382
383imx_diag(:,:)=0
384imy_diag(:,:)=0
385
386! a l'interieur du domaine
387!-------------------------
388
389where (flgzmx(:,:))
390   imx_diag(:,:)=2             ! shelf ou stream       
391elsewhere
392   imx_diag(:,:)=1             ! vitesse imposee
393end where
394
395where (flgzmy(:,:))     
396   imy_diag(:,:)=2             ! shelf ou stream
397elsewhere
398   imy_diag(:,:)=1             ! vitesse imposee
399end where
400
401! bord sud
402imx_diag(:,1)=-1
403imy_diag(:,2)=-1
404
405! bord nord
406imx_diag(:,ny)=-3
407imy_diag(:,ny)=-3
408
409! bord Est
410imx_diag(1,:)=0    ! hors domaine a cause des mailles alternees
411imx_diag(2,:)=-4
412imy_diag(1,:)=-4
413
414! bord West
415imx_diag(nx,:)=-2
416imy_diag(nx,:)=-2
417
418! Coins
419imx_diag(2,1)=-41       ! SW
420imy_diag(1,2)=-41
421
422imx_diag(nx,1)=-12      ! SE
423imy_diag(nx,2)=-12
424
425imx_diag(nx,ny)=-23     ! NE
426imy_diag(nx,ny)=-23
427
428imx_diag(2,ny)=-34      ! NW
429imy_diag(1,ny)=-34
430
431! hors domaine
432imx_diag(1,:)=0     ! hors domaine a cause des mailles alternees
433imy_diag(:,1)=0     ! hors domaine a cause des mailles alternees
434 
435end subroutine imx_imy_nx_ny
436!___________________________________________________________________________
437! pour imposer les conditions de mismip sur les bords du fleuve
438! a appeler apres imx_imy_nx_ny
439
440subroutine mismip_boundary_cond
441if (itracebug.eq.1)  call tracebug(' Subroutine mismip_boundray_cond')
442
443! Condition pas de flux sur les bords nord et sud
444
445 imy_diag(:,2)      = 1
446 imy_diag(:,3)      = 1
447 imy_diag(:,ny-1)   = 1
448 imy_diag(:,ny)     = 1
449
450
451 Uybar(:,2)    = 0.
452 Uybar(:,3)    = 0.
453 Uybar(:,ny-1) = 0.
454 Uybar(:,ny)   = 0.
455
456
457! condition pas de cisaillement sur les bords nord et sud
458imx_diag(:,2)    = -1
459imx_diag(:,ny-1) = -3
460
461! coins
462imx_diag(2,2)     =  -41
463imx_diag(nx,2)    =  -12
464imx_diag(2,ny-1)  =  -34
465imx_diag(nx,ny-1) =  -23
466imx_diag(1,:)     =    0     ! ces points sont hors grille
467
468end subroutine mismip_boundary_cond
469
470!end module diagno_L2_mod
471end module diagno_mod
Note: See TracBrowser for help on using the repository browser.