source: trunk/LMDZ.COMMON/libf/dyn3d_common/iniacademic.F90 @ 3552

Last change on this file since 3552 was 1508, checked in by emillour, 9 years ago

Common dynamics:
Updates in the dynamics (seq and ) to keep up with updates
in LMDZ5 (up to LMDZ5 trunk, rev 2325):
IMPORTANT: Modifications for isotopes are only done in dyn3d, not in dyn3dpar

as in LMDZ5 these modifications were done in dyn3dmem.
Related LMDZ5 revisions are r2270 and r2281

  • in dynlonlat_phylonlat:
  • add module "grid_atob_m.F90" (a regridding utility so far only used by phylmd/ce0l.F90, used to be dyn3d_common/grid_atob.F)
  • in misc:
  • follow up updates on wxios.F (add missing_val module variable)
  • in dyn3d_common:
  • pression.F => pression.F90
  • misc_mod.F90: moved from misc to dyn3d_common
  • added new iso_verif_dyn.F
  • covcont.F => covcont.F90
  • infotrac.F90 : add handling of isotopes (reading of corresponding traceur.def for planets not implemented)
  • dynetat0.F => dynetat0.F90 with some code factorization
  • dynredem.F => dynredem.F90 with some code factorization
  • added dynredem_mod.F90: routines used by dynredem
  • iniacademic.F90 : added isotopes-related initialization for Earth case
  • in dyn3d:
  • added check_isotopes.F
  • modified (isotopes) advtrac.F90, caladvtrac.F
  • guide_mod.F90: ported updates
  • leapfrog.F : (isotopes) updates (NB: call integrd with nqtot tracers)
  • qminimium.F : adaptations for isotopes (copied over, except that #include comvert.h is not needed).
  • vlsplt.F: adaptations for isotopes (copied over, except than #include logic.h, comvert.h not needed, and replace "include comconst.h" with use comconst_mod, ONLY: pi)
  • vlspltqs.F : same as vlsplt.F, but also keeping added modification for CP(T)
  • in dyn3dpar:
  • leapfrog_p.F: remove unecessary #ifdef CPP_EARTH cpp flag. and call integrd_p with nqtot tracers (only important for Earth)
  • dynredem_p.F => dynredem_p.F90 and some code factorization
  • and no isotopes-relates changes in dyn3dpar (since these changes have been made in LMDZ5 dyn3dmem).

EM

File size: 9.6 KB
Line 
1!
2! $Id: iniacademic.F90 1625 2012-05-09 13:14:48Z lguez $
3!
4SUBROUTINE iniacademic(vcov,ucov,teta,q,masse,ps,phis,time_0)
5
6  USE filtreg_mod, ONLY: inifilr
7  USE infotrac, ONLY : nqtot, ok_isotopes, iso_num, zone_num, &
8                       iqpere, tnat, alpha_ideal, iso_indnum, &
9                       phase_num, iqiso, ok_iso_verif
10  USE control_mod, ONLY: day_step,planet_type
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 exner_hyb_m, only: exner_hyb
19  use exner_milieu_m, only: exner_milieu
20  USE comvert_mod, ONLY: ap,bp,preff,presnivs,pressure_exner
21  USE comconst_mod, ONLY: im,jm,daysec,dtvr,kappa,cpp,g,pi
22  USE logic_mod, ONLY: iflag_phys,read_start
23  USE temps_mod, ONLY: annee_ref,day_ref,day_ini
24  USE ener_mod, ONLY: etot0,ptot0,ztot0,stot0,ang0
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
44  REAL,INTENT(OUT) :: time_0
45
46  !   fields
47  REAL,INTENT(OUT) :: vcov(ip1jm,llm) ! meridional covariant wind
48  REAL,INTENT(OUT) :: ucov(ip1jmp1,llm) ! zonal covariant wind
49  REAL,INTENT(OUT) :: teta(ip1jmp1,llm) ! potential temperature (K)
50  REAL,INTENT(OUT) :: q(ip1jmp1,llm,nqtot) ! advected tracers (.../kg_of_air)
51  REAL,INTENT(OUT) :: ps(ip1jmp1) ! surface pressure (Pa)
52  REAL,INTENT(OUT) :: masse(ip1jmp1,llm) ! air mass in grid cell (kg)
53  REAL,INTENT(OUT) :: phis(ip1jmp1) ! surface geopotential
54
55  !   Local:
56  !   ------
57
58  REAL p (ip1jmp1,llmp1  )               ! pression aux interfac.des couches
59  REAL pks(ip1jmp1)                      ! exner au  sol
60  REAL pk(ip1jmp1,llm)                   ! exner au milieu des couches
61  REAL phi(ip1jmp1,llm)                  ! geopotentiel
62  REAL ddsin,zsig,tetapv,w_pv  ! variables auxiliaires
63  real tetastrat ! potential temperature in the stratosphere, in K
64  real tetajl(jjp1,llm)
65  INTEGER i,j,l,lsup,ij
66
67  REAL teta0,ttp,delt_y,delt_z,eps ! Constantes pour profil de T
68  REAL k_f,k_c_a,k_c_s         ! Constantes de rappel
69  LOGICAL ok_geost             ! Initialisation vent geost. ou nul
70  LOGICAL ok_pv                ! Polar Vortex
71  REAL phi_pv,dphi_pv,gam_pv   ! Constantes pour polar vortex
72
73  real zz,ran1
74  integer idum
75
76  REAL alpha(ip1jmp1,llm),beta(ip1jmp1,llm),zdtvr
77 
78  character(len=*),parameter :: modname="iniacademic"
79  character(len=80) :: abort_message
80
81
82  ! Sanity check: verify that options selected by user are not incompatible
83  if ((iflag_phys==1).and. .not. read_start) then
84    write(lunout,*) trim(modname)," error: if read_start is set to ", &
85    " false then iflag_phys should not be 1"
86    write(lunout,*) "You most likely want an aquaplanet initialisation", &
87    " (iflag_phys >= 100)"
88    call abort_gcm(modname,"incompatible iflag_phys==1 and read_start==.false.",1)
89  endif
90 
91  !-----------------------------------------------------------------------
92  ! 1. Initializations for Earth-like case
93  ! --------------------------------------
94  !
95  ! initialize planet radius, rotation rate,...
96  call conf_planete
97
98  time_0=0.
99  day_ref=1
100  annee_ref=0
101
102  im         = iim
103  jm         = jjm
104  day_ini    = 1
105  dtvr    = daysec/REAL(day_step)
106  zdtvr=dtvr
107  etot0      = 0.
108  ptot0      = 0.
109  ztot0      = 0.
110  stot0      = 0.
111  ang0       = 0.
112
113  if (llm == 1) then
114     ! specific initializations for the shallow water case
115     kappa=1
116  endif
117
118  CALL iniconst
119  CALL inigeom
120  CALL inifilr
121
122  if (llm == 1) then
123     ! initialize fields for the shallow water case, if required
124     if (.not.read_start) then
125        phis(:)=0.
126        q(:,:,:)=0
127        CALL sw_case_williamson91_6(vcov,ucov,teta,masse,ps)
128     endif
129  endif
130
131  academic_case: if (iflag_phys >= 2) then
132     ! initializations
133
134     ! 1. local parameters
135     ! by convention, winter is in the southern hemisphere
136     ! Geostrophic wind or no wind?
137     ok_geost=.TRUE.
138     CALL getin('ok_geost',ok_geost)
139     ! Constants for Newtonian relaxation and friction
140     k_f=1.                !friction
141     CALL getin('k_j',k_f)
142     k_f=1./(daysec*k_f)
143     k_c_s=4.  !cooling surface
144     CALL getin('k_c_s',k_c_s)
145     k_c_s=1./(daysec*k_c_s)
146     k_c_a=40. !cooling free atm
147     CALL getin('k_c_a',k_c_a)
148     k_c_a=1./(daysec*k_c_a)
149     ! Constants for Teta equilibrium profile
150     teta0=315.     ! mean Teta (S.H. 315K)
151     CALL getin('teta0',teta0)
152     ttp=200.       ! Tropopause temperature (S.H. 200K)
153     CALL getin('ttp',ttp)
154     eps=0.         ! Deviation to N-S symmetry(~0-20K)
155     CALL getin('eps',eps)
156     delt_y=60.     ! Merid Temp. Gradient (S.H. 60K)
157     CALL getin('delt_y',delt_y)
158     delt_z=10.     ! Vertical Gradient (S.H. 10K)
159     CALL getin('delt_z',delt_z)
160     ! Polar vortex
161     ok_pv=.false.
162     CALL getin('ok_pv',ok_pv)
163     phi_pv=-50.            ! Latitude of edge of vortex
164     CALL getin('phi_pv',phi_pv)
165     phi_pv=phi_pv*pi/180.
166     dphi_pv=5.             ! Width of the edge
167     CALL getin('dphi_pv',dphi_pv)
168     dphi_pv=dphi_pv*pi/180.
169     gam_pv=4.              ! -dT/dz vortex (in K/km)
170     CALL getin('gam_pv',gam_pv)
171
172     ! 2. Initialize fields towards which to relax
173     ! Friction
174     knewt_g=k_c_a
175     DO l=1,llm
176        zsig=presnivs(l)/preff
177        knewt_t(l)=(k_c_s-k_c_a)*MAX(0.,(zsig-0.7)/0.3)
178        kfrict(l)=k_f*MAX(0.,(zsig-0.7)/0.3)
179     ENDDO
180     DO j=1,jjp1
181        clat4((j-1)*iip1+1:j*iip1)=cos(rlatu(j))**4
182     ENDDO
183
184     ! Potential temperature
185     DO l=1,llm
186        zsig=presnivs(l)/preff
187        tetastrat=ttp*zsig**(-kappa)
188        tetapv=tetastrat
189        IF ((ok_pv).AND.(zsig.LT.0.1)) THEN
190           tetapv=tetastrat*(zsig*10.)**(kappa*cpp*gam_pv/1000./g)
191        ENDIF
192        DO j=1,jjp1
193           ! Troposphere
194           ddsin=sin(rlatu(j))
195           tetajl(j,l)=teta0-delt_y*ddsin*ddsin+eps*ddsin &
196                -delt_z*(1.-ddsin*ddsin)*log(zsig)
197           if (planet_type=="giant") then
198             tetajl(j,l)=teta0+(delt_y*                   &
199                ((sin(rlatu(j)*3.14159*eps+0.0001))**2)   &
200                / ((rlatu(j)*3.14159*eps+0.0001)**2))     &
201                -delt_z*log(zsig)
202           endif
203           ! Profil stratospherique isotherme (+vortex)
204           w_pv=(1.-tanh((rlatu(j)-phi_pv)/dphi_pv))/2.
205           tetastrat=tetastrat*(1.-w_pv)+tetapv*w_pv             
206           tetajl(j,l)=MAX(tetajl(j,l),tetastrat) 
207        ENDDO
208     ENDDO
209
210     !          CALL writefield('theta_eq',tetajl)
211
212     do l=1,llm
213        do j=1,jjp1
214           do i=1,iip1
215              ij=(j-1)*iip1+i
216              tetarappel(ij,l)=tetajl(j,l)
217           enddo
218        enddo
219     enddo
220
221     ! 3. Initialize fields (if necessary)
222     IF (.NOT. read_start) THEN
223        ! surface pressure
224        if (iflag_phys>2) then
225           ! specific value for CMIP5 aqua/terra planets
226           ! "Specify the initial dry mass to be equivalent to
227           !  a global mean surface pressure (101325 minus 245) Pa."
228           ps(:)=101080. 
229        else
230           ! use reference surface pressure
231           ps(:)=preff
232        endif
233       
234        ! ground geopotential
235        phis(:)=0.
236
237        CALL pression ( ip1jmp1, ap, bp, ps, p       )
238        if (pressure_exner) then
239          CALL exner_hyb( ip1jmp1, ps, p, pks, pk)
240        else
241          call exner_milieu(ip1jmp1,ps,p,pks,pk)
242        endif
243        CALL massdair(p,masse)
244
245        ! bulk initialization of temperature
246        teta(:,:)=tetarappel(:,:)
247
248        ! geopotential
249        CALL geopot(ip1jmp1,teta,pk,pks,phis,phi)
250
251        ! winds
252        if (ok_geost) then
253           call ugeostr(phi,ucov)
254        else
255           ucov(:,:)=0.
256        endif
257        vcov(:,:)=0.
258
259        ! bulk initialization of tracers
260        if (planet_type=="earth") then
261           ! Earth: first two tracers will be water
262           do i=1,nqtot
263              if (i == 1) q(:,:,i)=1.e-10
264              if (i == 2) q(:,:,i)=1.e-15
265              if (i.gt.2) q(:,:,i)=0.
266
267              ! CRisi: init des isotopes
268              ! distill de Rayleigh très simplifiée
269              if (ok_isotopes) then
270                if ((iso_num(i).gt.0).and.(zone_num(i).eq.0)) then         
271                   q(:,:,i)=q(:,:,iqpere(i))             &
272      &                  *tnat(iso_num(i))               &
273      &                  *(q(:,:,iqpere(i))/30.e-3)      &
274      &                  **(alpha_ideal(iso_num(i))-1)
275                endif               
276                if ((iso_num(i).gt.0).and.(zone_num(i).eq.1)) then
277                  q(:,:,i)=q(:,:,iqiso(iso_indnum(i),phase_num(i)))
278                endif
279              endif !if (ok_isotopes) then
280
281           enddo
282        else
283           q(:,:,:)=0
284        endif ! of if (planet_type=="earth")
285
286        if (ok_iso_verif) then
287! Ehouarn: this will onyly work in serial mode
288!           call check_isotopes_seq(q,1,ip1jmp1,'iniacademic_loc')
289        endif !if (ok_iso_verif) then
290
291        ! add random perturbation to temperature
292        idum  = -1
293        zz = ran1(idum)
294        idum  = 0
295        do l=1,llm
296           do ij=iip2,ip1jm
297              teta(ij,l)=teta(ij,l)*(1.+0.005*ran1(idum))
298           enddo
299        enddo
300
301        ! maintain periodicity in longitude
302        do l=1,llm
303           do ij=1,ip1jmp1,iip1
304              teta(ij+iim,l)=teta(ij,l)
305           enddo
306        enddo
307
308     ENDIF ! of IF (.NOT. read_start)
309  endif academic_case
310
311END SUBROUTINE iniacademic
Note: See TracBrowser for help on using the repository browser.