source: trunk/LMDZ.COMMON/libf/dyn3d/inidissip.F90 @ 1000

Last change on this file since 1000 was 979, checked in by emillour, 12 years ago

Common dynamics:
Updates in the dynamics (seq and ) to keep up with updates
in LMDZ5 (up to LMDZ5 trunk, rev 1760):

  • General stuff (essentially to keep up with Earth model):
  • Updated makelmdz_fcm and makelmdz (more control on dimension.h, added option -mem, although it is not usefull for now)
  • Updated build_gcm with more control over fcm
  • Updated create_make_gcm (enable looking for code in subdirectories)
  • bibio:
  • updates (just renaming the files actually...) new_unit.F90 => new_unit_m.F90, pchsp_95.F90 => pchsp_95_m.F90 and pchfe_95.F90 => pchfe_95_m.F90
  • filtrez:
  • mod_fft.F90: use more baseline CPP directives for preprocessor compatibility
  • mod_filtre_fft_loc.F90: added this new file
  • filtreg_mod.f90: added calls to init_..._loc
  • filtreg.F: fixed calls to DGEMM into SGEMM (preprocessing does the switch)
  • dyn3d:
  • removed obsolete files: etat0_netcdf.F90 limit_netcdf.F90

pres2lev.F90

  • added new file : pres2lev_mod.F90 (module containing "old" pres2lev)
  • gcm.F: changed args to call to inidissip (added arg "vert_prof_dissip")
  • inidissip.F90: added arg "vert_prof_dissip" and the "earth model" discterizations (flagged with "planet_type=="earth")
  • comdissnew.h: added 'vert_prof_dissip' to the common block
  • guide_mod.F90: added the "use pres2lev_mod"
  • conf_gcm.F: cosmetics, and evaluation of vert_prof_dissip, (and also of dissip_* factors, for Earth model)
  • comconst.h : added dissip_factz,dissip_zref variables (for Earth mode dissip)
  • dyn3dpar:
  • removed obsolete files: etat0_netcdf.F90 limit_netcdf.F90

pres2lev.F90 mod_const_para.F90

  • added new files: pres2lev_mod.F90 (module containing "old" pres2lev)

mod_const_mpi

  • abort_gcm : better control of abort in parallel mode
  • gcm.F: changed args to call to inidissip (added arg "vert_prof_dissip")
  • inidissip.F90: added arg "vert_prof_dissip" and the "earth model" discterizations (flagged with "planet_type=="earth")
  • comdissnew.h: added 'vert_prof_dissip' to the common block
  • filtreg_p.F : bug correction (array bounds)
  • guide_p_mod.F90 : added the "use pres2lev_mod"
  • conf_gcm.F : cosmetics (and evaluation of vert_prof_dissip , and

also of dissip_* factors, for Earth model)
plus check if "adjust" is indeed not used in OpenMP

  • comconst.h : add dissip_factz,dissip_zref variables (for Earth mode dissip)

EM

File size: 7.4 KB
RevLine 
[1]1!
[270]2! $Id: inidissip.F90 1502 2011-03-21 16:07:54Z jghattas $
[1]3!
[270]4SUBROUTINE inidissip ( lstardis,nitergdiv,nitergrot,niterh  , &
[979]5     tetagdiv,tetagrot,tetatemp, vert_prof_dissip)
[270]6  !=======================================================================
7  !   initialisation de la dissipation horizontale
8  !=======================================================================
9  !-----------------------------------------------------------------------
10  !   declarations:
11  !   -------------
[1]12
[979]13  USE control_mod, only : dissip_period,iperiod,planet_type
[1]14
[270]15  IMPLICIT NONE
16  include "dimensions.h"
17  include "paramet.h"
18  include "comdissipn.h"
19  include "comconst.h"
20  include "comvert.h"
21  include "logic.h"
22  include "iniprint.h"
[1]23
[270]24  LOGICAL,INTENT(in) :: lstardis
25  INTEGER,INTENT(in) :: nitergdiv,nitergrot,niterh
26  REAL,INTENT(in) :: tetagdiv,tetagrot,tetatemp
[1]27
[979]28  integer, INTENT(in):: vert_prof_dissip ! for the Earth model !!
29  ! Vertical profile of horizontal dissipation
30  ! Allowed values:
31  ! 0: rational fraction, function of pressure
32  ! 1: tanh of altitude
33
[270]34! Local variables:
35  REAL fact,zvert(llm),zz
[776]36  REAL zh(ip1jmp1),zu(ip1jmp1), gx(ip1jmp1), divgra(ip1jmp1)
37  real zv(ip1jm), gy(ip1jm), deltap(ip1jmp1,llm)
[270]38  REAL ullm,vllm,umin,vmin,zhmin,zhmax
[776]39  REAL zllm
[1]40
[270]41  INTEGER l,ij,idum,ii
42  REAL tetamin
[979]43  REAL pseudoz
[270]44  REAL Pup
45  character (len=80) :: abort_message
[1]46
[270]47  REAL ran1
[1]48
49
[270]50  !-----------------------------------------------------------------------
51  !
52  !   calcul des valeurs propres des operateurs par methode iterrative:
53  !   -----------------------------------------------------------------
[1]54
[270]55  crot     = -1.
56  cdivu    = -1.
57  cdivh    = -1.
58
59  !   calcul de la valeur propre de divgrad:
60  !   --------------------------------------
61  idum = 0
62  DO l = 1, llm
63     DO ij = 1, ip1jmp1
[1]64        deltap(ij,l) = 1.
[270]65     ENDDO
66  ENDDO
[1]67
[270]68  idum  = -1
69  zh(1) = RAN1(idum)-.5
70  idum  = 0
71  DO ij = 2, ip1jmp1
72     zh(ij) = RAN1(idum) -.5
73  ENDDO
[1]74
[270]75  CALL filtreg (zh,jjp1,1,2,1,.TRUE.,1)
[1]76
[270]77  CALL minmax(iip1*jjp1,zh,zhmin,zhmax )
[1]78
[270]79  IF ( zhmin .GE. zhmax  )     THEN
80     write(lunout,*)'  Inidissip  zh min max  ',zhmin,zhmax
81     abort_message='probleme generateur alleatoire dans inidissip'
82     call abort_gcm('inidissip',abort_message,1)
83  ENDIF
[1]84
[270]85  zllm = ABS( zhmax )
86  DO l = 1,50
87     IF(lstardis) THEN
[776]88        CALL divgrad2(1,zh,deltap,niterh,divgra)
[270]89     ELSE
[776]90        CALL divgrad (1,zh,niterh,divgra)
[270]91     ENDIF
[1]92
[776]93     zllm  = ABS(maxval(divgra))
94     zh = divgra / zllm
[270]95  ENDDO
[1]96
[270]97  IF(lstardis) THEN
98     cdivh = 1./ zllm
99  ELSE
100     cdivh = zllm ** ( -1./niterh )
101  ENDIF
[1]102
[270]103  !   calcul des valeurs propres de gradiv (ii =1) et  nxgrarot(ii=2)
104  !   -----------------------------------------------------------------
105  write(lunout,*)'inidissip: calcul des valeurs propres'
[1]106
[270]107  DO    ii = 1, 2
108     !
109     DO ij = 1, ip1jmp1
110        zu(ij)  = RAN1(idum) -.5
111     ENDDO
112     CALL filtreg (zu,jjp1,1,2,1,.TRUE.,1)
113     DO ij = 1, ip1jm
114        zv(ij) = RAN1(idum) -.5
115     ENDDO
116     CALL filtreg (zv,jjm,1,2,1,.FALSE.,1)
[1]117
[270]118     CALL minmax(iip1*jjp1,zu,umin,ullm )
119     CALL minmax(iip1*jjm, zv,vmin,vllm )
[1]120
[270]121     ullm = ABS ( ullm )
122     vllm = ABS ( vllm )
[1]123
[270]124     DO    l = 1, 50
125        IF(ii.EQ.1) THEN
126           !cccc             CALL covcont( 1,zu,zv,zu,zv )
127           IF(lstardis) THEN
[776]128              CALL gradiv2( 1,zu,zv,nitergdiv,gx,gy )
[270]129           ELSE
[776]130              CALL gradiv ( 1,zu,zv,nitergdiv,gx,gy )
[270]131           ENDIF
132        ELSE
133           IF(lstardis) THEN
[776]134              CALL nxgraro2( 1,zu,zv,nitergrot,gx,gy )
[270]135           ELSE
[776]136              CALL nxgrarot( 1,zu,zv,nitergrot,gx,gy )
[270]137           ENDIF
138        ENDIF
[1]139
[776]140        zllm = max(abs(maxval(gx)), abs(maxval(gy)))
141        zu = gx / zllm
142        zv = gy / zllm
[270]143     end DO
[1]144
[270]145     IF ( ii.EQ.1 ) THEN
146        IF(lstardis) THEN
147           cdivu  = 1./zllm
148        ELSE
149           cdivu  = zllm **( -1./nitergdiv )
150        ENDIF
151     ELSE
152        IF(lstardis) THEN
153           crot   = 1./ zllm
154        ELSE
155           crot   = zllm **( -1./nitergrot )
156        ENDIF
157     ENDIF
[1]158
[270]159  end DO
[1]160
[270]161  !   petit test pour les operateurs non star:
162  !   ----------------------------------------
[1]163
[270]164  !     IF(.NOT.lstardis) THEN
165  fact    = rad*24./REAL(jjm)
166  fact    = fact*fact
167  write(lunout,*)'inidissip: coef u ', fact/cdivu, 1./cdivu
168  write(lunout,*)'inidissip: coef r ', fact/crot , 1./crot
169  write(lunout,*)'inidissip: coef h ', fact/cdivh, 1./cdivh
170  !     ENDIF
[1]171
[270]172  !-----------------------------------------------------------------------
173  !   variation verticale du coefficient de dissipation:
174  !   --------------------------------------------------
[979]175 
176  if (planet_type.eq."earth") then
[1]177
[979]178   if (vert_prof_dissip == 1) then
179     do l=1,llm
180        pseudoz=8.*log(preff/presnivs(l))
181        zvert(l)=1+ &
182             (tanh((pseudoz-dissip_zref)/dissip_deltaz)+1.)/2. &
183             *(dissip_factz-1.)
184     enddo
185   else
186     DO l=1,llm
187        zvert(l)=1.
188     ENDDO
189     fact=2.
190     DO l = 1, llm
191        zz      = 1. - preff/presnivs(l)
192        zvert(l)= fact -( fact-1.)/( 1.+zz*zz )
193     ENDDO
194   endif ! of if (vert_prof_dissip == 1)
195
196  else ! other planets
197 
[270]198! First step: going from 1 to dissip_fac_mid (in gcm.def)
199!============
[979]200   DO l=1,llm
[270]201     zz      = 1. - preff/presnivs(l)
202     zvert(l)= dissip_fac_mid -( dissip_fac_mid-1.)/( 1.+zz*zz )
[979]203   ENDDO
[108]204
[979]205   write(lunout,*) 'Dissipation : '
206   write(lunout,*) 'Multiplication de la dissipation en altitude :'
207   write(lunout,*) '  dissip_fac_mid =', dissip_fac_mid
[108]208
[270]209! Second step if ok_strato:  from dissip_fac_mid to dissip_fac_up (in gcm.def)
210!==========================
211! Utilisation de la fonction tangente hyperbolique pour augmenter
212! arbitrairement la dissipation et donc la stabilite du modele a
213! partir d'une certaine altitude.
[108]214
[270]215!   Le facteur multiplicatif de basse atmosphere etant deja pris
216!   en compte, il faut diviser le facteur multiplicatif de haute
217!   atmosphere par celui-ci.
[108]218
[979]219   if (ok_strato) then
[108]220
[270]221    Pup = dissip_pupstart*exp(-0.5*dissip_deltaz/dissip_hdelta)
222    do l=1,llm
223      zvert(l)= zvert(l)*(1.0+( (dissip_fac_up/dissip_fac_mid-1) &
224                *(1-(0.5*(1+tanh(-6./dissip_deltaz*              &
225               (-dissip_hdelta*log(presnivs(l)/Pup))  ))))  ))
226    enddo
[108]227
[270]228    write(*,*) '  dissip_fac_up =', dissip_fac_up
229    write(*,*) 'Transition mid /up:  Pupstart,delta =',           &
230                   dissip_pupstart,'Pa', dissip_deltaz , '(km)'
[108]231
[979]232   endif ! of if (ok_strato)
233 
234  endif ! of if (planet_type.eq."earth")
[1]235
236
[270]237  write(lunout,*)'inidissip: Constantes de temps de la diffusion horizontale'
[1]238
[270]239  tetamin =  1.e+6
[1]240
[270]241  DO l=1,llm
242     tetaudiv(l)   = zvert(l)/tetagdiv
243     tetaurot(l)   = zvert(l)/tetagrot
244     tetah(l)      = zvert(l)/tetatemp
[1]245
[270]246     IF( tetamin.GT. (1./tetaudiv(l)) ) tetamin = 1./ tetaudiv(l)
247     IF( tetamin.GT. (1./tetaurot(l)) ) tetamin = 1./ tetaurot(l)
248     IF( tetamin.GT. (1./   tetah(l)) ) tetamin = 1./    tetah(l)
249  ENDDO
[1]250
[270]251  ! If dissip_period=0 calculate value for dissipation period, else keep value read from gcm.def
252  IF (dissip_period == 0) THEN
253     dissip_period = INT( tetamin/( 2.*dtvr*iperiod) ) * iperiod
254     write(lunout,*)'inidissip: tetamin dtvr iperiod dissip_period(intermed) ',tetamin,dtvr,iperiod,dissip_period
255     dissip_period = MAX(iperiod,dissip_period)
256  END IF
[1]257
[270]258  dtdiss  = dissip_period * dtvr
259  write(lunout,*)'inidissip: dissip_period=',dissip_period,' dtdiss=',dtdiss,' dtvr=',dtvr
[1]260
[270]261  DO l = 1,llm
262     write(lunout,*)zvert(l),dtdiss*tetaudiv(l),dtdiss*tetaurot(l), &
263          dtdiss*tetah(l)
264  ENDDO
265
266END SUBROUTINE inidissip
Note: See TracBrowser for help on using the repository browser.