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

Last change on this file since 2393 was 2393, checked in by jyg, 9 years ago

Add various intializations of arrays in lmdz1d.F90
and in the convection scheme. Add output variables
for boundary layer splitting.

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