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

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

Modification pour le cas ayotte au format standard.
+ correction pour le cas forcing_les (correspondant a forcing_type=0)
partage par ayotte et rce_oce_les.
Dans la nouvelle version, ce dernier correspond a forcing_type=-1
Attention aux problemes de retrocomptabilite.

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