source: LMDZ6/branches/Ocean_skin/libf/phylmd/dyn1d/scm.F90 @ 3811

Last change on this file since 3811 was 3798, checked in by lguez, 4 years ago

Sync latest trunk changes to Ocean_skin

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