source: LMDZ6/trunk/libf/dyn3dmem/iniacademic_loc.F90 @ 4175

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