source: LMDZ5/trunk/libf/phylmd/dyn1d/lmdz1d.F90 @ 2601

Last change on this file since 2601 was 2601, checked in by Ehouarn Millour, 8 years ago

Cleanup in the dynamics: turn temps.h into module temps_mod.F90
EM

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