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

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

Improving the physics/dynamics interface:

  • added module callphysiq_mod.F90 in dynphy_lonlat/phy* which contains the routine "call_physiq" which is called by calfis* and calls the physics. This way different "physiq" routine from different physics packages may be called: The calfis* routines now exposes all available fields that might be transmitted to physiq but which is actually send (ie: expected/needed by physiq) is decided in call_physiq.
  • turned "physiq.F90" into module "physiq_mod.F90" for better control of "physiq" arguments. Extracted embeded "gr_fi_ecrit" as self-standing routine (but note that this routine actually only works in serial mode).

EM

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