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

Last change on this file since 140 was 140, checked in by aquiquet, 7 years ago

Additions for initMIP outputs

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