source: trunk/LMDZ.COMMON/libf/dyn3d_common/inidissip.F90 @ 3556

Last change on this file since 3556 was 2181, checked in by emillour, 5 years ago

Common dyn core:
Enforce using dissip_fac_mid and dissip_fac_up parameters read from .def files to be used by inidissip.
EM

File size: 11.7 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  !=======================================================================
[1021]7  !   Initialization for horizontal (lateral) dissipation
8  !  - in all cases, there is a multiplicative coefficient which increases
9  !    the dissipation in the higher levels of the atmosphere, but there
10  !    are different ways of seting the vertical profile of this coefficient
11  !    (see code below).
12  !  - the call to dissipation, every 'dissip_period' dynamical time step,
13  !    can be imposed via 'dissip_period=...' in run.def (or via the
14  !    'idissip=...' flag, but this flag should become obsolete, and is
15  !    overridden by the 'dissip_period' flag). Note that setting
16  !    dissip_period=0 has the special meaning of requesting an "optimal"
17  !    value for "dissip_period" (then taken as the largest possible value)
18  !  - the three characteristic time scales (relative to the smallest
19  !    longitudinal grid mesh size), which are privided in run.def, which
20  !    are used for the dissipations steps are:
21  !     tetagdiv : time scale for the gradient of the divergence of velocities
22  !     tetagrot : time scale for the curl of the curl of velocities
23  !     tetatemp : time scale for the laplacian of the potential temperature
[270]24  !=======================================================================
[1]25
[979]26  USE control_mod, only : dissip_period,iperiod,planet_type
[1422]27  USE comvert_mod, ONLY: preff,presnivs,scaleheight,pseudoalt
28  USE comconst_mod, ONLY: dtvr,dtdiss,rad,pi,dissip_zref,dissip_deltaz,         &
29                dissip_factz,dissip_fac_mid,dissip_fac_up,dissip_pupstart,      &
30                dissip_hdelta   
31  USE logic_mod, ONLY: ok_strato
[1]32
[270]33  IMPLICIT NONE
34  include "dimensions.h"
35  include "paramet.h"
36  include "comdissipn.h"
37  include "iniprint.h"
[1]38
[270]39  LOGICAL,INTENT(in) :: lstardis
40  INTEGER,INTENT(in) :: nitergdiv,nitergrot,niterh
41  REAL,INTENT(in) :: tetagdiv,tetagrot,tetatemp
[1]42
[1021]43  integer, INTENT(in):: vert_prof_dissip ! Vertical profile of horizontal dissipation
44  ! For the Earth model:
[979]45  ! Allowed values:
46  ! 0: rational fraction, function of pressure
47  ! 1: tanh of altitude
[1021]48  ! For planets:
49  ! 0: use fac_mid (read from run.def)
50  ! 1: use fac_mid, fac_up, startalt, delta (hard coded in inidissip.F90)
[270]51! Local variables:
52  REAL fact,zvert(llm),zz
[776]53  REAL zh(ip1jmp1),zu(ip1jmp1), gx(ip1jmp1), divgra(ip1jmp1)
54  real zv(ip1jm), gy(ip1jm), deltap(ip1jmp1,llm)
[270]55  REAL ullm,vllm,umin,vmin,zhmin,zhmax
[776]56  REAL zllm
[1]57
[270]58  INTEGER l,ij,idum,ii
59  REAL tetamin
[979]60  REAL pseudoz
[270]61  REAL Pup
62  character (len=80) :: abort_message
[1]63
[270]64  REAL ran1
[1021]65  logical,save :: firstcall=.true.
66  real,save :: fac_mid,fac_up,startalt,delta,middle
[1]67
[1021]68  if (firstcall) then
69    firstcall=.false.
70    if ((planet_type.ne."earth").and.(vert_prof_dissip.eq.1)) then
71      ! initialize values for dissipation variation along the vertical (Mars)
72      fac_mid=3 ! coefficient for lower/middle atmosphere
73      fac_up=30 ! coefficient for upper atmosphere
74      startalt=70. ! altitude (in km) for the transition from middle to upper atm.
75      delta=30.! Size (in km) of the transition region between middle/upper atm.
[2181]76     
77      if (ok_strato) then !overwrite defaults above with values from def file
78        fac_mid=dissip_fac_mid
79        fac_up=dissip_fac_up
80        delta=dissip_deltaz
81      endif
82     
[1021]83      middle=startalt+delta/2
84      write(lunout,*)"inidissip: multiplicative factors in altitude:", &
85        " fac_mid=",fac_mid," fac_up=",fac_up
86      write(lunout,*)" transition mid/up : startalt (km) =",startalt, &
87        " delta (km) =",delta
88    endif
89  endif !of if firstcall
90
[270]91  !-----------------------------------------------------------------------
92  !
93  !   calcul des valeurs propres des operateurs par methode iterrative:
94  !   -----------------------------------------------------------------
[1]95
[270]96  crot     = -1.
97  cdivu    = -1.
98  cdivh    = -1.
99
100  !   calcul de la valeur propre de divgrad:
101  !   --------------------------------------
102  idum = 0
103  DO l = 1, llm
104     DO ij = 1, ip1jmp1
[1]105        deltap(ij,l) = 1.
[270]106     ENDDO
107  ENDDO
[1]108
[270]109  idum  = -1
110  zh(1) = RAN1(idum)-.5
111  idum  = 0
112  DO ij = 2, ip1jmp1
113     zh(ij) = RAN1(idum) -.5
114  ENDDO
[1]115
[270]116  CALL filtreg (zh,jjp1,1,2,1,.TRUE.,1)
[1]117
[270]118  CALL minmax(iip1*jjp1,zh,zhmin,zhmax )
[1]119
[270]120  IF ( zhmin .GE. zhmax  )     THEN
121     write(lunout,*)'  Inidissip  zh min max  ',zhmin,zhmax
122     abort_message='probleme generateur alleatoire dans inidissip'
123     call abort_gcm('inidissip',abort_message,1)
124  ENDIF
[1]125
[270]126  zllm = ABS( zhmax )
127  DO l = 1,50
128     IF(lstardis) THEN
[776]129        CALL divgrad2(1,zh,deltap,niterh,divgra)
[270]130     ELSE
[776]131        CALL divgrad (1,zh,niterh,divgra)
[270]132     ENDIF
[1]133
[776]134     zllm  = ABS(maxval(divgra))
135     zh = divgra / zllm
[270]136  ENDDO
[1]137
[270]138  IF(lstardis) THEN
139     cdivh = 1./ zllm
140  ELSE
141     cdivh = zllm ** ( -1./niterh )
142  ENDIF
[1]143
[270]144  !   calcul des valeurs propres de gradiv (ii =1) et  nxgrarot(ii=2)
145  !   -----------------------------------------------------------------
146  write(lunout,*)'inidissip: calcul des valeurs propres'
[1]147
[270]148  DO    ii = 1, 2
149     !
150     DO ij = 1, ip1jmp1
151        zu(ij)  = RAN1(idum) -.5
152     ENDDO
153     CALL filtreg (zu,jjp1,1,2,1,.TRUE.,1)
154     DO ij = 1, ip1jm
155        zv(ij) = RAN1(idum) -.5
156     ENDDO
157     CALL filtreg (zv,jjm,1,2,1,.FALSE.,1)
[1]158
[270]159     CALL minmax(iip1*jjp1,zu,umin,ullm )
160     CALL minmax(iip1*jjm, zv,vmin,vllm )
[1]161
[270]162     ullm = ABS ( ullm )
163     vllm = ABS ( vllm )
[1]164
[270]165     DO    l = 1, 50
166        IF(ii.EQ.1) THEN
167           !cccc             CALL covcont( 1,zu,zv,zu,zv )
168           IF(lstardis) THEN
[776]169              CALL gradiv2( 1,zu,zv,nitergdiv,gx,gy )
[270]170           ELSE
[776]171              CALL gradiv ( 1,zu,zv,nitergdiv,gx,gy )
[270]172           ENDIF
173        ELSE
174           IF(lstardis) THEN
[776]175              CALL nxgraro2( 1,zu,zv,nitergrot,gx,gy )
[270]176           ELSE
[776]177              CALL nxgrarot( 1,zu,zv,nitergrot,gx,gy )
[270]178           ENDIF
179        ENDIF
[1]180
[776]181        zllm = max(abs(maxval(gx)), abs(maxval(gy)))
182        zu = gx / zllm
183        zv = gy / zllm
[270]184     end DO
[1]185
[270]186     IF ( ii.EQ.1 ) THEN
187        IF(lstardis) THEN
188           cdivu  = 1./zllm
189        ELSE
190           cdivu  = zllm **( -1./nitergdiv )
191        ENDIF
192     ELSE
193        IF(lstardis) THEN
194           crot   = 1./ zllm
195        ELSE
196           crot   = zllm **( -1./nitergrot )
197        ENDIF
198     ENDIF
[1]199
[270]200  end DO
[1]201
[270]202  !   petit test pour les operateurs non star:
203  !   ----------------------------------------
[1]204
[270]205  !     IF(.NOT.lstardis) THEN
206  fact    = rad*24./REAL(jjm)
207  fact    = fact*fact
208  write(lunout,*)'inidissip: coef u ', fact/cdivu, 1./cdivu
209  write(lunout,*)'inidissip: coef r ', fact/crot , 1./crot
210  write(lunout,*)'inidissip: coef h ', fact/cdivh, 1./cdivh
211  !     ENDIF
[1]212
[270]213  !-----------------------------------------------------------------------
214  !   variation verticale du coefficient de dissipation:
215  !   --------------------------------------------------
[979]216 
217  if (planet_type.eq."earth") then
[1]218
[979]219   if (vert_prof_dissip == 1) then
220     do l=1,llm
221        pseudoz=8.*log(preff/presnivs(l))
222        zvert(l)=1+ &
223             (tanh((pseudoz-dissip_zref)/dissip_deltaz)+1.)/2. &
224             *(dissip_factz-1.)
225     enddo
226   else
227     DO l=1,llm
228        zvert(l)=1.
229     ENDDO
230     fact=2.
231     DO l = 1, llm
232        zz      = 1. - preff/presnivs(l)
233        zvert(l)= fact -( fact-1.)/( 1.+zz*zz )
234     ENDDO
235   endif ! of if (vert_prof_dissip == 1)
236
237  else ! other planets
238 
[1021]239   if (vert_prof_dissip == 0) then
[270]240! First step: going from 1 to dissip_fac_mid (in gcm.def)
241!============
[1021]242    DO l=1,llm
[270]243     zz      = 1. - preff/presnivs(l)
244     zvert(l)= dissip_fac_mid -( dissip_fac_mid-1.)/( 1.+zz*zz )
[1021]245    ENDDO
[108]246
[1021]247    write(lunout,*) 'Dissipation : '
248    write(lunout,*) 'Multiplication de la dissipation en altitude :'
249    write(lunout,*) '  dissip_fac_mid =', dissip_fac_mid
[108]250
[270]251! Second step if ok_strato:  from dissip_fac_mid to dissip_fac_up (in gcm.def)
252!==========================
253! Utilisation de la fonction tangente hyperbolique pour augmenter
254! arbitrairement la dissipation et donc la stabilite du modele a
255! partir d'une certaine altitude.
[108]256
[270]257!   Le facteur multiplicatif de basse atmosphere etant deja pris
258!   en compte, il faut diviser le facteur multiplicatif de haute
259!   atmosphere par celui-ci.
[108]260
[1021]261    if (ok_strato) then
[108]262
[1021]263     Pup = dissip_pupstart*exp(-0.5*dissip_deltaz/dissip_hdelta)
264     do l=1,llm
[270]265      zvert(l)= zvert(l)*(1.0+( (dissip_fac_up/dissip_fac_mid-1) &
266                *(1-(0.5*(1+tanh(-6./dissip_deltaz*              &
267               (-dissip_hdelta*log(presnivs(l)/Pup))  ))))  ))
[1021]268     enddo
[108]269
[1021]270     write(*,*) '  dissip_fac_up =', dissip_fac_up
271     write(*,*) 'Transition mid /up:  Pupstart,delta =',           &
[270]272                   dissip_pupstart,'Pa', dissip_deltaz , '(km)'
[108]273
[1021]274    endif ! of if (ok_strato)
275   elseif (vert_prof_dissip==1) then
276    DO l=1,llm
277     zz      = 1. - preff/presnivs(l)
[1238]278!     zvert(l)= dissip_fac_mid -( dissip_fac_mid-1.)/( 1.+zz*zz )
279     zvert(l)= fac_mid -( fac_mid-1.)/( 1.+zz*zz )
[1021]280     
281     zvert(l)= zvert(l)*(1.0+((fac_up/fac_mid-1)*    &
282                (1-(0.5*(1+tanh(-6./                 &
283                delta*(scaleheight*(-log(presnivs(l)/preff))-middle))))) &
284                ))
285    ENDDO
[1238]286    write(lunout,*) "inidissip: vert_prof_disip=1, scaleheight=",scaleheight
287    write(lunout,*) "           fac_mid=",fac_mid,", fac_up=",fac_up
288   
[1021]289   else
290     write(lunout,*) 'wrong value for vert_prof_dissip:',vert_prof_dissip
291     abort_message='wrong value for vert_prof_dissip'
292     call abort_gcm('inidissip',abort_message,1)
293   endif ! of (vert_prof_dissip == 0)
[979]294  endif ! of if (planet_type.eq."earth")
[1]295
296
[1021]297  write(lunout,*)'inidissip: Time constants for lateral dissipation'
[1]298
[270]299  tetamin =  1.e+6
[1]300
[270]301  DO l=1,llm
302     tetaudiv(l)   = zvert(l)/tetagdiv
303     tetaurot(l)   = zvert(l)/tetagrot
304     tetah(l)      = zvert(l)/tetatemp
[1]305
[270]306     IF( tetamin.GT. (1./tetaudiv(l)) ) tetamin = 1./ tetaudiv(l)
307     IF( tetamin.GT. (1./tetaurot(l)) ) tetamin = 1./ tetaurot(l)
308     IF( tetamin.GT. (1./   tetah(l)) ) tetamin = 1./    tetah(l)
309  ENDDO
[1]310
[270]311  ! If dissip_period=0 calculate value for dissipation period, else keep value read from gcm.def
312  IF (dissip_period == 0) THEN
313     dissip_period = INT( tetamin/( 2.*dtvr*iperiod) ) * iperiod
314     write(lunout,*)'inidissip: tetamin dtvr iperiod dissip_period(intermed) ',tetamin,dtvr,iperiod,dissip_period
315     dissip_period = MAX(iperiod,dissip_period)
316  END IF
[1]317
[270]318  dtdiss  = dissip_period * dtvr
319  write(lunout,*)'inidissip: dissip_period=',dissip_period,' dtdiss=',dtdiss,' dtvr=',dtvr
[1]320
[1021]321  write(lunout,*)'pseudoZ(km)  zvert    dt(tetagdiv)   dt(tetagrot)   dt(divgrad)'
[270]322  DO l = 1,llm
[1021]323     write(lunout,'(f6.1,x,4(1pe14.7,x))') &
324     pseudoalt(l),zvert(l),dtdiss*tetaudiv(l),dtdiss*tetaurot(l),dtdiss*tetah(l)
325     ! test if disipation is not too strong (for Explicit Euler time marching)
326     if (dtdiss*tetaudiv(l).gt.1.9) then
327       write(lunout,*)"STOP : lateral dissipation is too intense and will"
328       write(lunout,*)"       generate instabilities in the model !"
329       write(lunout,*)" You must increase tetagdiv (or increase dissip_period"
330       write(lunout,*)"                             or increase day_stap)"
331     endif
332     if (dtdiss*tetaurot(l).gt.1.9) then
333       write(lunout,*)"STOP : lateral dissipation is too intense and will"
334       write(lunout,*)"       generate instabilities in the model !"
335       write(lunout,*)" You must increase tetaurot (or increase dissip_period"
336       write(lunout,*)"                             or increase day_stap)"
337     endif
338     if (dtdiss*tetah(l).gt.1.9) then
339       write(lunout,*)"STOP : lateral dissipation is too intense and will"
340       write(lunout,*)"       generate instabilities in the model !"
341       write(lunout,*)" You must increase tetah (or increase dissip_period"
342       write(lunout,*)"                          or increase day_stap)"
343     endif
344  ENDDO ! of DO l=1,llm
[270]345
346END SUBROUTINE inidissip
Note: See TracBrowser for help on using the repository browser.