source: LMDZ6/trunk/libf/phylmd/dyn1d/old_lmdz1d.F90 @ 3781

Last change on this file since 3781 was 3781, checked in by evignon, 4 years ago

Initialisation de la TKE pour les cas 1D (important pour GABLS1), Etienne

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