source: LMDZ6/branches/Amaury_dev/libf/phylmd/dyn1d/lmdz_old_lmdz1d.F90 @ 5142

Last change on this file since 5142 was 5142, checked in by abarral, 8 weeks ago

Put cvparam.h, fcg_gcssold.h, planete.h, tsoilnudge.h, YOECUMF.h into modules

File size: 55.4 KB
Line 
1! $Id: lmdz1d.F90 3540 2019-06-25 14:50:13Z fairhead $
2
3MODULE lmdz_old_lmdz1d
4  USE lmdz_old_1dconv, ONLY: copie, get_uvd2, get_uvd, rdgrads, spaces
5  USE lmdz_1dutils, ONLY: interp2_case_vertical, nudge_uv, nudge_rht, lstendh, nudge_uv_init, &
6          nudge_rht_init, disvert0
7   PRIVATE  ! -- We'd love to put IMPLICIT NONE;  here...
8  PUBLIC old_lmdz1d
9CONTAINS
10
11  SUBROUTINE old_lmdz1d
12
13    USE ioipsl, ONLY: ju2ymds, ymds2ju, ioconf_calendar, getin
14    USE phys_state_var_mod, ONLY: phys_state_var_init, phys_state_var_end, &
15            clwcon, detr_therm, &
16            qsol, fevap, z0m, z0h, agesno, &
17            du_gwd_rando, du_gwd_front, entr_therm, f0, fm_therm, &
18            falb_dir, falb_dif, &
19            ftsol, beta_aridity, pbl_tke, pctsrf, radsol, rain_fall, snow_fall, ratqs, &
20            rnebcon, rugoro, sig1, w01, solaire_etat0, sollw, sollwdown, &
21            solsw, solswfdiff, t_ancien, q_ancien, u_ancien, v_ancien, rneb_ancien, &
22            wake_delta_pbl_TKE, delta_tsurf, wake_fip, wake_pe, &
23            wake_deltaq, wake_deltat, wake_s, awake_s, wake_dens, &
24            awake_dens, cv_gen, wake_cstar, &
25            zgam, zmax0, zmea, zpic, zsig, &
26            zstd, zthe, zval, ale_bl, ale_bl_trig, alp_bl, ql_ancien, qs_ancien, &
27            prlw_ancien, prsw_ancien, prw_ancien, &
28            u10m, v10m, ale_wake, ale_bl_stat
29
30    USE dimphy
31    USE surface_data, ONLY: type_ocean, ok_veget
32    USE pbl_surface_mod, ONLY: ftsoil, pbl_surface_init, pbl_surface_final
33    USE fonte_neige_mod, ONLY: fonte_neige_init, fonte_neige_final
34
35    USE infotrac
36    USE control_mod
37    USE indice_sol_mod
38    USE phyaqua_mod
39    USE mod_1D_cases_read2
40    USE mod_1D_amma_read
41    USE lmdz_print_control, ONLY: lunout, prt_level
42    USE iniphysiq_mod, ONLY: iniphysiq
43    USE mod_const_mpi, ONLY: comm_lmdz
44    USE physiq_mod, ONLY: physiq
45    USE comvert_mod, ONLY: presnivs, ap, bp, dpres, nivsig, nivsigs, pa, &
46            preff, aps, bps, pseudoalt, scaleheight
47    USE temps_mod, ONLY: annee_ref, calend, day_end, day_ini, day_ref, &
48            itau_dyn, itau_phy, start_time, year_len
49    USE phys_cal_mod, ONLY: year_len_phys_cal_mod => year_len
50    USE mod_1D_cases_read, ONLY: interp_case_time ! used in included old_1D_read_forc_cases.h
51    USE lmdz_cppkeys_wrapper, ONLY: CPPKEY_OUTPUTPHYSSCM
52    USE lmdz_1dutils, ONLY: fq_sat, conf_unicol, dyn1deta0, dyn1dredem
53    USE lmdz_clesphys
54    USE lmdz_flux_arp, ONLY: fsens, flat, betaevap, ust, tg, ok_flux_surf, ok_prescr_ust, ok_prescr_beta, ok_forc_tsurf
55    USE lmdz_compbl, ONLY: iflag_pbl, iflag_pbl_split, iflag_order2_sollw, ifl_pbltree
56    USE lmdz_fcs_gcssold, ONLY: imp_fcg_gcssold, ts_fcg_gcssold, Tp_fcg_gcssold, Tp_ini_gcssold, xTurb_fcg_gcssold
57    USE lmdz_tsoilnudge, ONLY: nudge_tsoil, isoil_nudge, Tsoil_nudge, tau_soil_nudge
58
59    INCLUDE "dimensions.h"
60    INCLUDE "YOMCST.h"
61    INCLUDE "dimsoil.h"
62    INCLUDE "compar1d.h"
63    INCLUDE "date_cas.h"
64
65    !=====================================================================
66    ! DECLARATIONS
67    !=====================================================================
68    !---------------------------------------------------------------------
69    !  Arguments d' initialisations de la physique (USER DEFINE)
70    !---------------------------------------------------------------------
71
72    INTEGER, parameter :: ngrid = 1
73    REAL :: zcufi = 1.
74    REAL :: zcvfi = 1.
75
76    !-      REAL :: nat_surf
77    !-      LOGICAL :: ok_flux_surf
78    !-      REAL :: fsens
79    !-      REAL :: flat
80    !-      REAL :: tsurf
81    !-      REAL :: rugos
82    !-      REAL :: qsol(1:2)
83    !-      REAL :: qsurf
84    !-      REAL :: psurf
85    !-      REAL :: zsurf
86    !-      REAL :: albedo
87    !-
88    !-      REAL :: time     = 0.
89    !-      REAL :: time_ini
90    !-      REAL :: xlat
91    !-      REAL :: xlon
92    !-      REAL :: wtsurf
93    !-      REAL :: wqsurf
94    !-      REAL :: restart_runoff
95    !-      REAL :: xagesno
96    !-      REAL :: qsolinp
97    !-      REAL :: zpicinp
98    !-
99    REAL :: fnday
100    REAL :: day, daytime
101    REAL :: day1
102    REAL :: heure
103    INTEGER :: jour
104    INTEGER :: mois
105    INTEGER :: an
106
107    !---------------------------------------------------------------------
108    !  Declarations related to forcing and initial profiles
109    !---------------------------------------------------------------------
110
111    INTEGER :: kmax = llm
112    INTEGER llm700, nq1, nq2
113    INTEGER, PARAMETER :: nlev_max = 1000, nqmx = 1000
114    REAL timestep, frac
115    REAL height(nlev_max), tttprof(nlev_max), qtprof(nlev_max)
116    real  uprof(nlev_max), vprof(nlev_max), e12prof(nlev_max)
117    real  ugprof(nlev_max), vgprof(nlev_max), wfls(nlev_max)
118    real  dqtdxls(nlev_max), dqtdyls(nlev_max)
119    real  dqtdtls(nlev_max), thlpcar(nlev_max)
120    real  qprof(nlev_max, nqmx)
121
122    !        INTEGER :: forcing_type
123    LOGICAL :: forcing_les = .FALSE.
124    LOGICAL :: forcing_armcu = .FALSE.
125    LOGICAL :: forcing_rico = .FALSE.
126    LOGICAL :: forcing_radconv = .FALSE.
127    LOGICAL :: forcing_toga = .FALSE.
128    LOGICAL :: forcing_twpice = .FALSE.
129    LOGICAL :: forcing_amma = .FALSE.
130    LOGICAL :: forcing_dice = .FALSE.
131    LOGICAL :: forcing_gabls4 = .FALSE.
132
133    LOGICAL :: forcing_GCM2SCM = .FALSE.
134    LOGICAL :: forcing_GCSSold = .FALSE.
135    LOGICAL :: forcing_sandu = .FALSE.
136    LOGICAL :: forcing_astex = .FALSE.
137    LOGICAL :: forcing_fire = .FALSE.
138    LOGICAL :: forcing_case = .FALSE.
139    LOGICAL :: forcing_case2 = .FALSE.
140    LOGICAL :: forcing_SCM = .FALSE.
141    INTEGER :: type_ts_forcing ! 0 = SST constant; 1 = SST read from a file
142    !                                                            (cf read_tsurf1d.F)
143
144    REAL wwww
145    !vertical advection computation
146    !       real d_t_z(llm), d_q_z(llm)
147    !       real d_t_dyn_z(llm), dq_dyn_z(llm)
148    !       real zz(llm)
149    !       real zfact
150
151    !flag forcings
152    LOGICAL :: nudge_wind = .TRUE.
153    LOGICAL :: nudge_thermo = .FALSE.
154    LOGICAL :: cptadvw = .TRUE.
155    !=====================================================================
156    ! DECLARATIONS FOR EACH CASE
157    !=====================================================================
158
159    INCLUDE "old_1D_decl_cases.h"
160
161    !---------------------------------------------------------------------
162    !  Declarations related to nudging
163    !---------------------------------------------------------------------
164    INTEGER :: nudge_max
165    parameter (nudge_max = 9)
166    INTEGER :: inudge_RHT = 1
167    INTEGER :: inudge_UV = 2
168    LOGICAL :: nudge(nudge_max)
169    REAL :: t_targ(llm)
170    REAL :: rh_targ(llm)
171    REAL :: u_targ(llm)
172    REAL :: v_targ(llm)
173
174    !---------------------------------------------------------------------
175    !  Declarations related to vertical discretization:
176    !---------------------------------------------------------------------
177    REAL :: pzero = 1.e5
178    REAL :: play (llm), zlay (llm), sig_s(llm), plev(llm + 1)
179    REAL :: playd(llm), zlayd(llm), ap_amma(llm + 1), bp_amma(llm + 1)
180
181    !---------------------------------------------------------------------
182    !  Declarations related to variables
183    !---------------------------------------------------------------------
184
185    REAL :: phi(llm)
186    REAL :: teta(llm), tetal(llm), temp(llm), u(llm), v(llm), w(llm)
187    REAL rot(1, llm) ! relative vorticity, in s-1
188    REAL :: rlat_rad(1), rlon_rad(1)
189    REAL :: omega(llm + 1), omega2(llm), rho(llm + 1)
190    REAL :: ug(llm), vg(llm), fcoriolis
191    REAL :: sfdt, cfdt
192    REAL :: du_phys(llm), dv_phys(llm), dt_phys(llm)
193    REAL :: dt_dyn(llm)
194    REAL :: dt_cooling(llm), d_t_adv(llm), d_t_nudge(llm)
195    REAL :: d_u_nudge(llm), d_v_nudge(llm)
196    REAL :: du_adv(llm), dv_adv(llm)
197    REAL :: du_age(llm), dv_age(llm)
198    REAL :: alpha
199    REAL :: ttt
200
201    REAL, ALLOCATABLE, DIMENSION(:, :) :: q
202    REAL, ALLOCATABLE, DIMENSION(:, :) :: dq
203    REAL, ALLOCATABLE, DIMENSION(:, :) :: dq_dyn
204    REAL, ALLOCATABLE, DIMENSION(:, :) :: d_q_adv
205    REAL, ALLOCATABLE, DIMENSION(:, :) :: d_q_nudge
206    !      REAL, ALLOCATABLE, DIMENSION(:):: d_th_adv
207
208    !---------------------------------------------------------------------
209    !  Initialization of surface variables
210    !---------------------------------------------------------------------
211    REAL :: run_off_lic_0(1)
212    REAL :: fder(1), snsrf(1, nbsrf), qsurfsrf(1, nbsrf)
213    REAL :: tsoil(1, nsoilmx, nbsrf)
214    !     REAL :: agesno(1,nbsrf)
215
216    !---------------------------------------------------------------------
217    !  Call to phyredem
218    !---------------------------------------------------------------------
219    LOGICAL :: ok_writedem = .TRUE.
220    REAL :: sollw_in = 0.
221    REAL :: solsw_in = 0.
222
223    !---------------------------------------------------------------------
224    !  Call to physiq
225    !---------------------------------------------------------------------
226    LOGICAL :: firstcall = .TRUE.
227    LOGICAL :: lastcall = .FALSE.
228    REAL :: phis(1) = 0.0
229    REAL :: dpsrf(1)
230
231    !---------------------------------------------------------------------
232    !  Initializations of boundary conditions
233    !---------------------------------------------------------------------
234    REAL, ALLOCATABLE :: phy_nat (:)  ! 0=ocean libre,1=land,2=glacier,3=banquise
235    REAL, ALLOCATABLE :: phy_alb (:)  ! Albedo land only (old value condsurf_jyg=0.3)
236    REAL, ALLOCATABLE :: phy_sst (:)  ! SST (will not be used; cf read_tsurf1d.F)
237    REAL, ALLOCATABLE :: phy_bil (:)  ! Ne sert que pour les slab_ocean
238    REAL, ALLOCATABLE :: phy_rug (:) ! Longueur rugosite utilisee sur land only
239    REAL, ALLOCATABLE :: phy_ice (:) ! Fraction de glace
240    REAL, ALLOCATABLE :: phy_fter(:) ! Fraction de terre
241    REAL, ALLOCATABLE :: phy_foce(:) ! Fraction de ocean
242    REAL, ALLOCATABLE :: phy_fsic(:) ! Fraction de glace
243    REAL, ALLOCATABLE :: phy_flic(:) ! Fraction de glace
244
245    !---------------------------------------------------------------------
246    !  Fichiers et d'autres variables
247    !---------------------------------------------------------------------
248    INTEGER :: k, l, i, it = 1, mxcalc
249    INTEGER :: nsrf
250    INTEGER jcode
251    INTEGER read_climoz
252
253    INTEGER :: it_end ! iteration number of the last call
254    !Al1
255    INTEGER ecrit_slab_oc !1=ecrit,-1=lit,0=no file
256    data ecrit_slab_oc/-1/
257
258    !     if flag_inhib_forcing = 0, tendencies of forcing are added
259    !                           <> 0, tendencies of forcing are not added
260    INTEGER :: flag_inhib_forcing = 0
261
262    !=====================================================================
263    ! INITIALIZATIONS
264    !=====================================================================
265    du_phys(:) = 0.
266    dv_phys(:) = 0.
267    dt_phys(:) = 0.
268    dt_dyn(:) = 0.
269    dt_cooling(:) = 0.
270    d_t_adv(:) = 0.
271    d_t_nudge(:) = 0.
272    d_u_nudge(:) = 0.
273    d_v_nudge(:) = 0.
274    du_adv(:) = 0.
275    dv_adv(:) = 0.
276    du_age(:) = 0.
277    dv_age(:) = 0.
278
279    ! Initialization of Common turb_forcing
280    dtime_frcg = 0.
281    Turb_fcg_gcssold = .FALSE.
282    hthturb_gcssold = 0.
283    hqturb_gcssold = 0.
284
285
286
287
288    !---------------------------------------------------------------------
289    ! OPTIONS OF THE 1D SIMULATION (lmdz1d.def => unicol.def)
290    !---------------------------------------------------------------------
291    !Al1
292    CALL conf_unicol
293    !Al1 moves this gcssold var from common fcg_gcssold to
294    Turb_fcg_gcssold = xTurb_fcg_gcssold
295    ! --------------------------------------------------------------------
296    close(1)
297    !Al1
298    WRITE(*, *) 'lmdz1d.def lu => unicol.def'
299
300    ! forcing_type defines the way the SCM is forced:
301    !forcing_type = 0 ==> forcing_les = .TRUE.
302    !             initial profiles from file prof.inp.001
303    !             no forcing by LS convergence ;
304    !             surface temperature imposed ;
305    !             radiative cooling may be imposed (iflag_radia=0 in physiq.def)
306    !forcing_type = 1 ==> forcing_radconv = .TRUE.
307    !             idem forcing_type = 0, but the imposed radiative cooling
308    !             is set to 0 (hence, if iflag_radia=0 in physiq.def,
309    !             then there is no radiative cooling at all)
310    !forcing_type = 2 ==> forcing_toga = .TRUE.
311    !             initial profiles from TOGA-COARE IFA files
312    !             LS convergence and SST imposed from TOGA-COARE IFA files
313    !forcing_type = 3 ==> forcing_GCM2SCM = .TRUE.
314    !             initial profiles from the GCM output
315    !             LS convergence imposed from the GCM output
316    !forcing_type = 4 ==> forcing_twpice = .TRUE.
317    !             initial profiles from TWP-ICE cdf file
318    !             LS convergence, omega and SST imposed from TWP-ICE files
319    !forcing_type = 5 ==> forcing_rico = .TRUE.
320    !             initial profiles from RICO files
321    !             LS convergence imposed from RICO files
322    !forcing_type = 6 ==> forcing_amma = .TRUE.
323    !             initial profiles from AMMA nc file
324    !             LS convergence, omega and surface fluxes imposed from AMMA file
325    !forcing_type = 7 ==> forcing_dice = .TRUE.
326    !             initial profiles and large scale forcings in dice_driver.nc
327    !             Different stages: soil model alone, atm. model alone
328    !             then both models coupled
329    !forcing_type = 8 ==> forcing_gabls4 = .TRUE.
330    !             initial profiles and large scale forcings in gabls4_driver.nc
331    !forcing_type >= 100 ==> forcing_case = .TRUE.
332    !             initial profiles and large scale forcings in cas.nc
333    !             LS convergence, omega and SST imposed from CINDY-DYNAMO files
334    !             101=cindynamo
335    !             102=bomex
336    !forcing_type >= 100 ==> forcing_case2 = .TRUE.
337    !             temporary flag while all the 1D cases are not whith the same cas.nc forcing file
338    !             103=arm_cu2 ie arm_cu with new forcing format
339    !             104=rico2 ie rico with new forcing format
340    !forcing_type = 40 ==> forcing_GCSSold = .TRUE.
341    !             initial profile from GCSS file
342    !             LS convergence imposed from GCSS file
343    !forcing_type = 50 ==> forcing_fire = .TRUE.
344    !             forcing from fire.nc
345    !forcing_type = 59 ==> forcing_sandu = .TRUE.
346    !             initial profiles from sanduref file: see prof.inp.001
347    !             SST varying with time and divergence constante: see ifa_sanduref.txt file
348    !             Radiation has to be computed interactively
349    !forcing_type = 60 ==> forcing_astex = .TRUE.
350    !             initial profiles from file: see prof.inp.001
351    !             SST,divergence,ug,vg,ufa,vfa varying with time : see ifa_astex.txt file
352    !             Radiation has to be computed interactively
353    !forcing_type = 61 ==> forcing_armcu = .TRUE.
354    !             initial profiles from file: see prof.inp.001
355    !             sensible and latent heat flux imposed: see ifa_arm_cu_1.txt
356    !             large scale advective forcing & radiative tendencies applied below 1000m: see ifa_arm_cu_2.txt
357    !             use geostrophic wind ug=10m/s vg=0m/s. Duration of the case 53100s
358    !             Radiation to be switched off
359
360    IF (forcing_type <=0) THEN
361      forcing_les = .TRUE.
362    elseif (forcing_type ==1) THEN
363      forcing_radconv = .TRUE.
364    elseif (forcing_type ==2) THEN
365      forcing_toga = .TRUE.
366    elseif (forcing_type ==3) THEN
367      forcing_GCM2SCM = .TRUE.
368    elseif (forcing_type ==4) THEN
369      forcing_twpice = .TRUE.
370    elseif (forcing_type ==5) THEN
371      forcing_rico = .TRUE.
372    elseif (forcing_type ==6) THEN
373      forcing_amma = .TRUE.
374    elseif (forcing_type ==7) THEN
375      forcing_dice = .TRUE.
376    elseif (forcing_type ==8) THEN
377      forcing_gabls4 = .TRUE.
378    elseif (forcing_type ==101) THEN ! Cindynamo starts 1-10-2011 0h
379      forcing_case = .TRUE.
380      year_ini_cas = 2011
381      mth_ini_cas = 10
382      day_deb = 1
383      heure_ini_cas = 0.
384      pdt_cas = 3 * 3600.         ! forcing frequency
385    elseif (forcing_type ==102) THEN ! Bomex starts 24-6-1969 0h
386      forcing_case = .TRUE.
387      year_ini_cas = 1969
388      mth_ini_cas = 6
389      day_deb = 24
390      heure_ini_cas = 0.
391      pdt_cas = 1800.         ! forcing frequency
392    elseif (forcing_type ==103) THEN ! Arm_cu starts 21-6-1997 11h30
393      forcing_case2 = .TRUE.
394      year_ini_cas = 1997
395      mth_ini_cas = 6
396      day_deb = 21
397      heure_ini_cas = 11.5
398      pdt_cas = 1800.         ! forcing frequency
399    elseif (forcing_type ==104) THEN ! rico starts 16-12-2004 0h
400      forcing_case2 = .TRUE.
401      year_ini_cas = 2004
402      mth_ini_cas = 12
403      day_deb = 16
404      heure_ini_cas = 0.
405      pdt_cas = 1800.         ! forcing frequency
406    elseif (forcing_type ==105) THEN ! bomex starts 16-12-2004 0h
407      forcing_case2 = .TRUE.
408      year_ini_cas = 1969
409      mth_ini_cas = 6
410      day_deb = 24
411      heure_ini_cas = 0.
412      pdt_cas = 1800.         ! forcing frequency
413    elseif (forcing_type ==106) THEN ! ayotte_24SC starts 6-11-1992 0h
414      forcing_case2 = .TRUE.
415      year_ini_cas = 1992
416      mth_ini_cas = 11
417      day_deb = 6
418      heure_ini_cas = 10.
419      pdt_cas = 86400.        ! forcing frequency
420    elseif (forcing_type ==113) THEN ! Arm_cu starts 21-6-1997 11h30
421      forcing_SCM = .TRUE.
422      year_ini_cas = 1997
423      ! It is possible that those parameters are run twice.
424      CALL getin('anneeref', year_ini_cas)
425      CALL getin('dayref', day_deb)
426      mth_ini_cas = 1 ! pour le moment on compte depuis le debut de l'annee
427      CALL getin('time_ini', heure_ini_cas)
428    elseif (forcing_type ==40) THEN
429      forcing_GCSSold = .TRUE.
430    elseif (forcing_type ==50) THEN
431      forcing_fire = .TRUE.
432    elseif (forcing_type ==59) THEN
433      forcing_sandu = .TRUE.
434    elseif (forcing_type ==60) THEN
435      forcing_astex = .TRUE.
436    elseif (forcing_type ==61) THEN
437      forcing_armcu = .TRUE.
438      IF(llm/=19.AND.llm/=40) stop 'Erreur nombre de niveaux !!'
439    else
440      write (*, *) 'ERROR : unknown forcing_type ', forcing_type
441      stop 'Forcing_type should be 0,1,2,3,4,5,6 or 40,59,60,61'
442    ENDIF
443    PRINT*, "forcing type=", forcing_type
444
445    ! if type_ts_forcing=0, the surface temp of 1D simulation is constant in time
446    ! (specified by tsurf in lmdz1d.def); if type_ts_forcing=1, the surface temperature
447    ! varies in time according to a forcing (e.g. forcing_toga) and is passed to read_tsurf1d.F
448    ! through the common sst_forcing.
449
450    type_ts_forcing = 0
451    IF (forcing_toga.OR.forcing_sandu.OR.forcing_astex .OR. forcing_dice)                 &
452            type_ts_forcing = 1
453
454    ! Initialization of the LOGICAL switch for nudging
455    jcode = iflag_nudge
456    do i = 1, nudge_max
457      nudge(i) = mod(jcode, 10) >= 1
458      jcode = jcode / 10
459    enddo
460    !---------------------------------------------------------------------
461    !  Definition of the run
462    !---------------------------------------------------------------------
463
464    CALL conf_gcm(99, .TRUE.)
465
466    !-----------------------------------------------------------------------
467    allocate(phy_nat (year_len))  ! 0=ocean libre,1=land,2=glacier,3=banquise
468    phy_nat(:) = 0.0
469    allocate(phy_alb (year_len))  ! Albedo land only (old value condsurf_jyg=0.3)
470    allocate(phy_sst (year_len))  ! SST (will not be used; cf read_tsurf1d.F)
471    allocate(phy_bil (year_len))  ! Ne sert que pour les slab_ocean
472    phy_bil(:) = 1.0
473    allocate(phy_rug (year_len)) ! Longueur rugosite utilisee sur land only
474    allocate(phy_ice (year_len)) ! Fraction de glace
475    phy_ice(:) = 0.0
476    allocate(phy_fter(year_len)) ! Fraction de terre
477    phy_fter(:) = 0.0
478    allocate(phy_foce(year_len)) ! Fraction de ocean
479    phy_foce(:) = 0.0
480    allocate(phy_fsic(year_len)) ! Fraction de glace
481    phy_fsic(:) = 0.0
482    allocate(phy_flic(year_len)) ! Fraction de glace
483    phy_flic(:) = 0.0
484    !-----------------------------------------------------------------------
485    !   Choix du calendrier
486    !   -------------------
487
488    !      calend = 'earth_365d'
489    IF (calend == 'earth_360d') THEN
490      CALL ioconf_calendar('360_day')
491      WRITE(*, *)'CALENDRIER CHOISI: Terrestre a 360 jours/an'
492    ELSE IF (calend == 'earth_365d') THEN
493      CALL ioconf_calendar('noleap')
494      WRITE(*, *)'CALENDRIER CHOISI: Terrestre a 365 jours/an'
495    ELSE IF (calend == 'earth_366d') THEN
496      CALL ioconf_calendar('all_leap')
497      WRITE(*, *)'CALENDRIER CHOISI: Terrestre bissextile'
498    ELSE IF (calend == 'gregorian') THEN
499      stop 'gregorian calend should not be used by normal user'
500      CALL ioconf_calendar('gregorian') ! not to be used by normal users
501      WRITE(*, *)'CALENDRIER CHOISI: Gregorien'
502    else
503      write (*, *) 'ERROR : unknown calendar ', calend
504      stop 'calend should be 360d,earth_365d,earth_366d,gregorian'
505    endif
506    !-----------------------------------------------------------------------
507
508    !c Date :
509    !      La date est supposee donnee sous la forme [annee, numero du jour dans
510    !      l annee] ; l heure est donnee dans time_ini, lu dans lmdz1d.def.
511    !      On appelle ymds2ju pour convertir [annee, jour] en [jour Julien].
512    !      Le numero du jour est dans "day". L heure est traitee separement.
513    !      La date complete est dans "daytime" (l'unite est le jour).
514    IF (nday>0) THEN
515      fnday = nday
516    else
517      fnday = -nday / float(day_step)
518    endif
519    print *, 'fnday=', fnday
520    !     start_time doit etre en FRACTION DE JOUR
521    start_time = time_ini / 24.
522
523    ! Special case for arm_cu which lasts less than one day : 53100s !! (MPL 20111026)
524    IF(forcing_type == 61) fnday = 53100. / 86400.
525    IF(forcing_type == 103) fnday = 53100. / 86400.
526    ! Special case for amma which lasts less than one day : 64800s !! (MPL 20120216)
527    IF(forcing_type == 6) fnday = 64800. / 86400.
528    !     IF(forcing_type .EQ. 6) fnday=50400./86400.
529    IF(forcing_type == 8) fnday = 129600. / 86400.
530    annee_ref = anneeref
531    mois = 1
532    day_ref = dayref
533    heure = 0.
534    itau_dyn = 0
535    itau_phy = 0
536    CALL ymds2ju(annee_ref, mois, day_ref, heure, day)
537    day_ini = int(day)
538    day_end = day_ini + int(fnday)
539
540    IF (forcing_type ==2) THEN
541      ! Convert the initial date of Toga-Coare to Julian day
542      CALL ymds2ju                                                          &
543              (year_ini_toga, mth_ini_toga, day_ini_toga, heure, day_ju_ini_toga)
544
545    ELSEIF (forcing_type ==4) THEN
546      ! Convert the initial date of TWPICE to Julian day
547      CALL ymds2ju                                                          &
548              (year_ini_twpi, mth_ini_twpi, day_ini_twpi, heure_ini_twpi              &
549              , day_ju_ini_twpi)
550    ELSEIF (forcing_type ==6) THEN
551      ! Convert the initial date of AMMA to Julian day
552      CALL ymds2ju                                                          &
553              (year_ini_amma, mth_ini_amma, day_ini_amma, heure_ini_amma              &
554              , day_ju_ini_amma)
555    ELSEIF (forcing_type ==7) THEN
556      ! Convert the initial date of DICE to Julian day
557      CALL ymds2ju                                                         &
558              (year_ini_dice, mth_ini_dice, day_ini_dice, heure_ini_dice             &
559              , day_ju_ini_dice)
560    ELSEIF (forcing_type ==8) THEN
561      ! Convert the initial date of GABLS4 to Julian day
562      CALL ymds2ju                                                         &
563              (year_ini_gabls4, mth_ini_gabls4, day_ini_gabls4, heure_ini_gabls4     &
564              , day_ju_ini_gabls4)
565    ELSEIF (forcing_type >100) THEN
566      ! Convert the initial date to Julian day
567      day_ini_cas = day_deb
568      PRINT*, 'time case', year_ini_cas, mth_ini_cas, day_ini_cas
569      CALL ymds2ju                                                         &
570              (year_ini_cas, mth_ini_cas, day_ini_cas, heure_ini_cas * 3600            &
571              , day_ju_ini_cas)
572      PRINT*, 'time case 2', day_ini_cas, day_ju_ini_cas
573    ELSEIF (forcing_type ==59) THEN
574      ! Convert the initial date of Sandu case to Julian day
575      CALL ymds2ju                                                          &
576              (year_ini_sandu, mth_ini_sandu, day_ini_sandu, &
577              time_ini * 3600., day_ju_ini_sandu)
578
579    ELSEIF (forcing_type ==60) THEN
580      ! Convert the initial date of Astex case to Julian day
581      CALL ymds2ju                                                          &
582              (year_ini_astex, mth_ini_astex, day_ini_astex, &
583              time_ini * 3600., day_ju_ini_astex)
584
585    ELSEIF (forcing_type ==61) THEN
586      ! Convert the initial date of Arm_cu case to Julian day
587      CALL ymds2ju                                                          &
588              (year_ini_armcu, mth_ini_armcu, day_ini_armcu, heure_ini_armcu          &
589              , day_ju_ini_armcu)
590    ENDIF
591
592    IF (forcing_type >100) THEN
593      daytime = day + heure_ini_cas / 24. ! 1st day and initial time of the simulation
594    ELSE
595      daytime = day + time_ini / 24. ! 1st day and initial time of the simulation
596    ENDIF
597    ! Print out the actual date of the beginning of the simulation :
598    CALL ju2ymds(daytime, year_print, month_print, day_print, sec_print)
599    print *, ' Time of beginning : ', &
600            year_print, month_print, day_print, sec_print
601
602    !---------------------------------------------------------------------
603    ! Initialization of dimensions, geometry and initial state
604    !---------------------------------------------------------------------
605    !      CALL init_phys_lmdz(1,1,llm,1,(/1/)) ! job now done via iniphysiq
606    !     but we still need to initialize dimphy module (klon,klev,etc.)  here.
607    CALL init_dimphy1D(1, llm)
608    CALL suphel
609    CALL init_infotrac
610
611    IF (nqtot>nqmx) STOP 'Augmenter nqmx dans lmdz1d.F'
612    allocate(q(llm, nqtot)) ; q(:, :) = 0.
613    allocate(dq(llm, nqtot))
614    allocate(dq_dyn(llm, nqtot))
615    allocate(d_q_adv(llm, nqtot))
616    allocate(d_q_nudge(llm, nqtot))
617    !      allocate(d_th_adv(llm))
618
619    q(:, :) = 0.
620    dq(:, :) = 0.
621    dq_dyn(:, :) = 0.
622    d_q_adv(:, :) = 0.
623    d_q_nudge(:, :) = 0.
624
625    !   No ozone climatology need be read in this pre-initialization
626    !          (phys_state_var_init is called again in physiq)
627    read_climoz = 0
628    nsw = 6          ! EV et LF: sinon, falb_dir et falb_dif ne peuvent etre alloues
629
630    CALL phys_state_var_init(read_climoz)
631
632    IF (ngrid/=klon) THEN
633      PRINT*, 'stop in inifis'
634      PRINT*, 'Probleme de dimensions :'
635      PRINT*, 'ngrid = ', ngrid
636      PRINT*, 'klon  = ', klon
637      stop
638    endif
639    !!!=====================================================================
640    !!! Feedback forcing values for Gateaux differentiation (al1)
641    !!!=====================================================================
642    !!! Surface Planck forcing bracketing CALL radiation
643    !!      surf_Planck = 0.
644    !!      surf_Conv   = 0.
645    !!      WRITE(*,*) 'Gateaux-dif Planck,Conv:',surf_Planck,surf_Conv
646    !!! a mettre dans le lmdz1d.def ou autre
647    !!
648    !!
649    qsol = qsolinp
650    qsurf = fq_sat(tsurf, psurf / 100.)
651    beta_surf = 1.
652    beta_aridity(:, :) = beta_surf
653    day1 = day_ini
654    time = daytime - day
655    ts_toga(1) = tsurf ! needed by read_tsurf1d.F
656    rho(1) = psurf / (rd * tsurf * (1. + (rv / rd - 1.) * qsurf))
657
658    !! mpl et jyg le 22/08/2012 :
659    !!  pour que les cas a flux de surface imposes marchent
660    IF(.NOT.ok_flux_surf.OR.max(abs(wtsurf), abs(wqsurf))>0.) THEN
661      fsens = -wtsurf * rcpd * rho(1)
662      flat = -wqsurf * rlvtt * rho(1)
663      print *, 'Flux: ok_flux wtsurf wqsurf', ok_flux_surf, wtsurf, wqsurf
664    ENDIF
665    PRINT*, 'Flux sol ', fsens, flat
666    !!      ok_flux_surf=.FALSE.
667    !!      fsens=-wtsurf*rcpd*rho(1)
668    !!      flat=-wqsurf*rlvtt*rho(1)
669    !!!!
670
671    ! Vertical discretization and pressure levels at half and mid levels:
672
673    pa = 5e4
674    !!      preff= 1.01325e5
675    preff = psurf
676    IF (ok_old_disvert) THEN
677      CALL disvert0(pa, preff, ap, bp, dpres, presnivs, nivsigs, nivsig)
678      print *, 'On utilise disvert0'
679      aps(1:llm) = 0.5 * (ap(1:llm) + ap(2:llm + 1))
680      bps(1:llm) = 0.5 * (bp(1:llm) + bp(2:llm + 1))
681      scaleheight = 8.
682      pseudoalt(1:llm) = -scaleheight * log(presnivs(1:llm) / preff)
683    ELSE
684      CALL disvert()
685      print *, 'On utilise disvert'
686      !       Nouvelle version disvert permettant d imposer ap,bp (modif L.Guez) MPL 18092012
687      !       Dans ce cas, on lit ap,bp dans le fichier hybrid.txt
688    ENDIF
689
690    sig_s = presnivs / preff
691    plev = ap + bp * psurf
692    play = 0.5 * (plev(1:llm) + plev(2:llm + 1))
693    zlay = -rd * 300. * log(play / psurf) / rg ! moved after reading profiles
694
695    IF (forcing_type == 59) THEN
696      ! pour forcing_sandu, on cherche l'indice le plus proche de 700hpa#3000m
697      WRITE(*, *) '***********************'
698      do l = 1, llm
699        WRITE(*, *) 'l,play(l),presnivs(l): ', l, play(l), presnivs(l)
700        IF (trouve_700 .AND. play(l)<=70000) THEN
701          llm700 = l
702          print *, 'llm700,play=', llm700, play(l) / 100.
703          trouve_700 = .FALSE.
704        endif
705      enddo
706      WRITE(*, *) '***********************'
707    ENDIF
708
709    !=====================================================================
710    ! EVENTUALLY, READ FORCING DATA :
711    !=====================================================================
712
713    INCLUDE "old_1D_read_forc_cases.h"
714
715  IF (forcing_GCM2SCM) THEN
716  write (*, *) 'forcing_GCM2SCM not yet implemented'
717  stop 'in initialization'
718END IF ! forcing_GCM2SCM
719
720        PRINT*, 'mxcalc=', mxcalc
721        !     PRINT*,'zlay=',zlay(mxcalc)
722PRINT*, 'play=', play(mxcalc)
723
724                !Al1 pour SST forced, appell?? depuis ocean_forced_noice
725                ! EV tg instead of ts_cur
726
727                tg = tsurf ! SST used in read_tsurf1d
728                !=====================================================================
729                ! Initialisation de la physique :
730                !=====================================================================
731
732                !  Rq: conf_phys.F90 lit tous les flags de physiq.def; conf_phys appele depuis physiq.F
733
734                ! day_step, iphysiq lus dans gcm.def ci-dessus
735                ! timestep: calcule ci-dessous from rday et day_step
736                ! ngrid=1
737                ! llm: defini dans .../modipsl/modeles/LMDZ4/libf/grid/dimension
738                ! rday: defini dans suphel.F (86400.)
739                ! day_ini: lu dans run.def (dayref)
740                ! rlat_rad,rlon-rad: transformes en radian de rlat,rlon lus dans lmdz1d.def (en degres)
741                ! airefi,zcufi,zcvfi initialises au debut de ce programme
742                ! rday,ra,rg,rd,rcpd declares dans YOMCST.h et calcules dans suphel.F
743                day_step = float(nsplit_phys) * day_step / float(iphysiq)
744                write (*, *) 'Time step divided by nsplit_phys (=', nsplit_phys, ')'
745        timestep = rday / day_step
746        dtime_frcg = timestep
747
748        zcufi = airefi
749                zcvfi = airefi
750
751                rlat_rad(1) = xlat * rpi / 180.
752                rlon_rad(1) = xlon * rpi / 180.
753
754                ! iniphysiq will CALL iniaqua who needs year_len from phys_cal_mod
755                year_len_phys_cal_mod = year_len
756
757                ! Ehouarn: iniphysiq requires arrays related to (3D) dynamics grid,
758                ! e.g. for cell boundaries, which are meaningless in 1D; so pad these
759                ! with '0.' when necessary
760                CALL iniphysiq(iim, jjm, llm, &
761                1, comm_lmdz, &
762                rday, day_ini, timestep, &
763                (/rlat_rad(1), 0./), (/0./), &
764                (/0., 0./), (/rlon_rad(1), 0./), &
765                (/ (/airefi, 0./), (/0., 0./) /), &
766                (/zcufi, 0., 0., 0./), &
767                (/zcvfi, 0./), &
768                ra, rg, rd, rcpd, 1)
769                PRINT*, 'apres iniphysiq'
770
771                ! 2 PARAMETRES QUI DEVRAIENT ETRE LUS DANS run.def MAIS NE LE SONT PAS ICI:
772                co2_ppm = 330.0
773        solaire = 1370.0
774
775        ! Ecriture du startphy avant le premier appel a la physique.
776        ! On le met juste avant pour avoir acces a tous les champs
777
778        IF (ok_writedem) THEN
779        !--------------------------------------------------------------------------
780        ! pbl_surface_init (called here) and pbl_surface_final (called by phyredem)
781        ! need : qsol fder snow qsurf evap rugos agesno ftsoil
782        !--------------------------------------------------------------------------
783
784        type_ocean = "force"
785                run_off_lic_0(1) = restart_runoff
786        CALL fonte_neige_init(run_off_lic_0)
787
788                fder = 0.
789                snsrf(1, :) = snowmass ! masse de neige des sous surface
790        print *, 'snsrf', snsrf
791        qsurfsrf(1, :) = qsurf ! humidite de l'air des sous surface
792        fevap = 0.
793        z0m(1, :) = rugos     ! couverture de neige des sous surface
794        z0h(1, :) = rugosh    ! couverture de neige des sous surface
795        agesno = xagesno
796        tsoil(:, :, :) = tsurf
797        !------ AMMA 2e run avec modele sol et rayonnement actif (MPL 23052012)
798        !       tsoil(1,1,1)=299.18
799        !       tsoil(1,2,1)=300.08
800        !       tsoil(1,3,1)=301.88
801        !       tsoil(1,4,1)=305.48
802        !       tsoil(1,5,1)=308.00
803        !       tsoil(1,6,1)=308.00
804                !       tsoil(1,7,1)=308.00
805                !       tsoil(1,8,1)=308.00
806                !       tsoil(1,9,1)=308.00
807                !       tsoil(1,10,1)=308.00
808                !       tsoil(1,11,1)=308.00
809                !-----------------------------------------------------------------------
810                CALL pbl_surface_init(fder, snsrf, qsurfsrf, tsoil)
811
812        !------------------ prepare limit conditions for limit.nc -----------------
813        !--   Ocean force
814
815        PRINT*, 'avant phyredem'
816        pctsrf(1, :) = 0.
817        IF (nat_surf==0.) THEN
818        pctsrf(1, is_oce) = 1.
819        pctsrf(1, is_ter) = 0.
820        pctsrf(1, is_lic) = 0.
821        pctsrf(1, is_sic) = 0.
822        ELSE IF (nat_surf == 1) THEN
823        pctsrf(1, is_oce) = 0.
824        pctsrf(1, is_ter) = 1.
825        pctsrf(1, is_lic) = 0.
826        pctsrf(1, is_sic) = 0.
827        ELSE IF (nat_surf == 2) THEN
828        pctsrf(1, is_oce) = 0.
829        pctsrf(1, is_ter) = 0.
830        pctsrf(1, is_lic) = 1.
831        pctsrf(1, is_sic) = 0.
832        ELSE IF (nat_surf == 3) THEN
833        pctsrf(1, is_oce) = 0.
834        pctsrf(1, is_ter) = 0.
835        pctsrf(1, is_lic) = 0.
836        pctsrf(1, is_sic) = 1.
837
838        end if
839
840        PRINT*, 'nat_surf,pctsrf(1,is_oce),pctsrf(1,is_ter)', nat_surf         &
841        , pctsrf(1, is_oce), pctsrf(1, is_ter)
842
843                zmasq = pctsrf(1, is_ter) + pctsrf(1, is_lic)
844        zpic = zpicinp
845        ftsol = tsurf
846        falb_dir = albedo
847        falb_dif = albedo
848        rugoro = rugos
849                t_ancien(1, :) = temp(:)
850                q_ancien(1, :) = q(:, 1)
851        ql_ancien = 0.
852        qs_ancien = 0.
853        prlw_ancien = 0.
854        prsw_ancien = 0.
855        prw_ancien = 0.
856                !jyg<
857                !!        pbl_tke(:,:,:)=1.e-8
858                pbl_tke(:, :, :) = 0.
859                pbl_tke(:, 2, :) = 1.e-2
860                PRINT *, ' pbl_tke dans lmdz1d '
861                IF (prt_level >= 5) THEN
862                DO nsrf = 1, 4
863                PRINT *, 'pbl_tke(1,:,', nsrf, ') ', pbl_tke(1, :, nsrf)
864        ENDDO
865        end if
866
867        !>jyg
868
869        rain_fall = 0.
870        snow_fall = 0.
871        solsw = 0.
872        solswfdiff = 0.
873        sollw = 0.
874        sollwdown = rsigma * tsurf**4
875        radsol = 0.
876        rnebcon = 0.
877        ratqs = 0.
878        clwcon = 0.
879                zmax0 = 0.
880                zmea = zsurf
881                zstd = 0.
882        zsig = 0.
883        zgam = 0.
884                zval = 0.
885                zthe = 0.
886                sig1 = 0.
887        w01 = 0.
888
889        wake_deltaq = 0.
890                wake_deltat = 0.
891                wake_delta_pbl_TKE(:, :, :) = 0.
892        delta_tsurf = 0.
893        wake_fip = 0.
894                wake_pe = 0.
895                wake_s = 0.
896                awake_s = 0.
897                wake_dens = 0.
898                awake_dens = 0.
899                cv_gen = 0.
900                wake_cstar = 0.
901                ale_bl = 0.
902                ale_bl_trig = 0.
903                alp_bl = 0.
904                IF (ALLOCATED(du_gwd_rando)) du_gwd_rando = 0.
905                IF (ALLOCATED(du_gwd_front)) du_gwd_front = 0.
906                entr_therm = 0.
907                detr_therm = 0.
908        f0 = 0.
909        fm_therm = 0.
910        u_ancien(1, :) = u(:)
911                v_ancien(1, :) = v(:)
912                rneb_ancien(1, :) = 0.
913
914        u10m = 0.
915        v10m = 0.
916                ale_wake = 0.
917                ale_bl_stat = 0.
918
919                !------------------------------------------------------------------------
920                ! Make file containing restart for the physics (startphy.nc)
921
922                ! NB: List of the variables to be written by phyredem (via put_field):
923                ! rlon,rlat,zmasq,pctsrf(:,is_ter),pctsrf(:,is_lic),pctsrf(:,is_oce)
924                ! pctsrf(:,is_sic),ftsol(:,nsrf),tsoil(:,isoil,nsrf),qsurf(:,nsrf)
925                ! qsol,falb_dir(:,nsrf),falb_dif(:,nsrf),evap(:,nsrf),snow(:,nsrf)
926                ! radsol,solsw,solswfdiff,sollw, sollwdown,fder,rain_fall,snow_fall,frugs(:,nsrf)
927                ! agesno(:,nsrf),zmea,zstd,zsig,zgam,zthe,zpic,zval,rugoro
928                ! t_ancien,q_ancien,,frugs(:,is_oce),clwcon(:,1),rnebcon(:,1),ratqs(:,1)
929                ! run_off_lic_0,pbl_tke(:,1:klev,nsrf), zmax0,f0,sig1,w01
930                ! wake_deltat,wake_deltaq,wake_s,awake_s,wake_dens,awake_dens,cv_gen,wake_cstar,
931                ! wake_fip,wake_delta_pbl_tke(:,1:klev,nsrf)
932
933                ! NB2: The content of the startphy.nc file depends on some flags defined in
934                ! the ".def" files. However, since conf_phys is not called in lmdz1d.F90, these flags have
935                ! to be set at some arbitratry convenient values.
936                !------------------------------------------------------------------------
937                !Al1 =============== restart option ==========================
938                iflag_physiq = 0
939                CALL getin('iflag_physiq', iflag_physiq)
940
941                IF (.NOT.restart) THEN
942                iflag_pbl = 5
943                CALL phyredem ("startphy.nc")
944        else
945        ! (desallocations)
946        PRINT*, 'callin surf final'
947        CALL pbl_surface_final(fder, snsrf, qsurfsrf, tsoil)
948                PRINT*, 'after surf final'
949                CALL fonte_neige_final(run_off_lic_0)
950                endif
951
952                ok_writedem = .FALSE.
953                PRINT*, 'apres phyredem'
954
955        END IF ! ok_writedem
956
957        !------------------------------------------------------------------------
958        ! Make file containing boundary conditions (limit.nc) **Al1->restartdyn***
959        ! --------------------------------------------------
960        ! NB: List of the variables to be written in limit.nc
961        !     (by writelim.F, SUBROUTINE of 1DUTILS.h):
962                !        phy_nat,phy_alb,phy_sst,phy_bil,phy_rug,phy_ice,
963                !        phy_fter,phy_foce,phy_flic,phy_fsic)
964                !------------------------------------------------------------------------
965                DO i = 1, year_len
966                phy_nat(i) = nat_surf
967        phy_alb(i) = albedo
968        phy_sst(i) = tsurf ! read_tsurf1d will be used instead
969        phy_rug(i) = rugos
970        phy_fter(i) = pctsrf(1, is_ter)
971        phy_foce(i) = pctsrf(1, is_oce)
972                phy_fsic(i) = pctsrf(1, is_sic)
973                phy_flic(i) = pctsrf(1, is_lic)
974                END DO
975
976                ! fabrication de limit.nc
977                CALL writelim (1, phy_nat, phy_alb, phy_sst, phy_bil, phy_rug, &
978                phy_ice, phy_fter, phy_foce, phy_flic, phy_fsic)
979
980                CALL phys_state_var_end
981                !Al1
982                IF (restart) THEN
983                PRINT*, 'CALL to restart dyn 1d'
984                Call dyn1deta0("start1dyn.nc", plev, play, phi, phis, presnivs, &
985                u, v, temp, q, omega2)
986
987                PRINT*, 'fnday,annee_ref,day_ref,day_ini', &
988                fnday, annee_ref, day_ref, day_ini
989                !**      CALL ymds2ju(annee_ref,mois,day_ini,heure,day)
990                day = day_ini
991                day_end = day_ini + nday
992        daytime = day + time_ini / 24. ! 1st day and initial time of the simulation
993
994        ! Print out the actual date of the beginning of the simulation :
995        CALL ju2ymds(daytime, an, mois, jour, heure)
996                print *, ' Time of beginning : y m d h', an, mois, jour, heure / 3600.
997
998                day = int(daytime)
999                time = daytime - day
1000
1001                PRINT*, '****** intialised fields from restart1dyn *******'
1002                PRINT*, 'plev,play,phi,phis,presnivs,u,v,temp,q,omega2'
1003                PRINT*, 'temp(1),q(1,1),u(1),v(1),plev(1),phis :'
1004                PRINT*, temp(1), q(1, 1), u(1), v(1), plev(1), phis
1005        ! raz for safety
1006        do l = 1, llm
1007                dq_dyn(l, 1) = 0.
1008        enddo
1009        END IF
1010        !Al1 ================  end restart =================================
1011        IF (ecrit_slab_oc==1) THEN
1012        open(97, file = 'div_slab.dat', STATUS = 'UNKNOWN')
1013                elseif (ecrit_slab_oc==0) THEN
1014                open(97, file = 'div_slab.dat', STATUS = 'OLD')
1015                END IF
1016
1017                !---------------------------------------------------------------------
1018                !    Initialize target profile for RHT nudging if needed
1019                !---------------------------------------------------------------------
1020                IF (nudge(inudge_RHT)) THEN
1021        CALL nudge_RHT_init(plev, play, temp, q(:, 1), t_targ, rh_targ)
1022                END IF
1023                IF (nudge(inudge_UV)) THEN
1024                CALL nudge_UV_init(plev, play, u, v, u_targ, v_targ)
1025                END IF
1026
1027                !=====================================================================
1028                IF (CPPKEY_OUTPUTPHYSSCM) THEN
1029                CALL iophys_ini(timestep)
1030        END IF
1031        ! START OF THE TEMPORAL LOOP :
1032        !=====================================================================
1033
1034        it_end = nint(fnday * day_step)
1035                !test JLD     it_end = 10
1036                DO while(it<=it_end)
1037
1038                IF (prt_level>=1) THEN
1039        PRINT*, 'XXXXXXXXXXXXXXXXXXX ITAP,day,time=', &
1040        it, day, time, it_end, day_step
1041        PRINT*, 'PAS DE TEMPS ', timestep
1042        endif
1043        !Al1 demande de restartphy.nc
1044        IF (it==it_end) lastcall = .True.
1045
1046        !---------------------------------------------------------------------
1047        !  Geopotential :
1048        !---------------------------------------------------------------------
1049
1050        phi(1) = RD * temp(1) * (plev(1) - play(1)) / (.5 * (plev(1) + play(1)))
1051                do l = 1, llm - 1
1052                phi(l + 1) = phi(l) + RD * (temp(l) + temp(l + 1)) * &
1053        (play(l) - play(l + 1)) / (play(l) + play(l + 1))
1054                enddo
1055
1056                !---------------------------------------------------------------------
1057                ! Interpolation of forcings in time and onto model levels
1058                !---------------------------------------------------------------------
1059
1060                INCLUDE "old_1D_interp_cases.h"
1061
1062                IF (forcing_GCM2SCM) THEN
1063        write (*, *) 'forcing_GCM2SCM not yet implemented'
1064        stop 'in time loop'
1065        END IF ! forcing_GCM2SCM
1066
1067        !!!!---------------------------------------------------------------------
1068        !!!!  Geopotential :
1069                !!!!---------------------------------------------------------------------
1070                !!!
1071                !!!        phi(1)=RD*temp(1)*(plev(1)-play(1))/(.5*(plev(1)+play(1)))
1072                !!!        do l = 1, llm-1
1073                !!!          phi(l+1)=phi(l)+RD*(temp(l)+temp(l+1))*                           &
1074                !!!     &    (play(l)-play(l+1))/(play(l)+play(l+1))
1075                !!!        enddo
1076
1077                !---------------------------------------------------------------------
1078                ! Listing output for debug prt_level>=1
1079                !---------------------------------------------------------------------
1080                IF (prt_level>=1) THEN
1081                print *, ' avant physiq : -------- day time ', day, time
1082                WRITE(*, *) 'firstcall,lastcall,phis', &
1083                firstcall, lastcall, phis
1084                end if
1085                IF (prt_level>=5) THEN
1086        WRITE(*, '(a10,2a4,4a13)') 'BEFOR1 IT=', 'it', 'l', &
1087        'presniv', 'plev', 'play', 'phi'
1088        WRITE(*, '(a10,2i4,4f13.2)') ('BEFOR1 IT= ', it, l, &
1089        presnivs(l), plev(l), play(l), phi(l), l = 1, llm)
1090                WRITE(*, '(a11,2a4,a11,6a8)') 'BEFOR2', 'it', 'l', &
1091                'presniv', 'u', 'v', 'temp', 'q1', 'q2', 'omega2'
1092        WRITE(*, '(a11,2i4,f11.2,5f8.2,e10.2)') ('BEFOR2 IT= ', it, l, &
1093        presnivs(l), u(l), v(l), temp(l), q(l, 1), q(l, 2), omega2(l), l = 1, llm)
1094                END IF
1095
1096                !---------------------------------------------------------------------
1097                !   Call physiq :
1098                !---------------------------------------------------------------------
1099                CALL physiq(ngrid, llm, &
1100        firstcall, lastcall, timestep, &
1101        plev, play, phi, phis, presnivs, &
1102        u, v, rot, temp, q, omega2, &
1103        du_phys, dv_phys, dt_phys, dq, dpsrf)
1104                firstcall = .FALSE.
1105
1106                !---------------------------------------------------------------------
1107                ! Listing output for debug
1108                !---------------------------------------------------------------------
1109                IF (prt_level>=5) THEN
1110        WRITE(*, '(a11,2a4,4a13)') 'AFTER1 IT=', 'it', 'l', &
1111        'presniv', 'plev', 'play', 'phi'
1112        WRITE(*, '(a11,2i4,4f13.2)') ('AFTER1 it= ', it, l, &
1113        presnivs(l), plev(l), play(l), phi(l), l = 1, llm)
1114                WRITE(*, '(a11,2a4,a11,6a8)') 'AFTER2', 'it', 'l', &
1115                'presniv', 'u', 'v', 'temp', 'q1', 'q2', 'omega2'
1116                WRITE(*, '(a11,2i4,f11.2,5f8.2,e10.2)') ('AFTER2 it= ', it, l, &
1117                presnivs(l), u(l), v(l), temp(l), q(l, 1), q(l, 2), omega2(l), l = 1, llm)
1118        WRITE(*, '(a11,2a4,a11,5a8)') 'AFTER3', 'it', 'l', &
1119        'presniv', 'du_phys', 'dv_phys', 'dt_phys', 'dq1', 'dq2'
1120        WRITE(*, '(a11,2i4,f11.2,5f8.2)') ('AFTER3 it= ', it, l, &
1121        presnivs(l), 86400 * du_phys(l), 86400 * dv_phys(l), &
1122        86400 * dt_phys(l), 86400 * dq(l, 1), dq(l, 2), l = 1, llm)
1123                WRITE(*, *) 'dpsrf', dpsrf
1124                END IF
1125                !---------------------------------------------------------------------
1126                !   Add physical tendencies :
1127                !---------------------------------------------------------------------
1128
1129                fcoriolis = 2. * sin(rpi * xlat / 180.) * romega
1130        IF (forcing_radconv .OR. forcing_fire) THEN
1131        fcoriolis = 0.0
1132        dt_cooling = 0.0
1133                d_t_adv = 0.0
1134                d_q_adv = 0.0
1135                END IF
1136                !      PRINT*, 'calcul de fcoriolis ', fcoriolis
1137
1138                IF (forcing_toga .OR. forcing_GCSSold .OR. forcing_twpice            &
1139                .OR.forcing_amma .OR. forcing_type==101) THEN
1140                fcoriolis = 0.0 ; ug = 0. ; vg = 0.
1141                END IF
1142
1143                IF(forcing_rico) THEN
1144                dt_cooling = 0.
1145                END IF
1146
1147                !CRio:Attention modif sp??cifique cas de Caroline
1148IF (forcing_type==-1) THEN
1149        fcoriolis = 0.
1150        !Nudging
1151
1152        !on calcule dt_cooling
1153        do l = 1, llm
1154        IF (play(l)>=20000.) THEN
1155        dt_cooling(l) = -1.5 / 86400.
1156        elseif ((play(l)>=10000.).AND.((play(l)<20000.))) THEN
1157        dt_cooling(l) = -1.5 / 86400. * (play(l) - 10000.) / (10000.) - 1. / 86400. * (20000. - play(l)) / 10000. * (temp(l) - 200.)
1158                else
1159                dt_cooling(l) = -1. * (temp(l) - 200.) / 86400.
1160        endif
1161        enddo
1162
1163        END IF
1164                !RC
1165                IF (forcing_sandu) THEN
1166                ug(1:llm) = u_mod(1:llm)
1167                vg(1:llm) = v_mod(1:llm)
1168                END IF
1169
1170                IF (prt_level >= 5) PRINT*, 'fcoriolis, xlat,mxcalc ', &
1171                fcoriolis, xlat, mxcalc
1172
1173                !       print *,'u-ug=',u-ug
1174
1175                !!!!!!!!!!!!!!!!!!!!!!!!
1176                ! Geostrophic wind
1177                ! Le calcul ci dessous est insuffisamment precis
1178                !      du_age(1:mxcalc)=fcoriolis*(v(1:mxcalc)-vg(1:mxcalc))
1179                !      dv_age(1:mxcalc)=-fcoriolis*(u(1:mxcalc)-ug(1:mxcalc))
1180                !!!!!!!!!!!!!!!!!!!!!!!!
1181                sfdt = sin(0.5 * fcoriolis * timestep)
1182                cfdt = cos(0.5 * fcoriolis * timestep)
1183        !       print *,'fcoriolis,sfdt,cfdt,timestep',fcoriolis,sfdt,cfdt,timestep
1184
1185        du_age(1:mxcalc) = -2. * sfdt / timestep * &
1186        (sfdt * (u(1:mxcalc) - ug(1:mxcalc)) - &
1187                cfdt * (v(1:mxcalc) - vg(1:mxcalc)))
1188                !!     : fcoriolis*(v(1:mxcalc)-vg(1:mxcalc))
1189
1190                dv_age(1:mxcalc) = -2. * sfdt / timestep * &
1191        (cfdt * (u(1:mxcalc) - ug(1:mxcalc)) + &
1192                sfdt * (v(1:mxcalc) - vg(1:mxcalc)))
1193                !!     : -fcoriolis*(u(1:mxcalc)-ug(1:mxcalc))
1194
1195                !!!!!!!!!!!!!!!!!!!!!!!!
1196                !  Nudging
1197                !!!!!!!!!!!!!!!!!!!!!!!!
1198                d_t_nudge(:) = 0.
1199        d_q_nudge(:, :) = 0.
1200        d_u_nudge(:) = 0.
1201        d_v_nudge(:) = 0.
1202        IF (nudge(inudge_RHT)) THEN
1203        CALL nudge_RHT(timestep, plev, play, t_targ, rh_targ, temp, q(:, 1), &
1204        d_t_nudge, d_q_nudge(:, 1))
1205        END IF
1206        IF (nudge(inudge_UV)) THEN
1207        CALL nudge_UV(timestep, plev, play, u_targ, v_targ, u, v, &
1208        d_u_nudge, d_v_nudge)
1209                END IF
1210
1211                IF (forcing_fire) THEN
1212
1213        !let ww=if ( alt le 1100 ) then alt*-0.00001 else 0
1214        !let wt=if ( alt le 1100 ) then min( -3.75e-5 , -7.5e-8*alt)  else 0
1215        !let wq=if ( alt le 1100 ) then max( 1.5e-8 , 3e-11*alt)  else 0
1216        d_t_adv = 0.
1217        d_q_adv = 0.
1218        teta = temp * (pzero / play)**rkappa
1219        d_t_adv = 0.
1220        d_q_adv = 0.
1221        do l = 2, llm - 1
1222        IF (zlay(l)<=1100) THEN
1223        wwww = -0.00001 * zlay(l)
1224                d_t_adv(l) = -wwww * (teta(l) - teta(l + 1)) / (zlay(l) - zlay(l + 1)) / (pzero / play(l))**rkappa
1225        d_q_adv(l, 1:2) = -wwww * (q(l, 1:2) - q(l + 1, 1:2)) / (zlay(l) - zlay(l + 1))
1226                d_t_adv(l) = d_t_adv(l) + min(-3.75e-5, -7.5e-8 * zlay(l))
1227                d_q_adv(l, 1) = d_q_adv(l, 1) + max(1.5e-8, 3e-11 * zlay(l))
1228        endif
1229        enddo
1230
1231        END IF
1232
1233        !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
1234        !         call  writefield_phy('dv_age' ,dv_age,llm)
1235!         call  writefield_phy('du_age' ,du_age,llm)
1236        !         call  writefield_phy('du_phys' ,du_phys,llm)
1237        !         call  writefield_phy('u_tend' ,u,llm)
1238        !         call  writefield_phy('u_g' ,ug,llm)
1239
1240        !!!!!!!!!!!!!!!!!!!!!!!!!!!!!
1241        !! Increment state variables
1242        !!!!!!!!!!!!!!!!!!!!!!!!!!!!!
1243        IF (flag_inhib_forcing == 0) then ! if tendency of forcings should be added
1244
1245        ! pour les cas sandu et astex, on reclacule u,v,q,temp et teta dans 1D_nudge_sandu_astex.h
1246        ! au dessus de 700hpa, on relaxe vers les profils initiaux
1247        IF (forcing_sandu .OR. forcing_astex) THEN
1248        INCLUDE "1D_nudge_sandu_astex.h"
1249        else
1250        u(1:mxcalc) = u(1:mxcalc) + timestep * (&
1251        du_phys(1:mxcalc)                                       &
1252        + du_age(1:mxcalc) + du_adv(1:mxcalc)                       &
1253        + d_u_nudge(1:mxcalc))
1254                v(1:mxcalc) = v(1:mxcalc) + timestep * (&
1255                dv_phys(1:mxcalc)                                       &
1256                + dv_age(1:mxcalc) + dv_adv(1:mxcalc)                       &
1257                + d_v_nudge(1:mxcalc))
1258                        q(1:mxcalc, :) = q(1:mxcalc, :) + timestep * (&
1259        dq(1:mxcalc, :)                                        &
1260        + d_q_adv(1:mxcalc, :)                                   &
1261        + d_q_nudge(1:mxcalc, :))
1262
1263                IF (prt_level>=3) THEN
1264                print *, &
1265                'physiq-> temp(1),dt_phys(1),d_t_adv(1),dt_cooling(1) ', &
1266                temp(1), dt_phys(1), d_t_adv(1), dt_cooling(1)
1267                PRINT*, 'dv_phys=', dv_phys
1268                PRINT*, 'dv_age=', dv_age
1269                PRINT*, 'dv_adv=', dv_adv
1270                PRINT*, 'd_v_nudge=', d_v_nudge
1271                PRINT*, v
1272                PRINT*, vg
1273                endif
1274
1275                temp(1:mxcalc) = temp(1:mxcalc) + timestep * (&
1276        dt_phys(1:mxcalc)                                       &
1277        + d_t_adv(1:mxcalc)                                      &
1278        + d_t_nudge(1:mxcalc)                                      &
1279        + dt_cooling(1:mxcalc))  ! Taux de chauffage ou refroid.
1280
1281        IF (CPPKEY_OUTPUTPHYSSCM) THEN
1282        CALL iophys_ecrit('d_t_adv', klev, 'd_t_adv', 'm/s', d_t_adv)
1283                CALL iophys_ecrit('d_t_nudge', klev, 'd_t_nudge', 'm/s', d_t_nudge)
1284                END IF
1285
1286                endif  ! forcing_sandu or forcing_astex
1287
1288                teta = temp * (pzero / play)**rkappa
1289
1290        !---------------------------------------------------------------------
1291        !   Nudge soil temperature if requested
1292        !---------------------------------------------------------------------
1293
1294        IF (nudge_tsoil .AND. .NOT. lastcall) THEN
1295        ftsoil(1, isoil_nudge, :) = ftsoil(1, isoil_nudge, :)                     &
1296        - timestep / tau_soil_nudge * (ftsoil(1, isoil_nudge, :) - Tsoil_nudge)
1297                ENDIF
1298
1299                !---------------------------------------------------------------------
1300                !   Add large-scale tendencies (advection, etc) :
1301                !---------------------------------------------------------------------
1302
1303                !cc nrlmd
1304                !cc        tmpvar=teta
1305                !cc        CALL advect_vert(llm,omega,timestep,tmpvar,plev)
1306                !cc
1307                !cc        teta(1:mxcalc)=tmpvar(1:mxcalc)
1308                !cc        tmpvar(:)=q(:,1)
1309                !cc        CALL advect_vert(llm,omega,timestep,tmpvar,plev)
1310                !cc        q(1:mxcalc,1)=tmpvar(1:mxcalc)
1311                !cc        tmpvar(:)=q(:,2)
1312                !cc        CALL advect_vert(llm,omega,timestep,tmpvar,plev)
1313                !cc        q(1:mxcalc,2)=tmpvar(1:mxcalc)
1314
1315                END IF ! end if tendency of tendency should be added
1316
1317                !---------------------------------------------------------------------
1318                !   Air temperature :
1319                !---------------------------------------------------------------------
1320                IF (lastcall) THEN
1321                PRINT*, 'Pas de temps final ', it
1322                CALL ju2ymds(daytime, an, mois, jour, heure)
1323                        PRINT*, 'a la date : a m j h', an, mois, jour, heure / 3600.
1324                END IF
1325
1326                !  incremente day time
1327                        !        PRINT*,'daytime bef',daytime,1./day_step
1328                        daytime = daytime + 1. / day_step
1329                        !Al1dbg
1330                        day = int(daytime + 0.1 / day_step)
1331                !        time = max(daytime-day,0.0)
1332                !Al1&jyg: correction de bug
1333                !cc        time = real(mod(it,day_step))/day_step
1334                time = time_ini / 24. + real(mod(it, day_step)) / day_step
1335                !        PRINT*,'daytime nxt time',daytime,time
1336                it = it + 1
1337
1338                        END DO
1339
1340                        !Al1
1341                        IF (ecrit_slab_oc/=-1) close(97)
1342
1343                        !Al1 Call to 1D equivalent of dynredem (an,mois,jour,heure ?)
1344                        ! -------------------------------------
1345                        CALL dyn1dredem("restart1dyn.nc", &
1346                plev, play, phi, phis, presnivs, &
1347                u, v, temp, q, omega2)
1348
1349                CALL abort_gcm ('lmdz1d   ', 'The End  ', 0)
1350
1351                        END SUBROUTINE old_lmdz1d
1352
1353                        INCLUDE "old_1DUTILS_read_interp.h"
1354        END MODULE lmdz_old_lmdz1d
Note: See TracBrowser for help on using the repository browser.