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

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

Missing usage of comvert_mod in lmdzd1d since the creation of the former (r2600).
EM+MPL

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