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

Last change on this file since 5220 was 5186, checked in by abarral, 2 months ago

Encapsulate files in modules

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