source: LMDZ6/trunk/libf/phylmd/dyn1d/scm.F90 @ 5213

Last change on this file since 5213 was 5213, checked in by evignon, 8 weeks ago

petits ajustements du 1d pour inclure proprement la neige soufflee

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