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

Last change on this file since 5227 was 5223, checked in by abarral, 14 months ago

Merge r5200

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