source: LMDZ6/trunk/libf/dyn3d/iniacademic.F90 @ 4143

Last change on this file since 4143 was 4143, checked in by dcugnet, 2 years ago
  • Some variables are renamed or replaced by direct equivalents:
    • iso_indnum -> tracers(:)%iso_iName
    • niso_possibles -> niso
    • iqiso -> iqIsoPha ; index_trac -> itZonIso
    • ok_iso_verif -> isoCheck
    • ntraceurs_zone -> nzone ; ntraciso -> ntiso
    • qperemin -> min_qparent ; masseqmin -> min_qmass ; ratiomin -> min_ratio
  • Some renamed variables are only aliased with the older name (using USE <module>, ONLY: <oldName> => <newName>) in routines where they are repeated many times.
  • Few hard-coded indexes are now computed (examples: ilic, iso, ivap, irneb, iq_vap, iq_liq, iso_H2O, iso_HDO, iso_HTO, iso_O17, iso_O18).
  • The IF(isoCheck) test is now embedded in the check_isotopes_seq and check_isotopes_loc routines (lighter calling).
  • Property copyright set to
    Name of program: LMDZ
    Creation date: 1984
    Version: LMDZ5
    License: CeCILL version 2
    Holder: Laboratoire de m\'et\'eorologie dynamique, CNRS, UMR 8539
    See the license file in the root directory
  • Property svn:eol-style set to native
  • Property svn:keywords set to Author Date Id Revision
File size: 9.8 KB
Line 
1!
2! $Id: iniacademic.F90 4143 2022-05-09 10:35:40Z dcugnet $
3!
4SUBROUTINE iniacademic(vcov,ucov,teta,q,masse,ps,phis,time_0)
5
6  USE filtreg_mod, ONLY: inifilr
7  USE infotrac,    ONLY: nqtot, niso, tnat, alpha_ideal, iqIsoPha, tracers
8  USE control_mod, ONLY: day_step,planet_type
9  use exner_hyb_m, only: exner_hyb
10  use exner_milieu_m, only: exner_milieu
11#ifdef CPP_IOIPSL
12  USE IOIPSL, ONLY: getin
13#else
14  ! if not using IOIPSL, we still need to use (a local version of) getin
15  USE ioipsl_getincom, ONLY: getin
16#endif
17  USE Write_Field
18  USE comconst_mod, ONLY: cpp, kappa, g, daysec, dtvr, pi, im, jm
19  USE logic_mod, ONLY: iflag_phys, read_start
20  USE comvert_mod, ONLY: ap, bp, preff, presnivs, pressure_exner
21  USE temps_mod, ONLY: annee_ref, day_ini, day_ref
22  USE ener_mod, ONLY: etot0,ptot0,ztot0,stot0,ang0
23  USE readTracFiles_mod, ONLY: addPhase
24
25  !   Author:    Frederic Hourdin      original: 15/01/93
26  ! The forcing defined here is from Held and Suarez, 1994, Bulletin
27  ! of the American Meteorological Society, 75, 1825.
28
29  IMPLICIT NONE
30
31  !   Declararations:
32  !   ---------------
33
34  include "dimensions.h"
35  include "paramet.h"
36  include "comgeom.h"
37  include "academic.h"
38  include "iniprint.h"
39
40  !   Arguments:
41  !   ----------
42
43  REAL,INTENT(OUT) :: time_0
44
45  !   fields
46  REAL,INTENT(OUT) :: vcov(ip1jm,llm) ! meridional covariant wind
47  REAL,INTENT(OUT) :: ucov(ip1jmp1,llm) ! zonal covariant wind
48  REAL,INTENT(OUT) :: teta(ip1jmp1,llm) ! potential temperature (K)
49  REAL,INTENT(OUT) :: q(ip1jmp1,llm,nqtot) ! advected tracers (.../kg_of_air)
50  REAL,INTENT(OUT) :: ps(ip1jmp1) ! surface pressure (Pa)
51  REAL,INTENT(OUT) :: masse(ip1jmp1,llm) ! air mass in grid cell (kg)
52  REAL,INTENT(OUT) :: phis(ip1jmp1) ! surface geopotential
53
54  !   Local:
55  !   ------
56
57  REAL p (ip1jmp1,llmp1  )               ! pression aux interfac.des couches
58  REAL pks(ip1jmp1)                      ! exner au  sol
59  REAL pk(ip1jmp1,llm)                   ! exner au milieu des couches
60  REAL phi(ip1jmp1,llm)                  ! geopotentiel
61  REAL ddsin,zsig,tetapv,w_pv  ! variables auxiliaires
62  real tetastrat ! potential temperature in the stratosphere, in K
63  real tetajl(jjp1,llm)
64  INTEGER i,j,l,lsup,ij, iq, iName, iPhase, iqParent
65
66  REAL teta0,ttp,delt_y,delt_z,eps ! Constantes pour profil de T
67  REAL k_f,k_c_a,k_c_s         ! Constantes de rappel
68  LOGICAL ok_geost             ! Initialisation vent geost. ou nul
69  LOGICAL ok_pv                ! Polar Vortex
70  REAL phi_pv,dphi_pv,gam_pv,tetanoise   ! Constantes pour polar vortex
71
72  real zz,ran1
73  integer idum
74
75  REAL zdtvr
76 
77  character(len=*),parameter :: modname="iniacademic"
78  character(len=80) :: abort_message
79
80  ! Sanity check: verify that options selected by user are not incompatible
81  if ((iflag_phys==1).and. .not. read_start) then
82    write(lunout,*) trim(modname)," error: if read_start is set to ", &
83    " false then iflag_phys should not be 1"
84    write(lunout,*) "You most likely want an aquaplanet initialisation", &
85    " (iflag_phys >= 100)"
86    call abort_gcm(modname,"incompatible iflag_phys==1 and read_start==.false.",1)
87  endif
88 
89  !-----------------------------------------------------------------------
90  ! 1. Initializations for Earth-like case
91  ! --------------------------------------
92  !
93  ! initialize planet radius, rotation rate,...
94  call conf_planete
95
96  time_0=0.
97  day_ref=1
98  annee_ref=0
99
100  im         = iim
101  jm         = jjm
102  day_ini    = 1
103  dtvr    = daysec/REAL(day_step)
104  zdtvr=dtvr
105  etot0      = 0.
106  ptot0      = 0.
107  ztot0      = 0.
108  stot0      = 0.
109  ang0       = 0.
110
111  if (llm == 1) then
112     ! specific initializations for the shallow water case
113     kappa=1
114  endif
115
116  CALL iniconst
117  CALL inigeom
118  CALL inifilr
119
120  ! Initialize pressure and mass field if read_start=.false.
121  IF (.NOT. read_start) THEN
122     ! surface pressure
123     if (iflag_phys>2) then
124        ! specific value for CMIP5 aqua/terra planets
125        ! "Specify the initial dry mass to be equivalent to
126        !  a global mean surface pressure (101325 minus 245) Pa."
127        ps(:)=101080. 
128     else
129        ! use reference surface pressure
130        ps(:)=preff
131     endif
132
133     ! ground geopotential
134     phis(:)=0.
135     CALL pression ( ip1jmp1, ap, bp, ps, p       )
136
137     if (pressure_exner) then
138       CALL exner_hyb( ip1jmp1, ps, p, pks, pk)
139     else
140       call exner_milieu(ip1jmp1,ps,p,pks,pk)
141     endif
142     CALL massdair(p,masse)
143  ENDIF
144
145  if (llm == 1) then
146     ! initialize fields for the shallow water case, if required
147     if (.not.read_start) then
148        phis(:)=0.
149        q(:,:,:)=0
150        CALL sw_case_williamson91_6(vcov,ucov,teta,masse,ps)
151     endif
152  endif
153
154  academic_case: if (iflag_phys >= 2) then
155     ! initializations
156
157     ! 1. local parameters
158     ! by convention, winter is in the southern hemisphere
159     ! Geostrophic wind or no wind?
160     ok_geost=.TRUE.
161     CALL getin('ok_geost',ok_geost)
162     ! Constants for Newtonian relaxation and friction
163     k_f=1.                !friction
164     CALL getin('k_j',k_f)
165     k_f=1./(daysec*k_f)
166     k_c_s=4.  !cooling surface
167     CALL getin('k_c_s',k_c_s)
168     k_c_s=1./(daysec*k_c_s)
169     k_c_a=40. !cooling free atm
170     CALL getin('k_c_a',k_c_a)
171     k_c_a=1./(daysec*k_c_a)
172     ! Constants for Teta equilibrium profile
173     teta0=315.     ! mean Teta (S.H. 315K)
174     CALL getin('teta0',teta0)
175     ttp=200.       ! Tropopause temperature (S.H. 200K)
176     CALL getin('ttp',ttp)
177     eps=0.         ! Deviation to N-S symmetry(~0-20K)
178     CALL getin('eps',eps)
179     delt_y=60.     ! Merid Temp. Gradient (S.H. 60K)
180     CALL getin('delt_y',delt_y)
181     delt_z=10.     ! Vertical Gradient (S.H. 10K)
182     CALL getin('delt_z',delt_z)
183     ! Polar vortex
184     ok_pv=.false.
185     CALL getin('ok_pv',ok_pv)
186     phi_pv=-50.            ! Latitude of edge of vortex
187     CALL getin('phi_pv',phi_pv)
188     phi_pv=phi_pv*pi/180.
189     dphi_pv=5.             ! Width of the edge
190     CALL getin('dphi_pv',dphi_pv)
191     dphi_pv=dphi_pv*pi/180.
192     gam_pv=4.              ! -dT/dz vortex (in K/km)
193     CALL getin('gam_pv',gam_pv)
194     tetanoise=0.005
195     CALL getin('tetanoise',tetanoise)
196
197     ! 2. Initialize fields towards which to relax
198     ! Friction
199     knewt_g=k_c_a
200     DO l=1,llm
201        zsig=presnivs(l)/preff
202        knewt_t(l)=(k_c_s-k_c_a)*MAX(0.,(zsig-0.7)/0.3)
203        kfrict(l)=k_f*MAX(0.,(zsig-0.7)/0.3)
204     ENDDO
205     DO j=1,jjp1
206        clat4((j-1)*iip1+1:j*iip1)=cos(rlatu(j))**4
207     ENDDO
208
209     ! Potential temperature
210     DO l=1,llm
211        zsig=presnivs(l)/preff
212        tetastrat=ttp*zsig**(-kappa)
213        tetapv=tetastrat
214        IF ((ok_pv).AND.(zsig.LT.0.1)) THEN
215           tetapv=tetastrat*(zsig*10.)**(kappa*cpp*gam_pv/1000./g)
216        ENDIF
217        DO j=1,jjp1
218           ! Troposphere
219           ddsin=sin(rlatu(j))
220           tetajl(j,l)=teta0-delt_y*ddsin*ddsin+eps*ddsin &
221                -delt_z*(1.-ddsin*ddsin)*log(zsig)
222           if (planet_type=="giant") then
223             tetajl(j,l)=teta0+(delt_y*                   &
224                ((sin(rlatu(j)*3.14159*eps+0.0001))**2)   &
225                / ((rlatu(j)*3.14159*eps+0.0001)**2))     &
226                -delt_z*log(zsig)
227           endif
228           ! Profil stratospherique isotherme (+vortex)
229           w_pv=(1.-tanh((rlatu(j)-phi_pv)/dphi_pv))/2.
230           tetastrat=tetastrat*(1.-w_pv)+tetapv*w_pv             
231           tetajl(j,l)=MAX(tetajl(j,l),tetastrat) 
232        ENDDO
233     ENDDO
234
235     !          CALL writefield('theta_eq',tetajl)
236
237     do l=1,llm
238        do j=1,jjp1
239           do i=1,iip1
240              ij=(j-1)*iip1+i
241              tetarappel(ij,l)=tetajl(j,l)
242           enddo
243        enddo
244     enddo
245
246     ! 3. Initialize fields (if necessary)
247     IF (.NOT. read_start) THEN
248        ! bulk initialization of temperature
249        IF (iflag_phys>10000) THEN
250        ! Particular case to impose a constant temperature T0=0.01*iflag_physx
251           teta(:,:)= 0.01*iflag_phys/(pk(:,:)/cpp)
252        ELSE
253           teta(:,:)=tetarappel(:,:)
254        ENDIF
255        ! geopotential
256        CALL geopot(ip1jmp1,teta,pk,pks,phis,phi)
257
258        DO l=1,llm
259          print*,'presnivs,play,l',presnivs(l),(pk(1,l)/cpp)**(1./kappa)*preff
260         !pks(ij) = (cpp/preff) * ps(ij)
261         !pk(ij,1) = .5*pks(ij)
262         ! pk = cpp * (p/preff)^kappa
263        ENDDO
264
265        ! winds
266        if (ok_geost) then
267           call ugeostr(phi,ucov)
268        else
269           ucov(:,:)=0.
270        endif
271        vcov(:,:)=0.
272
273        ! bulk initialization of tracers
274        if (planet_type=="earth") then
275           ! Earth: first two tracers will be water
276           do iq=1,nqtot
277              q(:,:,iq)=0.
278              IF(tracers(iq)%name == addPhase('H2O', 'g')) q(:,:,iq)=1.e-10
279              IF(tracers(iq)%name == addPhase('H2O', 'l')) q(:,:,iq)=1.e-15
280
281              ! CRisi: init des isotopes
282              ! distill de Rayleigh très simplifiée
283              iName    = tracers(iq)%iso_iName
284              if (niso <= 0 .OR. iName <= 0) CYCLE
285              iPhase   = tracers(iq)%iso_iPhase
286              iqParent = tracers(iq)%iqParent
287              IF(tracers(iq)%iso_iZone == 0) THEN
288                 q(:,:,iq) = q(:,:,iqParent)*tnat(iName)*(q(:,:,iqParent)/30.e-3)**(alpha_ideal(iName)-1.)
289              ELSE
290                 q(:,:,iq) = q(:,:,iqIsoPha(iName,iPhase))
291              END IF
292           enddo
293        else
294           q(:,:,:)=0
295        endif ! of if (planet_type=="earth")
296
297        call check_isotopes_seq(q,1,ip1jmp1,'iniacademic_loc')
298
299        ! add random perturbation to temperature
300        idum  = -1
301        zz = ran1(idum)
302        idum  = 0
303        do l=1,llm
304           do ij=iip2,ip1jm
305              teta(ij,l)=teta(ij,l)*(1.+tetanoise*ran1(idum))
306           enddo
307        enddo
308
309        ! maintain periodicity in longitude
310        do l=1,llm
311           do ij=1,ip1jmp1,iip1
312              teta(ij+iim,l)=teta(ij,l)
313           enddo
314        enddo
315
316     ENDIF ! of IF (.NOT. read_start)
317  endif academic_case
318
319END SUBROUTINE iniacademic
Note: See TracBrowser for help on using the repository browser.