source: LMDZ6/branches/Amaury_dev/libf/dyn3d/iniacademic.F90 @ 5112

Last change on this file since 5112 was 5106, checked in by abarral, 4 months ago

Turn coefils.h into lmdz_coefils.f90
Put filtreg.F90 inside lmdz_filtreg.F90
Turn mod_filtreg_p.F90 into lmdz_filtreg_p.F90
Delete obsolete parafilt.h*
(lint) remove spaces between routine name and args

  • 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: 11.8 KB
Line 
1
2! $Id: iniacademic.F90 5106 2024-07-23 20:21:18Z abarral $
3
4SUBROUTINE iniacademic(vcov,ucov,teta,q,masse,ps,phis,time_0)
5
6  USE lmdz_filtreg, ONLY: inifilr
7  USE infotrac,    ONLY: nqtot, niso, iqIsoPha, tracers, getKey, isoName
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  USE IOIPSL, ONLY: getin
12  USE Write_Field
13  USE comconst_mod, ONLY: cpp, kappa, g, daysec, dtvr, pi, im, jm
14  USE logic_mod, ONLY: iflag_phys, read_start
15  USE comvert_mod, ONLY: ap, bp, preff, pa, presnivs, pressure_exner
16  USE temps_mod, ONLY: annee_ref, day_ini, day_ref
17  USE ener_mod, ONLY: etot0,ptot0,ztot0,stot0,ang0
18  USE readTracFiles_mod, ONLY: addPhase
19  use netcdf, ONLY: nf90_nowrite,nf90_open,nf90_noerr,nf90_inq_varid,nf90_close,nf90_get_var
20
21  !   Author:    Frederic Hourdin      original: 15/01/93
22  ! The forcing defined here is from Held and Suarez, 1994, Bulletin
23  ! of the American Meteorological Society, 75, 1825.
24
25  IMPLICIT NONE
26
27  !   Declararations:
28  !   ---------------
29
30  include "dimensions.h"
31  include "paramet.h"
32  include "comgeom.h"
33  include "academic.h"
34  include "iniprint.h"
35
36  !   Arguments:
37  !   ----------
38
39  REAL,INTENT(OUT) :: time_0
40
41  !   fields
42  REAL,INTENT(OUT) :: vcov(ip1jm,llm) ! meridional covariant wind
43  REAL,INTENT(OUT) :: ucov(ip1jmp1,llm) ! zonal covariant wind
44  REAL,INTENT(OUT) :: teta(ip1jmp1,llm) ! potential temperature (K)
45  REAL,INTENT(OUT) :: q(ip1jmp1,llm,nqtot) ! advected tracers (.../kg_of_air)
46  REAL,INTENT(OUT) :: ps(ip1jmp1) ! surface pressure (Pa)
47  REAL,INTENT(OUT) :: masse(ip1jmp1,llm) ! air mass in grid cell (kg)
48  REAL,INTENT(OUT) :: phis(ip1jmp1) ! surface geopotential
49
50  !   Local:
51  !   ------
52
53  REAL p (ip1jmp1,llmp1  )               ! pression aux interfac.des couches
54  REAL pks(ip1jmp1)                      ! exner au  sol
55  REAL pk(ip1jmp1,llm)                   ! exner au milieu des couches
56  REAL phi(ip1jmp1,llm)                  ! geopotentiel
57  REAL ddsin,zsig,tetapv,w_pv  ! variables auxiliaires
58  real tetastrat ! potential temperature in the stratosphere, in K
59  real tetajl(jjp1,llm)
60  INTEGER i,j,l,lsup,ij, iq, iName, iPhase, iqParent
61
62  integer :: nid_relief,varid,ierr
63  real, dimension(iip1,jjp1) :: relief
64
65  REAL teta0,ttp,delt_y,delt_z,eps ! Constantes pour profil de T
66  REAL k_f,k_c_a,k_c_s         ! Constantes de rappel
67  LOGICAL ok_geost             ! Initialisation vent geost. ou nul
68  LOGICAL ok_pv                ! Polar Vortex
69  REAL phi_pv,dphi_pv,gam_pv,tetanoise   ! Constantes pour polar vortex
70
71  real zz,ran1
72  integer idum
73
74  REAL zdtvr, tnat, alpha_ideal
75  LOGICAL,PARAMETER :: tnat1=.TRUE.
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
121  !------------------------------------------------------------------
122  ! Initialize pressure and mass field if read_start=.FALSE.
123  !------------------------------------------------------------------
124
125  IF (.NOT. read_start) THEN
126
127     !------------------------------------------------------------------
128     ! Lecture eventuelle d'un fichier de relief interpollee sur la grille
129     ! du modele.
130     ! On suppose que le fichier relief_in.nc est stoké sur une grille
131     ! iim*jjp1
132     ! Facile a créer à partir de la commande
133     ! cdo remapcon,fichier_output_phys.nc Relief.nc relief_in.nc
134     !------------------------------------------------------------------
135
136     relief=0.
137     ierr = nf90_open ('relief_in.nc', nf90_nowrite,nid_relief)
138     if (ierr==nf90_noerr) THEN
139         ierr=nf90_inq_varid(nid_relief,'RELIEF',varid)
140         if (ierr==nf90_noerr) THEN
141              ierr=nf90_get_var(nid_relief,varid,relief(1:iim,1:jjp1))
142              relief(iip1,:)=relief(1,:)
143         else
144              CALL abort_gcm ('iniacademic','variable RELIEF pas la',1)
145         endif
146     endif
147     ierr = nf90_close (nid_relief)
148
149     !------------------------------------------------------------------
150     ! Initialisation du geopotentiel au sol et de la pression
151     !------------------------------------------------------------------
152
153     PRINT*,'relief=',minval(relief),maxval(relief),'g=',g
154     do j=1,jjp1
155        do i=1,iip1
156           phis((j-1)*iip1+i)=g*relief(i,j)
157        enddo
158     enddo
159     PRINT*,'phis=',minval(phis),maxval(phis),'g=',g
160
161     ! ground geopotential
162     !phis(:)=0.
163     ps(:)=preff
164     CALL pression ( ip1jmp1, ap, bp, ps, p       )
165
166     if (pressure_exner) then
167       CALL exner_hyb( ip1jmp1, ps, p, pks, pk)
168     else
169       CALL exner_milieu(ip1jmp1,ps,p,pks,pk)
170     endif
171     CALL massdair(p,masse)
172  ENDIF
173
174  if (llm == 1) then
175     ! initialize fields for the shallow water case, if required
176     if (.not.read_start) then
177        phis(:)=0.
178        q(:,:,:)=0
179        CALL sw_case_williamson91_6(vcov,ucov,teta,masse,ps)
180     endif
181  endif
182
183  academic_case: if (iflag_phys >= 2) then
184     ! initializations
185
186     ! 1. local parameters
187     ! by convention, winter is in the southern hemisphere
188     ! Geostrophic wind or no wind?
189     ok_geost=.TRUE.
190     CALL getin('ok_geost',ok_geost)
191     ! Constants for Newtonian relaxation and friction
192     k_f=1.                !friction
193     CALL getin('k_j',k_f)
194     k_f=1./(daysec*k_f)
195     k_c_s=4.  !cooling surface
196     CALL getin('k_c_s',k_c_s)
197     k_c_s=1./(daysec*k_c_s)
198     k_c_a=40. !cooling free atm
199     CALL getin('k_c_a',k_c_a)
200     k_c_a=1./(daysec*k_c_a)
201     ! Constants for Teta equilibrium profile
202     teta0=315.     ! mean Teta (S.H. 315K)
203     CALL getin('teta0',teta0)
204     ttp=200.       ! Tropopause temperature (S.H. 200K)
205     CALL getin('ttp',ttp)
206     eps=0.         ! Deviation to N-S symmetry(~0-20K)
207     CALL getin('eps',eps)
208     delt_y=60.     ! Merid Temp. Gradient (S.H. 60K)
209     CALL getin('delt_y',delt_y)
210     delt_z=10.     ! Vertical Gradient (S.H. 10K)
211     CALL getin('delt_z',delt_z)
212     ! Polar vortex
213     ok_pv=.FALSE.
214     CALL getin('ok_pv',ok_pv)
215     phi_pv=-50.            ! Latitude of edge of vortex
216     CALL getin('phi_pv',phi_pv)
217     phi_pv=phi_pv*pi/180.
218     dphi_pv=5.             ! Width of the edge
219     CALL getin('dphi_pv',dphi_pv)
220     dphi_pv=dphi_pv*pi/180.
221     gam_pv=4.              ! -dT/dz vortex (in K/km)
222     CALL getin('gam_pv',gam_pv)
223     tetanoise=0.005
224     CALL getin('tetanoise',tetanoise)
225
226     ! 2. Initialize fields towards which to relax
227     ! Friction
228     knewt_g=k_c_a
229     DO l=1,llm
230        zsig=presnivs(l)/preff
231        knewt_t(l)=(k_c_s-k_c_a)*MAX(0.,(zsig-0.7)/0.3)
232        kfrict(l)=k_f*MAX(0.,(zsig-0.7)/0.3)
233     ENDDO
234     DO j=1,jjp1
235        clat4((j-1)*iip1+1:j*iip1)=cos(rlatu(j))**4
236     ENDDO
237
238     ! Potential temperature
239     DO l=1,llm
240        zsig=presnivs(l)/preff
241        tetastrat=ttp*zsig**(-kappa)
242        tetapv=tetastrat
243        IF ((ok_pv).AND.(zsig<0.1)) THEN
244           tetapv=tetastrat*(zsig*10.)**(kappa*cpp*gam_pv/1000./g)
245        ENDIF
246        DO j=1,jjp1
247           ! Troposphere
248           ddsin=sin(rlatu(j))
249           tetajl(j,l)=teta0-delt_y*ddsin*ddsin+eps*ddsin &
250                -delt_z*(1.-ddsin*ddsin)*log(zsig)
251           if (planet_type=="giant") then
252             tetajl(j,l)=teta0+(delt_y*                   &
253                ((sin(rlatu(j)*3.14159*eps+0.0001))**2)   &
254                / ((rlatu(j)*3.14159*eps+0.0001)**2))     &
255                -delt_z*log(zsig)
256           endif
257           ! Profil stratospherique isotherme (+vortex)
258           w_pv=(1.-tanh((rlatu(j)-phi_pv)/dphi_pv))/2.
259           tetastrat=tetastrat*(1.-w_pv)+tetapv*w_pv             
260           tetajl(j,l)=MAX(tetajl(j,l),tetastrat) 
261        ENDDO
262     ENDDO
263
264     !          CALL writefield('theta_eq',tetajl)
265
266     do l=1,llm
267        do j=1,jjp1
268           do i=1,iip1
269              ij=(j-1)*iip1+i
270              tetarappel(ij,l)=tetajl(j,l)
271           enddo
272        enddo
273     enddo
274
275     ! 3. Initialize fields (if necessary)
276     IF (.NOT. read_start) THEN
277        ! bulk initialization of temperature
278        IF (iflag_phys>10000) THEN
279        ! Particular case to impose a constant temperature T0=0.01*iflag_physx
280           teta(:,:)= 0.01*iflag_phys/(pk(:,:)/cpp)
281        ELSE
282           teta(:,:)=tetarappel(:,:)
283        ENDIF
284        ! geopotential
285        CALL geopot(ip1jmp1,teta,pk,pks,phis,phi)
286
287        DO l=1,llm
288          PRINT*,'presnivs,play,l',presnivs(l),(pk(1,l)/cpp)**(1./kappa)*preff
289         !pks(ij) = (cpp/preff) * ps(ij)
290         !pk(ij,1) = .5*pks(ij)
291         ! pk = cpp * (p/preff)^kappa
292        ENDDO
293
294        ! winds
295        if (ok_geost) then
296           CALL ugeostr(phi,ucov)
297        else
298           ucov(:,:)=0.
299        endif
300        vcov(:,:)=0.
301
302        ! bulk initialization of tracers
303        if (planet_type=="earth") then
304           ! Earth: first two tracers will be water
305           do iq=1,nqtot
306              q(:,:,iq)=0.
307              IF(tracers(iq)%name == addPhase('H2O', 'g')) q(:,:,iq)=1.e-10
308              IF(tracers(iq)%name == addPhase('H2O', 'l')) q(:,:,iq)=1.e-15
309
310              ! CRisi: init des isotopes
311              ! distill de Rayleigh très simplifiée
312              iName    = tracers(iq)%iso_iName
313              if (niso <= 0 .OR. iName <= 0) CYCLE
314              iPhase   = tracers(iq)%iso_iPhase
315              iqParent = tracers(iq)%iqParent
316              IF(tracers(iq)%iso_iZone == 0) THEN
317                 if (tnat1) then
318                         tnat=1.0
319                         alpha_ideal=1.0
320                         write(*,*) 'Attention dans iniacademic: alpha_ideal=1'
321                 else
322                    IF(getKey('tnat', tnat, isoName(iName)) .OR. getKey('alpha', alpha_ideal, isoName(iName))) &
323                    CALL abort_gcm(TRIM(modname), 'missing isotopic parameters', 1)
324                 endif
325                 q(:,:,iq) = q(:,:,iqParent)*tnat*(q(:,:,iqParent)/30.e-3)**(alpha_ideal-1.)
326              ELSE !IF(tracers(iq)%iso_iZone == 0) THEN
327                 IF(tracers(iq)%iso_iZone == 1) THEN
328                    ! correction le 14 mai 2024 pour que tous les traceurs soient de la couleur 1.
329                    ! Sinon, on va avoir des porblèmes de conservation de masse de traceurs.
330                    q(:,:,iq) = q(:,:,iqIsoPha(iName,iPhase))
331                 else !IF(tracers(iq)%iso_iZone == 1) THEN
332                    q(:,:,iq) = 0.
333                 endif !IF(tracers(iq)%iso_iZone == 1) THEN
334              END IF !IF(tracers(iq)%iso_iZone == 0) THEN
335           enddo
336        else
337           q(:,:,:)=0
338        endif ! of if (planet_type=="earth")
339
340        CALL check_isotopes_seq(q,1,ip1jmp1,'iniacademic_loc')
341
342        ! add random perturbation to temperature
343        idum  = -1
344        zz = ran1(idum)
345        idum  = 0
346        do l=1,llm
347           do ij=iip2,ip1jm
348              teta(ij,l)=teta(ij,l)*(1.+tetanoise*ran1(idum))
349           enddo
350        enddo
351
352        ! maintain periodicity in longitude
353        do l=1,llm
354           do ij=1,ip1jmp1,iip1
355              teta(ij+iim,l)=teta(ij,l)
356           enddo
357        enddo
358
359     ENDIF ! of IF (.NOT. read_start)
360  endif academic_case
361
362END SUBROUTINE iniacademic
Note: See TracBrowser for help on using the repository browser.