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

Last change on this file since 5218 was 5186, checked in by abarral, 3 months ago

Encapsulate files in modules

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