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

Last change on this file since 5182 was 5182, checked in by abarral, 10 days ago

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