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

Last change on this file since 2716 was 2716, checked in by fhourdin, 7 years ago

Inclusion du cas arm_cu2, avec les nouveaux formats de forçage 1D
(Marie-Pierre Lefebvre)

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