source: LMDZ5/branches/LMDZ_tree_FC/libf/phylmd/dyn1d/lmdz1d.F90 @ 5409

Last change on this file since 5409 was 2925, checked in by fcheruy, 7 years ago

Update tree branch to trunk version

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