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

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

Further 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  !$OMP BARRIER
276
277  ! afq -- initMIP outputs:
278  !$OMP DO
279  do j=2,ny-1 
280     do i=2,nx-1
281        taub(i,j) = sqrt( ((tobmx(i,j)+tobmx(i+1,j))*0.5)**2 &
282             + ((tobmy(i,j)+tobmy(i,j+1))*0.5)**2 )
283     end do
284  end do
285  !$OMP END DO
286  !$OMP WORKSHARE
287  debug_3d(:,:,117) = taub(:,:)
288  !$OMP END WORKSHARE
289
290 
291  ! Mise ne réserve des vitesses flgzmx et flgzmy
292  !$OMP WORKSHARE
293  where (flgzmx(:,:))
294     uxflgz(:,:)=uxbar(:,:)
295  elsewhere
296     uxflgz(:,:)=0.
297  endwhere
298
299  where (flgzmy(:,:))
300     uyflgz(:,:)=uybar(:,:)
301  elsewhere
302     uyflgz(:,:)=0.
303  endwhere
304  !$OMP END WORKSHARE
305  !$OMP END PARALLEL
306  !i=92
307  !j=152
308  !write(6,*) 'time',time, uxbar(92,152),gzmx(92,152),ilemx(92,152),flotmx(92,152), flgzmx(92,152)
309
310  return
311end subroutine diagnoshelf
312
313
314!-------------------------------------------------------------------
315subroutine calc_pvi
316
317! calcule les viscosites integrees pvi et pvm
318! loi polynomiale + couplage thermomécanique
319!
320!     Attention ne marche que si la loi est la loi en n=3 + n=1
321!     y compris le pur glen (n=3) ou le pur Newtonien (n=1)
322! --------------------------------------------------------------------
323
324!     les deformations sont supposées indépendantes de la profondeur
325!     et sont calculées dans strain_rate (appelé par main)       
326
327!     eps(i,j) -> eps
328!     ttau -> tau (2eme invariant du deviateur des contraintes)
329!     BT2 loi en n=3, phiphi loi en n=1 calculés dans flowlaw
330
331!
332     
333
334! La viscosité est calculée partout y compris pour la glace posée (ou elle est moins
335! précise qu'un calcul direct avec la loi de déformation.)
336! Pour les noeuds posés mais ayant un voisin stream ou flottant, on calcule
337! la viscosité avec stream/shelves
338! le calcul se fait sur les noeuds majeurs
339
340!$  integer :: rang ,nb_taches
341!$  logical :: paral
342
343  integer :: t1,t2,ir
344  real                           :: temps, t_cpu_0, t_cpu_1, t_cpu, norme 
345
346if (itracebug.eq.1)  call tracebug(' Calc pvi')
347
348!$OMP PARALLEL PRIVATE(rang,iplus1,jplus1,sf3,sf1,BT2,phiphi,ttau,d02,discr)
349!$ paral = OMP_IN_PARALLEL()
350!$ rang=OMP_GET_THREAD_NUM()
351!$ nb_taches=OMP_GET_NUM_THREADS()
352
353!$OMP WORKSHARE
354pvi(:,:)  = pvimin
355Abar(:,:) = 0.
356!$OMP END WORKSHARE
357
358!$OMP DO
359do j=1,ny
360   do i=1,nx
361      iplus1=min(i+1,nx)
362      jplus1=min(j+1,ny)
363
364
365      if (flot(i,j)) then    ! noeuds flottants
366         sf3=sf03*sf_shelf     
367         sf1=sf01*sf_shelf
368     
369      else if (gzmx(i,j).or.gzmx(iplus1,j).or.gzmy(i,j).or.gzmy(i,jplus1)) then
370         sf1=sf01
371         sf3=max(sf03,0.01)   ! pour les fleuves de glace, un peu de Glen
372     
373      else if (ilemx(i,j).or.ilemx(iplus1,j).or.ilemy(i,j).or.ilemy(i,jplus1)) then
374         sf1=sf01
375         sf3=max(sf03,0.01)   ! pour les iles aussi
376
377      else
378!         sf1=1
379!         sf3=1
380         sf1=sf01             ! pour la viscosite anisotrope (ici en longitudinal)
381         sf3=sf03
382      endif
383
384
385      do k=1,nz 
386
387         BT2=BTT(i,j,k,1)*sf3    ! changement du sf
388         phiphi=BTT(i,j,k,2)*sf1 !  changement du sf
389
390
391
392         if (BT2.lt.1.e-25) then      ! pur newtonien
393            visc(i,j,k)=1./phiphi
394            ttau=2.*visc(i,j,k)*eps(i,j)
395         else                         ! polynomial
396
397
398!  en mettant Bt2 en facteur
399            d02=eps(i,j) 
400            discr=((phiphi/3.)**3.)/Bt2+d02**2
401            discr=discr**0.5
402
403            ttau=(d02+discr)**(1./3.)-(discr-d02)**(1./3.)
404
405
406            ttau=ttau*Bt2**(-1./3.)
407
408
409            visc(i,j,k)=Bt2*ttau*ttau+phiphi
410
411            if (visc(i,j,k).gt.1.e-15) then
412               visc(i,j,k)=1./visc(i,j,k)
413            else
414               visc(i,j,k)=1.e15
415            endif
416         endif
417         pvi(i,j)=pvi(i,j)+visc(i,j,k)
418         Abar(i,j) =(Bt2/2.)**(-1./3.) + Abar(i,j) 
419
420         Taushelf(i,j)=Taushelf(i,j)+ttau
421
422
423      end do
424     
425
426   
427      pvi(i,j)  = pvi(i,j)*H(i,j)/nz
428      Abar(i,j) = (Abar(i,j) /nz)**(-3.)
429     
430
431      Taushelf(i,j)=Taushelf(i,j)/nz
432   
433   
434   end do
435end do
436!$OMP END DO
437
438! cas des noeuds fictifs, si l'épaisseur est très petite
439! pvimin est très petit
440!$OMP WORKSHARE
441where (H(:,:).le.1.)
442   pvi(:,:) = pvimin
443end where
444
445where (ramollo(:,:).ge..5)
446   pvi(:,:) = pvimin
447end where
448
449
450debug_3D(:,:,27)=pvi(:,:)
451!$OMP END WORKSHARE
452! attention run 35
453!--------------------
454!!$if (time.gt.10.) then
455!!$   where (flot(:,:))
456!!$      pvi(:,:)=pvimin
457!!$   end where
458!!$end if
459
460!  calcul de la viscosite integree au milieu des mailles (pvm)
461!$OMP DO
462do i=2,nx
463   do j=2,ny
464
465! les lignes suivantes pour un pvm moyenne des pvi
466      pvm(i,j)=0.25*((pvi(i,j)+pvi(i-1,j-1))+    &   
467           (pvi(i,j-1)+pvi(i-1,j)))
468
469   end do
470end do
471!$OMP END DO
472!$OMP END PARALLEL
473
474end subroutine calc_pvi
475!------------------------------------------------------------------
476
477subroutine imx_imy_nx_ny
478
479! definition des masques
480! pour rempli_L2 : calcule les masques imx et imy qui
481! donnent les cas de conditions aux limites
482! version pour travailler sur tout le domaine nx ny
483!----------------------------------------------------
484
485
486
487!   -34                -3   Nord              -23   
488!     !----------------------------------------!
489!     !                                        !
490!     !              1 (prescrite)             !
491!  -4 !                   ou                   ! -2
492! West!              2 (L2)                    !  Est
493!     !                                        !
494!     !                                        !
495!     !----------------------------------------!
496!    -41               -1  Sud                -12
497
498!$OMP PARALLEL
499!$OMP WORKSHARE
500imx_diag(:,:)=0
501imy_diag(:,:)=0
502
503! a l'interieur du domaine
504!-------------------------
505
506where (flgzmx(:,:))
507   imx_diag(:,:)=2             ! shelf ou stream       
508elsewhere
509   imx_diag(:,:)=1             ! vitesse imposee
510end where
511
512where (flgzmy(:,:))     
513   imy_diag(:,:)=2             ! shelf ou stream
514elsewhere
515   imy_diag(:,:)=1             ! vitesse imposee
516end where
517
518! bord sud
519imx_diag(:,1)=-1
520imy_diag(:,2)=-1
521
522! bord nord
523imx_diag(:,ny)=-3
524imy_diag(:,ny)=-3
525
526! bord Est
527imx_diag(1,:)=0    ! hors domaine a cause des mailles alternees
528imx_diag(2,:)=-4
529imy_diag(1,:)=-4
530
531! bord West
532imx_diag(nx,:)=-2
533imy_diag(nx,:)=-2
534
535! Coins
536imx_diag(2,1)=-41       ! SW
537imy_diag(1,2)=-41
538
539imx_diag(nx,1)=-12      ! SE
540imy_diag(nx,2)=-12
541
542imx_diag(nx,ny)=-23     ! NE
543imy_diag(nx,ny)=-23
544
545imx_diag(2,ny)=-34      ! NW
546imy_diag(1,ny)=-34
547
548! hors domaine
549imx_diag(1,:)=0     ! hors domaine a cause des mailles alternees
550imy_diag(:,1)=0     ! hors domaine a cause des mailles alternees
551!$OMP END WORKSHARE
552!$OMP END PARALLEL
553 
554end subroutine imx_imy_nx_ny
555!___________________________________________________________________________
556! pour imposer les conditions de mismip sur les bords du fleuve
557! a appeler apres imx_imy_nx_ny
558
559subroutine mismip_boundary_cond
560if (itracebug.eq.1)  call tracebug(' Subroutine mismip_boundray_cond')
561
562! Condition pas de flux sur les bords nord et sud
563
564 imy_diag(:,2)      = 1
565 imy_diag(:,3)      = 1
566 imy_diag(:,ny-1)   = 1
567 imy_diag(:,ny)     = 1
568
569
570 Uybar(:,2)    = 0.
571 Uybar(:,3)    = 0.
572 Uybar(:,ny-1) = 0.
573 Uybar(:,ny)   = 0.
574
575
576! condition pas de cisaillement sur les bords nord et sud
577imx_diag(:,2)    = -1
578imx_diag(:,ny-1) = -3
579
580! coins
581imx_diag(2,2)     =  -41
582imx_diag(nx,2)    =  -12
583imx_diag(2,ny-1)  =  -34
584imx_diag(nx,ny-1) =  -23
585imx_diag(1,:)     =    0     ! ces points sont hors grille
586
587end subroutine mismip_boundary_cond
588
589!end module diagno_L2_mod
590end module diagno_mod
Note: See TracBrowser for help on using the repository browser.