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

Last change on this file since 2483 was 2465, checked in by fhourdin, 9 years ago

Correction pour axe des temps en 1D. Correction d'une correction.
Bug fixing for 1D time axis

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