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

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

rename modules properly lmdz_*
move some unused files to obsolete/
(lint) uppercase fortran keywords

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