source: LMDZ5/branches/IPSLCM5A2.1_ISO/libf/phylmd/dyn1d/lmdz1d.F90 @ 5456

Last change on this file since 5456 was 2611, checked in by jyg, 8 years ago

Introduction of a new option inhibiting the
evolution of the state variables, while calling
all parametrizations. The option is controlled by
the parameter flag_inhib_tend.

For the time being the flag is set to 0 (= no

inhibition of tendencies).

jyg for jld.

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