source: LMDZ5/trunk/libf/dyn3dmem/iniacademic_loc.F90 @ 2597

Last change on this file since 2597 was 2597, checked in by Ehouarn Millour, 8 years ago

Cleanup in the dynamics: get rid of comconst.h, make it a module comconst_mod.
EM

  • 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.6 KB
Line 
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
6  USE filtreg_mod, ONLY: inifilr
7  use exner_hyb_m, only: exner_hyb
8  use exner_milieu_m, only: exner_milieu
9  USE infotrac, ONLY: nqtot,niso_possibles,ok_isotopes,iqpere,ok_iso_verif,tnat,alpha_ideal, &
10        & iqiso,phase_num,iso_indnum,iso_num,zone_num
11  USE control_mod, ONLY: day_step,planet_type
12  USE parallel_lmdz, ONLY: ijb_u, ije_u, ijb_v, ije_v
13#ifdef CPP_IOIPSL
14  USE IOIPSL, ONLY: getin
15#else
16  ! if not using IOIPSL, we still need to use (a local version of) getin
17  USE ioipsl_getincom, ONLY: getin
18#endif
19  USE Write_Field
20  USE comconst_mod, ONLY: cpp, kappa, g, daysec, dtvr, pi, im, jm
21
22  !   Author:    Frederic Hourdin      original: 15/01/93
23  ! The forcing defined here is from Held and Suarez, 1994, Bulletin
24  ! of the American Meteorological Society, 75, 1825.
25
26  IMPLICIT NONE
27
28  !   Declararations:
29  !   ---------------
30
31  include "dimensions.h"
32  include "paramet.h"
33  include "comvert.h"
34  include "comgeom.h"
35  include "academic.h"
36  include "ener.h"
37  include "temps.h"
38  include "iniprint.h"
39  include "logic.h"
40
41  !   Arguments:
42  !   ----------
43
44  REAL,INTENT(OUT) :: time_0
45
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
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)
68  INTEGER i,j,l,lsup,ij
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
74  REAL phi_pv,dphi_pv,gam_pv   ! Constantes pour polar vortex
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
84  ! Sanity check: verify that options selected by user are not incompatible
85  if ((iflag_phys==1).and. .not. read_start) then
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 
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
102  annee_ref=0
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.
113  ang0       = 0.     
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
124  if (llm == 1) then
125     ! initialize fields for the shallow water case, if required
126     if (.not.read_start) then
127        phis(ijb_u:ije_u)=0.
128        q(ijb_u:ije_u,1:llm,1:nqtot)=0
129        CALL sw_case_williamson91_6_loc(vcov,ucov,teta,masse,ps)
130     endif
131  endif
132
133  academic_case: if (iflag_phys >= 2) then
134     ! initializations
135
136     ! 1. local parameters
137     ! by convention, winter is in the southern hemisphere
138     ! Geostrophic wind or no wind?
139     ok_geost=.TRUE.
140     CALL getin('ok_geost',ok_geost)
141     ! Constants for Newtonian relaxation and friction
142     k_f=1.                !friction
143     CALL getin('k_j',k_f)
144     k_f=1./(daysec*k_f)
145     k_c_s=4.  !cooling surface
146     CALL getin('k_c_s',k_c_s)
147     k_c_s=1./(daysec*k_c_s)
148     k_c_a=40. !cooling free atm
149     CALL getin('k_c_a',k_c_a)
150     k_c_a=1./(daysec*k_c_a)
151     ! Constants for Teta equilibrium profile
152     teta0=315.     ! mean Teta (S.H. 315K)
153     CALL getin('teta0',teta0)
154     ttp=200.       ! Tropopause temperature (S.H. 200K)
155     CALL getin('ttp',ttp)
156     eps=0.         ! Deviation to N-S symmetry(~0-20K)
157     CALL getin('eps',eps)
158     delt_y=60.     ! Merid Temp. Gradient (S.H. 60K)
159     CALL getin('delt_y',delt_y)
160     delt_z=10.     ! Vertical Gradient (S.H. 10K)
161     CALL getin('delt_z',delt_z)
162     ! Polar vortex
163     ok_pv=.false.
164     CALL getin('ok_pv',ok_pv)
165     phi_pv=-50.            ! Latitude of edge of vortex
166     CALL getin('phi_pv',phi_pv)
167     phi_pv=phi_pv*pi/180.
168     dphi_pv=5.             ! Width of the edge
169     CALL getin('dphi_pv',dphi_pv)
170     dphi_pv=dphi_pv*pi/180.
171     gam_pv=4.              ! -dT/dz vortex (in K/km)
172     CALL getin('gam_pv',gam_pv)
173
174     ! 2. Initialize fields towards which to relax
175     ! Friction
176     knewt_g=k_c_a
177     DO l=1,llm
178        zsig=presnivs(l)/preff
179        knewt_t(l)=(k_c_s-k_c_a)*MAX(0.,(zsig-0.7)/0.3)
180        kfrict(l)=k_f*MAX(0.,(zsig-0.7)/0.3)
181     ENDDO
182     DO j=1,jjp1
183        clat4((j-1)*iip1+1:j*iip1)=cos(rlatu(j))**4
184     ENDDO
185
186     ! Potential temperature
187     DO l=1,llm
188        zsig=presnivs(l)/preff
189        tetastrat=ttp*zsig**(-kappa)
190        tetapv=tetastrat
191        IF ((ok_pv).AND.(zsig.LT.0.1)) THEN
192           tetapv=tetastrat*(zsig*10.)**(kappa*cpp*gam_pv/1000./g)
193        ENDIF
194        DO j=1,jjp1
195           ! Troposphere
196           ddsin=sin(rlatu(j))
197           tetajl(j,l)=teta0-delt_y*ddsin*ddsin+eps*ddsin &
198                -delt_z*(1.-ddsin*ddsin)*log(zsig)
199           if (planet_type=="giant") then
200             tetajl(j,l)=teta0+(delt_y*                   &
201                ((sin(rlatu(j)*3.14159*eps+0.0001))**2)   &
202                / ((rlatu(j)*3.14159*eps+0.0001)**2))     &
203                -delt_z*log(zsig)
204           endif
205           ! Profil stratospherique isotherme (+vortex)
206           w_pv=(1.-tanh((rlatu(j)-phi_pv)/dphi_pv))/2.
207           tetastrat=tetastrat*(1.-w_pv)+tetapv*w_pv             
208           tetajl(j,l)=MAX(tetajl(j,l),tetastrat) 
209        ENDDO
210     ENDDO
211
212     !          CALL writefield('theta_eq',tetajl)
213
214     do l=1,llm
215        do j=1,jjp1
216           do i=1,iip1
217              ij=(j-1)*iip1+i
218              tetarappel(ij,l)=tetajl(j,l)
219           enddo
220        enddo
221     enddo
222
223     ! 3. Initialize fields (if necessary)
224     IF (.NOT. read_start) THEN
225       ! allocate global fields:
226!       allocate(vcov_glo(ip1jm,llm))
227       allocate(ucov_glo(ip1jmp1,llm))
228       allocate(teta_glo(ip1jmp1,llm))
229       allocate(ps_glo(ip1jmp1))
230       allocate(masse_glo(ip1jmp1,llm))
231       allocate(phis_glo(ip1jmp1))
232
233        ! surface pressure
234        if (iflag_phys>2) then
235           ! specific value for CMIP5 aqua/terra planets
236           ! "Specify the initial dry mass to be equivalent to
237           !  a global mean surface pressure (101325 minus 245) Pa."
238           ps_glo(:)=101080. 
239        else
240           ! use reference surface pressure
241           ps_glo(:)=preff
242        endif
243       
244        ! ground geopotential
245        phis_glo(:)=0.
246
247        CALL pression ( ip1jmp1, ap, bp, ps_glo, p       )
248        if (pressure_exner) then
249          CALL exner_hyb( ip1jmp1, ps_glo, p, pks, pk )
250        else
251          call exner_milieu(ip1jmp1,ps_glo,p,pks,pk)
252        endif
253        CALL massdair(p,masse_glo)
254
255        ! bulk initialization of temperature
256        teta_glo(:,:)=tetarappel(:,:)
257
258        ! geopotential
259        CALL geopot(ip1jmp1,teta_glo,pk,pks,phis_glo,phi)
260
261        ! winds
262        if (ok_geost) then
263           call ugeostr(phi,ucov_glo)
264        else
265           ucov_glo(:,:)=0.
266        endif
267        vcov(ijb_v:ije_v,1:llm)=0.
268
269        ! bulk initialization of tracers
270        if (planet_type=="earth") then
271           ! Earth: first two tracers will be water
272
273           do i=1,nqtot
274              if (i == 1) q(ijb_u:ije_u,:,i)=1.e-10
275              if (i == 2) q(ijb_u:ije_u,:,i)=1.e-15
276              if (i.gt.2) q(ijb_u:ije_u,:,i)=0.
277
278              ! CRisi: init des isotopes
279              ! distill de Rayleigh très simplifiée
280              if (ok_isotopes) then
281                if ((iso_num(i).gt.0).and.(zone_num(i).eq.0)) then         
282                   q(ijb_u:ije_u,:,i)=q(ijb_u:ije_u,:,iqpere(i))       &
283      &                  *tnat(iso_num(i))                             &
284      &                  *(q(ijb_u:ije_u,:,iqpere(i))/30.e-3)                              &
285     &                   **(alpha_ideal(iso_num(i))-1)
286                endif               
287                if ((iso_num(i).gt.0).and.(zone_num(i).eq.1)) then
288                  q(ijb_u:ije_u,:,i)=q(ijb_u:ije_u,:,iqiso(iso_indnum(i),phase_num(i)))
289                endif
290              endif !if (ok_isotopes) then
291
292           enddo
293        else
294           q(ijb_u:ije_u,:,:)=0
295        endif ! of if (planet_type=="earth")
296
297        if (ok_iso_verif) then
298           call check_isotopes(q,ijb_u,ije_u,'iniacademic_loc')
299        endif !if (ok_iso_verif) then
300
301        ! add random perturbation to temperature
302        idum  = -1
303        zz = ran1(idum)
304        idum  = 0
305        do l=1,llm
306           do ij=iip2,ip1jm
307              teta_glo(ij,l)=teta_glo(ij,l)*(1.+0.005*ran1(idum))
308           enddo
309        enddo
310
311        ! maintain periodicity in longitude
312        do l=1,llm
313           do ij=1,ip1jmp1,iip1
314              teta_glo(ij+iim,l)=teta_glo(ij,l)
315           enddo
316        enddo
317
318        ! copy data from global array to local array:
319        teta(ijb_u:ije_u,:)=teta_glo(ijb_u:ije_u,:)
320        ucov(ijb_u:ije_u,:)=ucov_glo(ijb_u:ije_u,:)
321!        vcov(ijb_v:ije_v,:)=vcov_glo(ijb_v:ije_v,:)
322        masse(ijb_u:ije_u,:)=masse_glo(ijb_u:ije_u,:)
323        ps(ijb_u:ije_u)=ps_glo(ijb_u:ije_u)
324        phis(ijb_u:ije_u)=phis_glo(ijb_u:ije_u)
325
326        deallocate(teta_glo)
327        deallocate(ucov_glo)
328!        deallocate(vcov_glo)
329        deallocate(masse_glo)
330        deallocate(ps_glo)
331        deallocate(phis_glo)
332     ENDIF ! of IF (.NOT. read_start)
333  endif academic_case
334
335END SUBROUTINE iniacademic_loc
Note: See TracBrowser for help on using the repository browser.