source: LMDZ6/branches/contrails/libf/phylmd/dyn1d/scm.f90 @ 5428

Last change on this file since 5428 was 5302, checked in by abarral, 2 months ago

Turn compar1d.h date_cas.h into module
Move fcg_racmo.h to obsolete

File size: 40.7 KB
Line 
1SUBROUTINE scm
2
3USE flux_arp_mod_h
4      USE compbl_mod_h
5         USE clesphys_mod_h
6      USE ioipsl, only: ju2ymds, ymds2ju, ioconf_calendar,getin
7   USE phys_state_var_mod, ONLY : phys_state_var_init, phys_state_var_end, &
8       clwcon, detr_therm, &
9       qsol, fevap, z0m, z0h, agesno, &
10       du_gwd_rando, du_gwd_front, entr_therm, f0, fm_therm, &
11       falb_dir, falb_dif, &
12       ftsol, beta_aridity, pbl_tke, pctsrf, radsol, rain_fall, snow_fall, ratqs, &
13       rnebcon, rugoro, sig1, w01, solaire_etat0, sollw, sollwdown, &
14       solsw, solswfdiff, t_ancien, q_ancien, u_ancien, v_ancien, &
15       wake_delta_pbl_TKE, delta_tsurf, wake_fip, wake_pe, &
16       wake_deltaq, wake_deltat, wake_s, awake_s, wake_dens, &
17       awake_dens, cv_gen, wake_cstar, &
18       zgam, zmax0, zmea, zpic, zsig, &
19       zstd, zthe, zval, ale_bl, ale_bl_trig, alp_bl, &
20       ql_ancien, qs_ancien, qbs_ancien, cf_ancien, rvc_ancien, &
21       prlw_ancien, prsw_ancien, prbsw_ancien, prw_ancien, &
22       u10m,v10m,ale_wake,ale_bl_stat, ratqs_inter_
23
24
25   USE dimphy
26   USE surface_data, only : type_ocean,ok_veget
27   USE pbl_surface_mod, only : ftsoil, pbl_surface_init, &
28                                 pbl_surface_final
29   USE fonte_neige_mod, only : fonte_neige_init, fonte_neige_final
30
31   USE infotrac ! new
32   USE control_mod
33   USE indice_sol_mod
34   USE phyaqua_mod
35!  USE mod_1D_cases_read
36   USE mod_1D_cases_read_std
37   !USE mod_1D_amma_read
38   USE print_control_mod, ONLY: lunout, prt_level
39   USE iniphysiq_mod, ONLY: iniphysiq
40   USE mod_const_mpi, ONLY: comm_lmdz
41   USE physiq_mod, ONLY: physiq
42   USE comvert_mod, ONLY: presnivs, ap, bp, dpres,nivsig, nivsigs, pa, &
43                          preff, aps, bps, pseudoalt, scaleheight
44   USE temps_mod, ONLY: annee_ref, calend, day_end, day_ini, day_ref, &
45                        itau_dyn, itau_phy, start_time, year_len
46   USE phys_cal_mod, ONLY : year_len_phys_cal_mod => year_len
47   USE lmdz_cppkeys_wrapper, ONLY: CPPKEY_OUTPUTPHYSSCM
48
49   USE dimensions_mod, ONLY: iim, jjm, llm, ndm
50   USE dimsoil_mod_h, ONLY: nsoilmx
51  USE yomcst_mod_h
52  USE tsoilnudge_mod_h
53  USE fcg_gcssold_mod_h
54  USE compar1d_mod_h
55  USE date_cas_mod_h
56implicit none
57
58!=====================================================================
59! DECLARATIONS
60!=====================================================================
61
62!---------------------------------------------------------------------
63!  Externals
64!---------------------------------------------------------------------
65      external fq_sat
66      real fq_sat
67
68!---------------------------------------------------------------------
69!  Arguments d' initialisations de la physique (USER DEFINE)
70!---------------------------------------------------------------------
71
72      integer, parameter :: ngrid=1
73      real :: zcufi    = 1.
74      real :: zcvfi    = 1.
75      real :: fnday
76      real :: day, daytime
77      real :: day1
78      real :: heure
79      integer :: jour
80      integer :: mois
81      integer :: an
82 
83!---------------------------------------------------------------------
84!  Declarations related to forcing and initial profiles
85!---------------------------------------------------------------------
86
87        integer :: kmax = llm
88        integer llm700,nq1,nq2
89        INTEGER, PARAMETER :: nlev_max=1000, nqmx=1000
90        real timestep, frac
91        real height(nlev_max),tttprof(nlev_max),qtprof(nlev_max)
92        real  uprof(nlev_max),vprof(nlev_max),e12prof(nlev_max)
93        real  ugprof(nlev_max),vgprof(nlev_max),wfls(nlev_max)
94        real  dqtdxls(nlev_max),dqtdyls(nlev_max)
95        real  dqtdtls(nlev_max),thlpcar(nlev_max)
96        real  qprof(nlev_max,nqmx)
97
98!        integer :: forcing_type
99        logical :: forcing_les     = .false.
100        logical :: forcing_armcu   = .false.
101        logical :: forcing_rico    = .false.
102        logical :: forcing_radconv = .false.
103        logical :: forcing_toga    = .false.
104        logical :: forcing_twpice  = .false.
105        logical :: forcing_amma    = .false.
106        logical :: forcing_dice    = .false.
107        logical :: forcing_gabls4  = .false.
108
109        logical :: forcing_GCM2SCM = .false.
110        logical :: forcing_GCSSold = .false.
111        logical :: forcing_sandu   = .false.
112        logical :: forcing_astex   = .false.
113        logical :: forcing_fire    = .false.
114        logical :: forcing_case    = .false.
115        logical :: forcing_case2   = .false.
116        logical :: forcing_SCM   = .false.
117
118!flag forcings
119        logical :: nudge_wind=.true.
120        logical :: nudge_thermo=.false.
121        logical :: cptadvw=.true.
122
123
124!=====================================================================
125! DECLARATIONS FOR EACH CASE
126!=====================================================================
127!
128      INCLUDE "1D_decl_cases.h"
129!
130!---------------------------------------------------------------------
131!  Declarations related to nudging
132!---------------------------------------------------------------------
133     integer :: nudge_max
134     parameter (nudge_max=9)
135     integer :: inudge_RHT=1
136     integer :: inudge_UV=2
137     logical :: nudge(nudge_max)
138     real :: t_targ(llm)
139     real :: rh_targ(llm)
140     real :: u_targ(llm)
141     real :: v_targ(llm)
142!
143!---------------------------------------------------------------------
144!  Declarations related to vertical discretization:
145!---------------------------------------------------------------------
146      real :: pzero=1.e5
147      real :: play (llm),zlay (llm),sig_s(llm),plev(llm+1)
148      real :: playd(llm),zlayd(llm),ap_amma(llm+1),bp_amma(llm+1)
149
150!---------------------------------------------------------------------
151!  Declarations related to variables
152!---------------------------------------------------------------------
153
154      real :: phi(llm)
155      real :: teta(llm),tetal(llm),temp(llm),u(llm),v(llm),w(llm)
156      REAL rot(1, llm) ! relative vorticity, in s-1
157      real :: rlat_rad(1),rlon_rad(1)
158      real :: omega(llm),omega2(llm),rho(llm+1)
159      real :: ug(llm),vg(llm),fcoriolis
160      real :: sfdt, cfdt
161      real :: du_phys(llm),dv_phys(llm),dt_phys(llm)
162      real :: w_adv(llm),z_adv(llm)
163      real :: d_t_vert_adv(llm),d_u_vert_adv(llm),d_v_vert_adv(llm)
164      real :: dt_cooling(llm),d_t_adv(llm),d_t_nudge(llm)
165      real :: d_u_nudge(llm),d_v_nudge(llm)
166!      real :: d_u_adv(llm),d_v_adv(llm)
167      real :: d_u_age(llm),d_v_age(llm)
168      real :: alpha
169      real :: ttt
170
171      REAL, ALLOCATABLE, DIMENSION(:,:):: q
172      REAL, ALLOCATABLE, DIMENSION(:,:):: dq
173      REAL, ALLOCATABLE, DIMENSION(:,:):: d_q_vert_adv
174      REAL, ALLOCATABLE, DIMENSION(:,:):: d_q_adv
175      REAL, ALLOCATABLE, DIMENSION(:,:):: d_q_nudge
176!      REAL, ALLOCATABLE, DIMENSION(:):: d_th_adv
177
178!---------------------------------------------------------------------
179!  Initialization of surface variables
180!---------------------------------------------------------------------
181      real :: run_off_lic_0(1)
182      real :: fder(1),snsrf(1,nbsrf),qsurfsrf(1,nbsrf)
183      real :: tsoil(1,nsoilmx,nbsrf)
184!     real :: agesno(1,nbsrf)
185
186!---------------------------------------------------------------------
187!  Call to phyredem
188!---------------------------------------------------------------------
189      logical :: ok_writedem =.true.
190      real :: sollw_in = 0.
191      real :: solsw_in = 0.
192     
193!---------------------------------------------------------------------
194!  Call to physiq
195!---------------------------------------------------------------------
196      logical :: firstcall=.true.
197      logical :: lastcall=.false.
198      real :: phis(1)    = 0.0
199      real :: dpsrf(1)
200
201!---------------------------------------------------------------------
202!  Initializations of boundary conditions
203!---------------------------------------------------------------------
204      real, allocatable :: phy_nat (:)  ! 0=ocean libre,1=land,2=glacier,3=banquise
205      real, allocatable :: phy_alb (:)  ! Albedo land only (old value condsurf_jyg=0.3)
206      real, allocatable :: phy_sst (:)  ! SST (will not be used; cf read_tsurf1d.F)
207      real, allocatable :: phy_bil (:)  ! Ne sert que pour les slab_ocean
208      real, allocatable :: phy_rug (:) ! Longueur rugosite utilisee sur land only
209      real, allocatable :: phy_ice (:) ! Fraction de glace
210      real, allocatable :: phy_fter(:) ! Fraction de terre
211      real, allocatable :: phy_foce(:) ! Fraction de ocean
212      real, allocatable :: phy_fsic(:) ! Fraction de glace
213      real, allocatable :: phy_flic(:) ! Fraction de glace
214
215!---------------------------------------------------------------------
216!  Fichiers et d'autres variables
217!---------------------------------------------------------------------
218      integer :: k,l,i,it=1,mxcalc
219      integer :: nsrf
220      integer jcode
221      INTEGER read_climoz
222!
223      integer :: it_end ! iteration number of the last call
224!Al1,plev,play,phi,phis,presnivs,
225      integer ecrit_slab_oc !1=ecrit,-1=lit,0=no file
226      data ecrit_slab_oc/-1/
227!
228!     if flag_inhib_forcing = 0, tendencies of forcing are added
229!                           <> 0, tendencies of forcing are not added
230      INTEGER :: flag_inhib_forcing = 0
231
232
233      print*,'VOUS ENTREZ DANS LE 1D FORMAT STANDARD'
234
235!=====================================================================
236! INITIALIZATIONS
237!=====================================================================
238      du_phys(:)=0.
239      dv_phys(:)=0.
240      dt_phys(:)=0.
241      d_t_vert_adv(:)=0.
242      d_u_vert_adv(:)=0.
243      d_v_vert_adv(:)=0.
244      dt_cooling(:)=0.
245      d_t_adv(:)=0.
246      d_t_nudge(:)=0.
247      d_u_nudge(:)=0.
248      d_v_nudge(:)=0.
249      d_u_adv(:)=0.
250      d_v_adv(:)=0.
251      d_u_age(:)=0.
252      d_v_age(:)=0.
253     
254     
255! Initialization of Common turb_forcing
256       dtime_frcg = 0.
257       Turb_fcg_gcssold=.false.
258       hthturb_gcssold = 0.
259       hqturb_gcssold = 0.
260
261
262
263
264!---------------------------------------------------------------------
265! OPTIONS OF THE 1D SIMULATION (lmdz1d.def => unicol.def)
266!---------------------------------------------------------------------
267        call conf_unicol
268!Al1 moves this gcssold var from common fcg_gcssold to
269        Turb_fcg_gcssold = xTurb_fcg_gcssold
270! --------------------------------------------------------------------
271        close(1)
272        write(*,*) 'lmdz1d.def lu => unicol.def'
273
274       forcing_SCM = .true.
275       year_ini_cas=1997
276       ! It is possible that those parameters are run twice.
277       ! A REVOIR : LIRE PEUT ETRE AN MOIS JOUR DIRECETEMENT
278
279
280       call getin('anneeref',year_ini_cas)
281       call getin('dayref',day_deb)
282       mth_ini_cas=1 ! pour le moment on compte depuis le debut de l'annee
283       call getin('time_ini',heure_ini_cas)
284
285      print*,'NATURE DE LA SURFACE ',nat_surf
286!
287! Initialization of the logical switch for nudging
288
289     jcode = iflag_nudge
290     do i = 1,nudge_max
291       nudge(i) = mod(jcode,10) .ge. 1
292       jcode = jcode/10
293     enddo
294!-----------------------------------------------------------------------
295!  Definition of the run
296!-----------------------------------------------------------------------
297
298      call conf_gcm( 99, .TRUE. )
299     
300!-----------------------------------------------------------------------
301      allocate( phy_nat (year_len))  ! 0=ocean libre,1=land,2=glacier,3=banquise
302      phy_nat(:)=0.0
303      allocate( phy_alb (year_len))  ! Albedo land only (old value condsurf_jyg=0.3)
304      allocate( phy_sst (year_len))  ! SST (will not be used; cf read_tsurf1d.F)
305      allocate( phy_bil (year_len))  ! Ne sert que pour les slab_ocean
306      phy_bil(:)=1.0
307      allocate( phy_rug (year_len)) ! Longueur rugosite utilisee sur land only
308      allocate( phy_ice (year_len)) ! Fraction de glace
309      phy_ice(:)=0.0
310      allocate( phy_fter(year_len)) ! Fraction de terre
311      phy_fter(:)=0.0
312      allocate( phy_foce(year_len)) ! Fraction de ocean
313      phy_foce(:)=0.0
314      allocate( phy_fsic(year_len)) ! Fraction de glace
315      phy_fsic(:)=0.0
316      allocate( phy_flic(year_len)) ! Fraction de glace
317      phy_flic(:)=0.0
318
319
320!-----------------------------------------------------------------------
321!   Choix du calendrier
322!   -------------------
323
324!      calend = 'earth_365d'
325      if (calend == 'earth_360d') then
326        call ioconf_calendar('360_day')
327        write(*,*)'CALENDRIER CHOISI: Terrestre a 360 jours/an'
328      else if (calend == 'earth_365d') then
329        call ioconf_calendar('noleap')
330        write(*,*)'CALENDRIER CHOISI: Terrestre a 365 jours/an'
331      else if (calend == 'earth_366d') then
332        call ioconf_calendar('all_leap')
333        write(*,*)'CALENDRIER CHOISI: Terrestre bissextile'
334      else if (calend == 'gregorian') then
335        stop 'gregorian calend should not be used by normal user'
336        call ioconf_calendar('gregorian') ! not to be used by normal users
337        write(*,*)'CALENDRIER CHOISI: Gregorien'
338      else
339        write (*,*) 'ERROR : unknown calendar ', calend
340        stop 'calend should be 360d,earth_365d,earth_366d,gregorian'
341      endif
342!-----------------------------------------------------------------------
343!
344!c Date :
345!      La date est supposee donnee sous la forme [annee, numero du jour dans
346!      l annee] ; l heure est donnee dans time_ini, lu dans lmdz1d.def.
347!      On appelle ymds2ju pour convertir [annee, jour] en [jour Julien].
348!      Le numero du jour est dans "day". L heure est traitee separement.
349!      La date complete est dans "daytime" (l'unite est le jour).
350
351
352      if (nday>0) then
353         fnday=nday
354      else
355         fnday=-nday/float(day_step)
356      endif
357      print *,'fnday=',fnday
358!     start_time doit etre en FRACTION DE JOUR
359      start_time=time_ini/24.
360
361      annee_ref = anneeref
362      mois = 1
363      day_ref = dayref
364      heure = 0.
365      itau_dyn = 0
366      itau_phy = 0
367      call ymds2ju(annee_ref,mois,day_ref,heure,day)
368      day_ini = int(day)
369      day_end = day_ini + int(fnday)
370
371! Convert the initial date to Julian day
372      day_ini_cas=day_deb
373      print*,'time case',year_ini_cas,mth_ini_cas,day_ini_cas
374      call ymds2ju                                                         &
375     & (year_ini_cas,mth_ini_cas,day_ini_cas,heure_ini_cas*3600            &
376     & ,day_ju_ini_cas)
377      print*,'time case 2',day_ini_cas,day_ju_ini_cas
378      daytime = day + heure_ini_cas/24. ! 1st day and initial time of the simulation
379
380! Print out the actual date of the beginning of the simulation :
381      call ju2ymds(daytime,year_print, month_print,day_print,sec_print)
382      print *,' Time of beginning : ',                                      &
383     &        year_print, month_print, day_print, sec_print
384
385!---------------------------------------------------------------------
386! Initialization of dimensions, geometry and initial state
387!---------------------------------------------------------------------
388!     call init_phys_lmdz(1,1,llm,1,(/1/)) ! job now done via iniphysiq
389!     but we still need to initialize dimphy module (klon,klev,etc.)  here.
390      call init_dimphy1D(1,llm)
391      call suphel
392      call init_infotrac
393
394      if (nqtot>nqmx) STOP 'Augmenter nqmx dans lmdz1d.F'
395      allocate(q(llm,nqtot)) ; q(:,:)=0.
396      allocate(dq(llm,nqtot))
397      allocate(d_q_vert_adv(llm,nqtot))
398      allocate(d_q_adv(llm,nqtot))
399      allocate(d_q_nudge(llm,nqtot))
400!      allocate(d_th_adv(llm))
401
402      q(:,:) = 0.
403      dq(:,:) = 0.
404      d_q_vert_adv(:,:) = 0.
405      d_q_adv(:,:) = 0.
406      d_q_nudge(:,:) = 0.
407
408!
409!   No ozone climatology need be read in this pre-initialization
410!          (phys_state_var_init is called again in physiq)
411      read_climoz = 0
412      nsw=6
413
414      call phys_state_var_init(read_climoz)
415
416      if (ngrid.ne.klon) then
417         print*,'stop in inifis'
418         print*,'Probleme de dimensions :'
419         print*,'ngrid = ',ngrid
420         print*,'klon  = ',klon
421         stop
422      endif
423!!!=====================================================================
424!!! Feedback forcing values for Gateaux differentiation (al1)
425!!!=====================================================================
426!!
427      qsol = qsolinp
428      qsurf = fq_sat(tsurf,psurf/100.)
429      beta_aridity(:,:) = beta_surf
430      day1= day_ini
431      time=daytime-day
432      ts_toga(1)=tsurf ! needed by read_tsurf1d.F
433      rho(1)=psurf/(rd*tsurf*(1.+(rv/rd-1.)*qsurf))
434
435!
436!! mpl et jyg le 22/08/2012 :
437!!  pour que les cas a flux de surface imposes marchent
438      IF(.NOT.ok_flux_surf.or.max(abs(wtsurf),abs(wqsurf))>0.) THEN
439       fsens=-wtsurf*rcpd*rho(1)
440       flat=-wqsurf*rlvtt*rho(1)
441       print *,'Flux: ok_flux wtsurf wqsurf',ok_flux_surf,wtsurf,wqsurf
442      ENDIF
443      print*,'Flux sol ',fsens,flat
444
445! Vertical discretization and pressure levels at half and mid levels:
446
447      pa   = 5e4
448!!      preff= 1.01325e5
449      preff = psurf
450      IF (ok_old_disvert) THEN
451        call disvert0(pa,preff,ap,bp,dpres,presnivs,nivsigs,nivsig)
452        print *,'On utilise disvert0'
453        aps(1:llm)=0.5*(ap(1:llm)+ap(2:llm+1))
454        bps(1:llm)=0.5*(bp(1:llm)+bp(2:llm+1))
455        scaleheight=8.
456        pseudoalt(1:llm)=-scaleheight*log(presnivs(1:llm)/preff)
457      ELSE
458        call disvert()
459        print *,'On utilise disvert'
460!       Nouvelle version disvert permettant d imposer ap,bp (modif L.Guez) MPL 18092012
461!       Dans ce cas, on lit ap,bp dans le fichier hybrid.txt
462      ENDIF
463
464      sig_s=presnivs/preff
465      plev =ap+bp*psurf
466      play = 0.5*(plev(1:llm)+plev(2:llm+1))
467      zlay=-rd*300.*log(play/psurf)/rg ! moved after reading profiles.
468
469      IF (forcing_type .eq. 59) THEN
470! pour forcing_sandu, on cherche l'indice le plus proche de 700hpa#3000m
471      write(*,*) '***********************'
472      do l = 1, llm
473       write(*,*) 'l,play(l),presnivs(l): ',l,play(l),presnivs(l)
474       if (trouve_700 .and. play(l).le.70000) then
475         llm700=l
476         print *,'llm700,play=',llm700,play(l)/100.
477         trouve_700= .false.
478       endif
479      enddo
480      write(*,*) '***********************'
481      ENDIF
482
483!
484!=====================================================================
485! EVENTUALLY, READ FORCING DATA :
486!=====================================================================
487
488      INCLUDE "1D_read_forc_cases.h"
489   print*,'A d_t_adv ',d_t_adv(1:20)*86400
490
491      if (forcing_GCM2SCM) then
492        write (*,*) 'forcing_GCM2SCM not yet implemented'
493        stop 'in initialization'
494      endif ! forcing_GCM2SCM
495
496
497!=====================================================================
498! Initialisation de la physique :
499!=====================================================================
500
501!  Rq: conf_phys.F90 lit tous les flags de physiq.def; conf_phys appele depuis physiq.F
502!
503! day_step, iphysiq lus dans gcm.def ci-dessus
504! timestep: calcule ci-dessous from rday et day_step
505! ngrid=1
506! llm: defini dans .../modipsl/modeles/LMDZ4/libf/grid/dimension
507! rday: defini dans suphel.F (86400.)
508! day_ini: lu dans run.def (dayref)
509! rlat_rad,rlon-rad: transformes en radian de rlat,rlon lus dans lmdz1d.def (en degres)
510! airefi,zcufi,zcvfi initialises au debut de ce programme
511! rday,ra,rg,rd,rcpd declares dans YOMCST.h et calcules dans suphel.F
512
513
514      day_step = float(nsplit_phys)*day_step/float(iphysiq)
515      write (*,*) 'Time step divided by nsplit_phys (=',nsplit_phys,')'
516      timestep =rday/day_step
517      dtime_frcg = timestep
518!
519      zcufi=airefi
520      zcvfi=airefi
521!
522      rlat_rad(1)=xlat*rpi/180.
523      rlon_rad(1)=xlon*rpi/180.
524
525     ! iniphysiq will call iniaqua who needs year_len from phys_cal_mod
526     year_len_phys_cal_mod=year_len
527           
528     ! Ehouarn: iniphysiq requires arrays related to (3D) dynamics grid,
529     ! e.g. for cell boundaries, which are meaningless in 1D; so pad these
530     ! with '0.' when necessary
531
532      call iniphysiq(iim,jjm,llm, &
533           1,comm_lmdz, &
534           rday,day_ini,timestep,  &
535           (/rlat_rad(1),0./),(/0./), &
536           (/0.,0./),(/rlon_rad(1),0./),  &
537           (/ (/airefi,0./),(/0.,0./) /), &
538           (/zcufi,0.,0.,0./), &
539           (/zcvfi,0./), &
540           ra,rg,rd,rcpd,1)
541      print*,'apres iniphysiq'
542
543! 2 PARAMETRES QUI DEVRAIENT ETRE LUS DANS run.def MAIS NE LE SONT PAS ICI:
544      co2_ppm= 330.0
545      solaire=1370.0
546
547! Ecriture du startphy avant le premier appel a la physique.
548! On le met juste avant pour avoir acces a tous les champs
549
550      if (ok_writedem) then
551
552!--------------------------------------------------------------------------
553! pbl_surface_init (called here) and pbl_surface_final (called by phyredem)
554! need : qsol fder snow qsurf evap rugos agesno ftsoil
555!--------------------------------------------------------------------------
556
557        type_ocean = "force"
558        run_off_lic_0(1) = restart_runoff
559        call fonte_neige_init(run_off_lic_0)
560
561        fder=0.
562        snsrf(1,:)=snowmass ! masse de neige des sous surface
563        qsurfsrf(1,:)=qsurf ! humidite de l'air des sous surface
564        fevap=0.
565        z0m(1,:)=rugos     ! couverture de neige des sous surface
566        z0h(1,:)=rugosh    ! couverture de neige des sous surface
567        agesno  = xagesno
568        tsoil(:,:,:)=tsurf
569!------ AMMA 2e run avec modele sol et rayonnement actif (MPL 23052012)
570!       tsoil(1,1,1)=299.18
571!       tsoil(1,2,1)=300.08
572!       tsoil(1,3,1)=301.88
573!       tsoil(1,4,1)=305.48
574!       tsoil(1,5,1)=308.00
575!       tsoil(1,6,1)=308.00
576!       tsoil(1,7,1)=308.00
577!       tsoil(1,8,1)=308.00
578!       tsoil(1,9,1)=308.00
579!       tsoil(1,10,1)=308.00
580!       tsoil(1,11,1)=308.00
581!-----------------------------------------------------------------------
582        call pbl_surface_init(fder, snsrf, qsurfsrf, tsoil)
583
584!------------------ prepare limit conditions for limit.nc -----------------
585!--   Ocean force
586
587        print*,'avant phyredem'
588        pctsrf(1,:)=0.
589          if (nat_surf.eq.0.) then
590          pctsrf(1,is_oce)=1.
591          pctsrf(1,is_ter)=0.
592          pctsrf(1,is_lic)=0.
593          pctsrf(1,is_sic)=0.
594        else if (nat_surf .eq. 1) then
595          pctsrf(1,is_oce)=0.
596          pctsrf(1,is_ter)=1.
597          pctsrf(1,is_lic)=0.
598          pctsrf(1,is_sic)=0.
599        else if (nat_surf .eq. 2) then
600          pctsrf(1,is_oce)=0.
601          pctsrf(1,is_ter)=0.
602          pctsrf(1,is_lic)=1.
603          pctsrf(1,is_sic)=0.
604        else if (nat_surf .eq. 3) then
605          pctsrf(1,is_oce)=0.
606          pctsrf(1,is_ter)=0.
607          pctsrf(1,is_lic)=0.
608          pctsrf(1,is_sic)=1.
609
610     end if
611
612
613        print*,'nat_surf,pctsrf(1,is_oce),pctsrf(1,is_ter)',nat_surf         &
614     &        ,pctsrf(1,is_oce),pctsrf(1,is_ter)
615
616        zmasq=pctsrf(1,is_ter)+pctsrf(1,is_lic)
617        zpic = zpicinp
618        ftsol=tsurf
619        falb_dir=albedo
620        falb_dif=albedo
621        rugoro=rugos
622        t_ancien(1,:)=temp(:)
623        q_ancien(1,:)=q(:,1)
624        ql_ancien = 0.
625        qs_ancien = 0.
626        prlw_ancien = 0.
627        prsw_ancien = 0.
628        prw_ancien = 0.
629        IF ( ok_bs ) THEN
630           qbs_ancien = 0.
631           prbsw_ancien = 0.
632        ENDIF
633        IF ( ok_ice_supersat ) THEN
634          cf_ancien = 0.
635          rvc_ancien = 0.
636        ENDIF
637!jyg<
638! Etienne: comment those lines since now the TKE is inialized in 1D_read_forc_cases
639!!      pbl_tke(:,:,:)=1.e-8
640!        pbl_tke(:,:,:)=0.
641!        pbl_tke(:,2,:)=1.e-2
642!>jyg
643        rain_fall=0.
644        snow_fall=0.
645        solsw=0.
646        solswfdiff=0.
647        sollw=0.
648        sollwdown=rsigma*tsurf**4
649        radsol=0.
650        rnebcon=0.
651        ratqs=0.
652        clwcon=0.
653        zmax0 = 0.
654        zmea=zsurf
655        zstd=0.
656        zsig=0.
657        zgam=0.
658        zval=0.
659        zthe=0.
660        sig1=0.
661        w01=0.
662!
663        wake_deltaq = 0.
664        wake_deltat = 0.
665        wake_delta_pbl_TKE(:,:,:) = 0.
666        delta_tsurf = 0.
667        wake_fip = 0.
668        wake_pe = 0.
669        wake_s = 0.
670        awake_s = 0.
671        wake_dens = 0.
672        awake_dens = 0.
673        cv_gen = 0.
674        wake_cstar = 0.
675        ale_bl = 0.
676        ale_bl_trig = 0.
677        alp_bl = 0.
678        IF (ALLOCATED(du_gwd_rando)) du_gwd_rando = 0.
679        IF (ALLOCATED(du_gwd_front)) du_gwd_front = 0.
680        entr_therm = 0.
681        detr_therm = 0.
682        f0 = 0.
683        fm_therm = 0.
684        u_ancien(1,:)=u(:)
685        v_ancien(1,:)=v(:)
686 
687        u10m=0.
688        v10m=0.
689        ale_wake=0.
690        ale_bl_stat=0.
691        ratqs_inter_(:,:)= 0.002
692
693!------------------------------------------------------------------------
694! Make file containing restart for the physics (startphy.nc)
695!
696! NB: List of the variables to be written by phyredem (via put_field):
697! rlon,rlat,zmasq,pctsrf(:,is_ter),pctsrf(:,is_lic),pctsrf(:,is_oce)
698! pctsrf(:,is_sic),ftsol(:,nsrf),tsoil(:,isoil,nsrf),qsurf(:,nsrf)
699! qsol,falb_dir(:,nsrf),falb_dif(:,nsrf),evap(:,nsrf),snow(:,nsrf)
700! radsol,solsw,solswfdiff,sollw, sollwdown,fder,rain_fall,snow_fall,frugs(:,nsrf)
701! agesno(:,nsrf),zmea,zstd,zsig,zgam,zthe,zpic,zval,rugoro
702! t_ancien,q_ancien,,frugs(:,is_oce),clwcon(:,1),rnebcon(:,1),ratqs(:,1)
703! run_off_lic_0,pbl_tke(:,1:klev,nsrf), zmax0,f0,sig1,w01
704! wake_deltat,wake_deltaq,wake_s,awake_s,wake_dens,awake_dens,cv_gen,wake_cstar,
705! wake_fip,wake_delta_pbl_tke(:,1:klev,nsrf)
706!
707! NB2: The content of the startphy.nc file depends on some flags defined in
708! the ".def" files. However, since conf_phys is not called in lmdz1d.F90, these flags have
709! to be set at some arbitratry convenient values.
710!------------------------------------------------------------------------
711!Al1 =============== restart option ======================================
712        iflag_physiq=0
713        call getin('iflag_physiq',iflag_physiq)
714
715        if (.not.restart) then
716          iflag_pbl = 5
717          call phyredem ("startphy.nc")
718        else
719! (desallocations)
720        print*,'callin surf final'
721          call pbl_surface_final( fder, snsrf, qsurfsrf, tsoil)
722        print*,'after surf final'
723          CALL fonte_neige_final(run_off_lic_0)
724        endif
725
726        ok_writedem=.false.
727        print*,'apres phyredem'
728
729      endif ! ok_writedem
730     
731!------------------------------------------------------------------------
732! Make file containing boundary conditions (limit.nc) **Al1->restartdyn***
733! --------------------------------------------------
734! NB: List of the variables to be written in limit.nc
735!     (by writelim.F, subroutine of 1DUTILS.h):
736!        phy_nat,phy_alb,phy_sst,phy_bil,phy_rug,phy_ice,
737!        phy_fter,phy_foce,phy_flic,phy_fsic)
738!------------------------------------------------------------------------
739      do i=1,year_len
740        phy_nat(i)  = nat_surf
741        phy_alb(i)  = albedo
742        phy_sst(i)  = tsurf ! read_tsurf1d will be used instead
743        phy_rug(i)  = rugos
744        phy_fter(i) = pctsrf(1,is_ter)
745        phy_foce(i) = pctsrf(1,is_oce)
746        phy_fsic(i) = pctsrf(1,is_sic)
747        phy_flic(i) = pctsrf(1,is_lic)
748      enddo
749
750! fabrication de limit.nc
751      call writelim (1,phy_nat,phy_alb,phy_sst,phy_bil,phy_rug,             &
752     &               phy_ice,phy_fter,phy_foce,phy_flic,phy_fsic)
753
754
755      call phys_state_var_end
756!Al1
757      if (restart) then
758        print*,'call to restart dyn 1d'
759        Call dyn1deta0("start1dyn.nc",plev,play,phi,phis,presnivs,          &
760     &              u,v,temp,q,omega2)
761
762       print*,'fnday,annee_ref,day_ref,day_ini',                            &
763     &     fnday,annee_ref,day_ref,day_ini
764!**      call ymds2ju(annee_ref,mois,day_ini,heure,day)
765       day = day_ini
766       day_end = day_ini + nday
767       daytime = day + time_ini/24. ! 1st day and initial time of the simulation
768
769! Print out the actual date of the beginning of the simulation :
770       call ju2ymds(daytime, an, mois, jour, heure)
771       print *,' Time of beginning : y m d h',an, mois,jour,heure/3600.
772
773       day = int(daytime)
774       time=daytime-day
775 
776       print*,'****** intialised fields from restart1dyn *******'
777       print*,'plev,play,phi,phis,presnivs,u,v,temp,q,omega2'
778       print*,'temp(1),q(1,1),u(1),v(1),plev(1),phis :'
779       print*,temp(1),q(1,1),u(1),v(1),plev(1),phis(1)
780! raz for safety
781       do l=1,llm
782         d_q_vert_adv(l,1) = 0.
783       enddo
784      endif
785!======================  end restart =================================
786      IF (ecrit_slab_oc.eq.1) then
787         open(97,file='div_slab.dat',STATUS='UNKNOWN')
788       elseif (ecrit_slab_oc.eq.0) then
789         open(97,file='div_slab.dat',STATUS='OLD')
790       endif
791!
792!=====================================================================
793IF (CPPKEY_OUTPUTPHYSSCM) THEN
794       CALL iophys_ini(timestep)
795END IF
796
797!=====================================================================
798! START OF THE TEMPORAL LOOP :
799!=====================================================================
800           
801      it_end = nint(fnday*day_step)
802      do while(it.le.it_end)
803
804       if (prt_level.ge.1) then
805         print*,'XXXXXXXXXXXXXXXXXXX ITAP,day,time=',                       &
806     &             it,day,time,it_end,day_step
807         print*,'PAS DE TEMPS ',timestep
808       endif
809       if (it.eq.it_end) lastcall=.True.
810
811!---------------------------------------------------------------------
812! Interpolation of forcings in time and onto model levels
813!---------------------------------------------------------------------
814
815      INCLUDE "1D_interp_cases.h"
816
817!---------------------------------------------------------------------
818!  Geopotential :
819!---------------------------------------------------------------------
820        phis(1)=zsurf*RG
821!        phi(1)=phis(1)+RD*temp(1)*(plev(1)-play(1))/(.5*(plev(1)+play(1)))
822
823        ! Calculate geopotential from the ground surface since phi and phis are added in physiq_mod
824        phi(1)=RD*temp(1)*(plev(1)-play(1))/(.5*(plev(1)+play(1)))
825
826        do l = 1, llm-1
827          phi(l+1)=phi(l)+RD*(temp(l)+temp(l+1))*                           &
828     &    (play(l)-play(l+1))/(play(l)+play(l+1))
829        enddo
830
831!---------------------------------------------------------------------
832!  Vertical advection
833!---------------------------------------------------------------------
834
835   IF ( forc_w+forc_omega > 0 ) THEN
836
837      IF ( forc_w == 1 ) THEN
838         w_adv=w_mod_cas
839         z_adv=phi/RG
840      ELSE
841         w_adv=omega
842         z_adv=play
843      ENDIF
844
845      teta=temp*(pzero/play)**rkappa
846      do l=2,llm-1
847        ! vertical tendencies computed as d X / d t = -W  d X / d z
848        d_u_vert_adv(l)=-w_adv(l)*(u(l+1)-u(l-1))/(z_adv(l+1)-z_adv(l-1))
849        d_v_vert_adv(l)=-w_adv(l)*(v(l+1)-v(l-1))/(z_adv(l+1)-z_adv(l-1))
850        ! d theta / dt = -W d theta / d z, transformed into d temp / d t dividing by (pzero/play(l))**rkappa
851        d_t_vert_adv(l)=-w_adv(l)*(teta(l+1)-teta(l-1))/(z_adv(l+1)-z_adv(l-1)) / (pzero/play(l))**rkappa
852        d_q_vert_adv(l,1)=-w_adv(l)*(q(l+1,1)-q(l-1,1))/(z_adv(l+1)-z_adv(l-1))
853      enddo
854      d_u_adv(:)=d_u_adv(:)+d_u_vert_adv(:)
855      d_v_adv(:)=d_v_adv(:)+d_v_vert_adv(:)
856      d_t_adv(:)=d_t_adv(:)+d_t_vert_adv(:)
857      d_q_adv(:,1)=d_q_adv(:,1)+d_q_vert_adv(:,1)
858   
859   ENDIF
860
861!---------------------------------------------------------------------
862! Listing output for debug prt_level>=1
863!---------------------------------------------------------------------
864       if (prt_level>=1) then
865         print *,' avant physiq : -------- day time ',day,time
866         write(*,*) 'firstcall,lastcall,phis',                               &
867     &               firstcall,lastcall,phis
868       end if
869       if (prt_level>=5) then
870         write(*,'(a10,2a4,4a13)') 'BEFOR1 IT=','it','l',                   &
871     &        'presniv','plev','play','phi'
872         write(*,'(a10,2i4,4f13.2)') ('BEFOR1 IT= ',it,l,                   &
873     &         presnivs(l),plev(l),play(l),phi(l),l=1,llm)
874         write(*,'(a11,2a4,a11,6a8)') 'BEFOR2','it','l',                    &
875     &         'presniv','u','v','temp','q1','q2','omega2'
876         write(*,'(a11,2i4,f11.2,5f8.2,e10.2)') ('BEFOR2 IT= ',it,l,         &
877     &   presnivs(l),u(l),v(l),temp(l),q(l,1),q(l,2),omega2(l),l=1,llm)
878       endif
879
880!---------------------------------------------------------------------
881!   Call physiq :
882!---------------------------------------------------------------------
883       call physiq(ngrid,llm, &
884                    firstcall,lastcall,timestep, &
885                    plev,play,phi,phis,presnivs, &
886                    u,v, rot, temp,q,omega2, &
887                    du_phys,dv_phys,dt_phys,dq,dpsrf)
888                firstcall=.false.
889
890!---------------------------------------------------------------------
891! Listing output for debug
892!---------------------------------------------------------------------
893        if (prt_level>=5) then
894          write(*,'(a11,2a4,4a13)') 'AFTER1 IT=','it','l',                  &
895     &        'presniv','plev','play','phi'
896          write(*,'(a11,2i4,4f13.2)') ('AFTER1 it= ',it,l,                  &
897     &    presnivs(l),plev(l),play(l),phi(l),l=1,llm)
898          write(*,'(a11,2a4,a11,6a8)') 'AFTER2','it','l',                   &
899     &         'presniv','u','v','temp','q1','q2','omega2'
900          write(*,'(a11,2i4,f11.2,5f8.2,e10.2)') ('AFTER2 it= ',it,l,       &
901     &    presnivs(l),u(l),v(l),temp(l),q(l,1),q(l,2),omega2(l),l=1,llm)
902          write(*,'(a11,2a4,a11,5a8)') 'AFTER3','it','l',                   &
903     &         'presniv','du_phys','dv_phys','dt_phys','dq1','dq2'   
904           write(*,'(a11,2i4,f11.2,5f8.2)') ('AFTER3 it= ',it,l,            &
905     &      presnivs(l),86400*du_phys(l),86400*dv_phys(l),                   &
906     &       86400*dt_phys(l),86400*dq(l,1),dq(l,2),l=1,llm)
907          write(*,*) 'dpsrf',dpsrf
908        endif
909!---------------------------------------------------------------------
910!   Add physical tendencies :
911!---------------------------------------------------------------------
912
913       fcoriolis=2.*sin(rpi*xlat/180.)*romega
914
915      IF (prt_level >= 5) print*, 'fcoriolis, xlat,mxcalc ', &
916                                   fcoriolis, xlat,mxcalc
917
918!---------------------------------------------------------------------
919! Geostrophic forcing
920!---------------------------------------------------------------------
921
922      IF ( forc_geo == 0 ) THEN
923              d_u_age(1:mxcalc)=0.
924              d_v_age(1:mxcalc)=0.
925      ELSE
926       sfdt = sin(0.5*fcoriolis*timestep)
927       cfdt = cos(0.5*fcoriolis*timestep)
928
929        d_u_age(1:mxcalc)= -2.*sfdt/timestep*                                &
930     &          (sfdt*(u(1:mxcalc)-ug(1:mxcalc)) -                          &
931     &           cfdt*(v(1:mxcalc)-vg(1:mxcalc))  )
932!!     : fcoriolis*(v(1:mxcalc)-vg(1:mxcalc))
933!
934       d_v_age(1:mxcalc)= -2.*sfdt/timestep*                                 &
935     &          (cfdt*(u(1:mxcalc)-ug(1:mxcalc)) +                           &
936     &           sfdt*(v(1:mxcalc)-vg(1:mxcalc))  )
937!!     : -fcoriolis*(u(1:mxcalc)-ug(1:mxcalc))
938      ENDIF
939!
940!---------------------------------------------------------------------
941!  Nudging
942!  EV: rewrite the section to account for a time- and height-varying
943!  nudging
944!---------------------------------------------------------------------
945      d_t_nudge(:) = 0.
946      d_u_nudge(:) = 0.
947      d_v_nudge(:) = 0.
948      d_q_nudge(:,:) = 0.
949
950      DO l=1,llm
951
952         IF (nudging_u .LT. 0) THEN
953             
954             d_u_nudge(l)=(u_nudg_mod_cas(l)-u(l))*invtau_u_nudg_mod_cas(l)
955       
956         ELSE
957
958             IF ( play(l) < p_nudging_u .AND. nint(nudging_u) /= 0 ) &
959             & d_u_nudge(l)=(u_nudg_mod_cas(l)-u(l))/nudging_u
960
961         ENDIF
962
963
964         IF (nudging_v .LT. 0) THEN
965             
966             d_v_nudge(l)=(v_nudg_mod_cas(l)-v(l))*invtau_v_nudg_mod_cas(l)
967       
968         ELSE
969
970
971             IF ( play(l) < p_nudging_v .AND. nint(nudging_v) /= 0 ) &
972             & d_v_nudge(l)=(v_nudg_mod_cas(l)-v(l))/nudging_v
973
974         ENDIF
975
976
977         IF (nudging_t .LT. 0) THEN
978             
979             d_t_nudge(l)=(temp_nudg_mod_cas(l)-temp(l))*invtau_temp_nudg_mod_cas(l)
980       
981         ELSE
982
983
984             IF ( play(l) < p_nudging_t .AND. nint(nudging_t) /= 0 ) &
985             & d_t_nudge(l)=(temp_nudg_mod_cas(l)-temp(l))/nudging_t
986
987          ENDIF
988
989
990         IF (nudging_qv .LT. 0) THEN
991             
992             d_q_nudge(l,1)=(qv_nudg_mod_cas(l)-q(l,1))*invtau_qv_nudg_mod_cas(l)
993       
994         ELSE
995
996             IF ( play(l) < p_nudging_qv .AND. nint(nudging_qv) /= 0 ) &
997             & d_q_nudge(l,1)=(qv_nudg_mod_cas(l)-q(l,1))/nudging_qv
998
999         ENDIF
1000
1001      ENDDO
1002
1003!---------------------------------------------------------------------
1004!  Optional outputs
1005!---------------------------------------------------------------------
1006
1007IF (CPPKEY_OUTPUTPHYSSCM) THEN
1008      CALL iophys_ecrit('w_adv',klev,'w_adv','K/day',w_adv)
1009      CALL iophys_ecrit('z_adv',klev,'z_adv','K/day',z_adv)
1010      CALL iophys_ecrit('dtadv',klev,'dtadv','K/day',86400*d_t_adv)
1011      CALL iophys_ecrit('dtdyn',klev,'dtdyn','K/day',86400*d_t_vert_adv)
1012      CALL iophys_ecrit('qv',klev,'qv','g/kg',1000*q(:,1))
1013      CALL iophys_ecrit('qvnud',klev,'qvnud','g/kg',1000*u_nudg_mod_cas)
1014      CALL iophys_ecrit('u',klev,'u','m/s',u)
1015      CALL iophys_ecrit('unud',klev,'unud','m/s',u_nudg_mod_cas)
1016      CALL iophys_ecrit('v',klev,'v','m/s',v)
1017      CALL iophys_ecrit('vnud',klev,'vnud','m/s',v_nudg_mod_cas)
1018      CALL iophys_ecrit('temp',klev,'temp','K',temp)
1019      CALL iophys_ecrit('tempnud',klev,'temp_nudg_mod_cas','K',temp_nudg_mod_cas)
1020      CALL iophys_ecrit('dtnud',klev,'dtnud','K/day',86400*d_t_nudge)
1021      CALL iophys_ecrit('dqnud',klev,'dqnud','K/day',1000*86400*d_q_nudge(:,1))
1022END IF
1023
1024    IF (flag_inhib_forcing == 0) then ! if tendency of forcings should be added
1025
1026        u(1:mxcalc)=u(1:mxcalc) + timestep*(                                &
1027     &              du_phys(1:mxcalc)                                       &
1028     &             +d_u_age(1:mxcalc)+d_u_adv(1:mxcalc)                       &
1029     &             +d_u_nudge(1:mxcalc) )           
1030        v(1:mxcalc)=v(1:mxcalc) + timestep*(                                 &
1031     &              dv_phys(1:mxcalc)                                       &
1032     &             +d_v_age(1:mxcalc)+d_v_adv(1:mxcalc)                       &
1033     &             +d_v_nudge(1:mxcalc) )
1034        q(1:mxcalc,:)=q(1:mxcalc,:)+timestep*(                              &
1035     &                dq(1:mxcalc,:)                                        &
1036     &               +d_q_adv(1:mxcalc,:)                                   &
1037     &               +d_q_nudge(1:mxcalc,:) )
1038
1039        if (prt_level.ge.3) then
1040          print *,                                                          &
1041     &    'physiq-> temp(1),dt_phys(1),d_t_adv(1),dt_cooling(1) ',         &
1042     &              temp(1),dt_phys(1),d_t_adv(1),dt_cooling(1)
1043           print* ,'dv_phys=',dv_phys
1044           print* ,'d_v_age=',d_v_age
1045           print* ,'d_v_adv=',d_v_adv
1046           print* ,'d_v_nudge=',d_v_nudge
1047           print*, v
1048           print*, vg
1049        endif
1050
1051        temp(1:mxcalc)=temp(1:mxcalc)+timestep*(                            &
1052     &              dt_phys(1:mxcalc)                                       &
1053     &             +d_t_adv(1:mxcalc)                                       &
1054     &             +d_t_nudge(1:mxcalc)                                     &
1055     &             +dt_cooling(1:mxcalc))  ! Taux de chauffage ou refroid.
1056
1057
1058!=======================================================================
1059!! CONSERVE EN ATTENDANT QUE LE CAS EN QUESTION FONCTIONNE EN STD !!
1060!=======================================================================
1061
1062        teta=temp*(pzero/play)**rkappa
1063
1064!---------------------------------------------------------------------
1065!   Nudge soil temperature if requested
1066!---------------------------------------------------------------------
1067
1068      IF (nudge_tsoil .AND. .NOT. lastcall) THEN
1069       ftsoil(1,isoil_nudge,:) = ftsoil(1,isoil_nudge,:)                     &
1070     &  -timestep/tau_soil_nudge*(ftsoil(1,isoil_nudge,:)-Tsoil_nudge)
1071      ENDIF
1072
1073!---------------------------------------------------------------------
1074!   Add large-scale tendencies (advection, etc) :
1075!---------------------------------------------------------------------
1076
1077!cc nrlmd
1078!cc        tmpvar=teta
1079!cc        call advect_vert(llm,omega,timestep,tmpvar,plev)
1080!cc
1081!cc        teta(1:mxcalc)=tmpvar(1:mxcalc)
1082!cc        tmpvar(:)=q(:,1)
1083!cc        call advect_vert(llm,omega,timestep,tmpvar,plev)
1084!cc        q(1:mxcalc,1)=tmpvar(1:mxcalc)
1085!cc        tmpvar(:)=q(:,2)
1086!cc        call advect_vert(llm,omega,timestep,tmpvar,plev)
1087!cc        q(1:mxcalc,2)=tmpvar(1:mxcalc)
1088
1089   END IF ! end if tendency of tendency should be added
1090
1091!---------------------------------------------------------------------
1092!   Air temperature :
1093!---------------------------------------------------------------------       
1094        if (lastcall) then
1095          print*,'Pas de temps final ',it
1096          call ju2ymds(daytime, an, mois, jour, heure)
1097          print*,'a la date : a m j h',an, mois, jour ,heure/3600.
1098        endif
1099
1100!  incremente day time
1101        daytime = daytime+1./day_step
1102        day = int(daytime+0.1/day_step)
1103!        time = max(daytime-day,0.0)
1104!Al1&jyg: correction de bug
1105!cc        time = real(mod(it,day_step))/day_step
1106        time = time_ini/24.+real(mod(it,day_step))/day_step
1107        it=it+1
1108
1109      enddo
1110
1111      if (ecrit_slab_oc.ne.-1) close(97)
1112
1113!Al1 Call to 1D equivalent of dynredem (an,mois,jour,heure ?)
1114! ---------------------------------------------------------------------------
1115       call dyn1dredem("restart1dyn.nc",                                    &
1116     &              plev,play,phi,phis,presnivs,                            &
1117     &              u,v,temp,q,omega2)
1118
1119        CALL abort_gcm ('lmdz1d   ','The End  ',0)
1120
1121END SUBROUTINE scm
Note: See TracBrowser for help on using the repository browser.