source: LMDZ6/trunk/libf/phylmd/dyn1d/lmdz1d.F90 @ 3511

Last change on this file since 3511 was 3442, checked in by Laurent Fairhead, 6 years ago

Following merging of DYNAMICO/LMDZ branches, modifications needed for the 1D LMDZ model
Added a 1D special version of dimphy which is called before the real dimphy in lmdz1d
MPL

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