source: LMDZ6/trunk/libf/phylmd/dyn1d/scm.f90 @ 5627

Last change on this file since 5627 was 5627, checked in by amaison, 2 months ago

Representation of heterogeneous continental subsurfaces with parameter or flux aggregation in the simplified surface model (bucket) for 1D case studies.

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