source: LMDZ6/trunk/libf/phylmd/dyn1d/old_lmdz1d.f90 @ 5452

Last change on this file since 5452 was 5368, checked in by abarral, 4 weeks ago

(WIP) Turn implicit into explicit declarations
Turn 1dconv.h into a module to that end

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