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

Last change on this file since 3679 was 3595, checked in by fhourdin, 5 years ago

Details SCM standard format

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