source: LMDZ5/trunk/libf/phylmd/lmdz1d_sub.F @ 2017

Last change on this file since 2017 was 2017, checked in by fhourdin, 10 years ago

Reintegration de phy1d dans phylmd

File size: 35.3 KB
Line 
1
2      SUBROUTINE lmdz1d
3
4      USE ioipsl, only: ju2ymds, ymds2ju, ioconf_calendar
5      use phys_state_var_mod
6      use comgeomphy
7      use dimphy
8      use surface_data, only : type_ocean,ok_veget
9      use pbl_surface_mod, only : ftsoil, pbl_surface_init,
10     $                            pbl_surface_final
11      use fonte_neige_mod, only : fonte_neige_init, fonte_neige_final
12
13      use infotrac ! new
14      use control_mod
15      USE indice_sol_mod
16      USE phyaqua_mod
17
18      implicit none
19#include "dimensions.h"
20#include "YOMCST.h"
21#include "temps.h"
22!!#include "control.h"
23#include "iniprint.h"
24#include "clesphys.h"
25#include "dimsoil.h"
26!#include "indicesol.h"
27
28#include "comvert.h"
29#include "compar1d.h"
30#include "flux_arp.h"
31#include "tsoilnudge.h"
32#include "fcg_gcssold.h"
33!!!#include "fbforcing.h"
34
35!=====================================================================
36! DECLARATIONS
37!=====================================================================
38
39!---------------------------------------------------------------------
40!  Externals
41!---------------------------------------------------------------------
42      external fq_sat
43      real fq_sat
44
45!---------------------------------------------------------------------
46!  Arguments d' initialisations de la physique (USER DEFINE)
47!---------------------------------------------------------------------
48
49      integer, parameter :: ngrid=1
50      real :: zcufi    = 1.
51      real :: zcvfi    = 1.
52
53!-      real :: nat_surf
54!-      logical :: ok_flux_surf
55!-      real :: fsens
56!-      real :: flat
57!-      real :: tsurf
58!-      real :: rugos
59!-      real :: qsol(1:2)
60!-      real :: qsurf
61!-      real :: psurf
62!-      real :: zsurf
63!-      real :: albedo
64!-
65!-      real :: time     = 0.
66!-      real :: time_ini
67!-      real :: xlat
68!-      real :: xlon
69!-      real :: wtsurf
70!-      real :: wqsurf
71!-      real :: restart_runoff
72!-      real :: xagesno
73!-      real :: qsolinp
74!-      real :: zpicinp
75!-
76      real :: fnday
77      real :: day, daytime
78      real :: day1
79      real :: heure
80      integer :: jour
81      integer :: mois
82      integer :: an
83 
84!---------------------------------------------------------------------
85!  Declarations related to forcing and initial profiles
86!---------------------------------------------------------------------
87
88        integer :: kmax = llm
89        integer llm700,nq1,nq2
90        INTEGER, PARAMETER :: nlev_max=1000, nqmx=1000
91        real timestep, frac
92        real height(nlev_max),tttprof(nlev_max),qtprof(nlev_max),
93     .              uprof(nlev_max),vprof(nlev_max),e12prof(nlev_max),
94     .              ugprof(nlev_max),vgprof(nlev_max),wfls(nlev_max),
95     .              dqtdxls(nlev_max),dqtdyls(nlev_max),
96     .              dqtdtls(nlev_max),thlpcar(nlev_max),
97     .              qprof(nlev_max,nqmx)
98
99c        integer :: forcing_type
100        logical :: forcing_les     = .false.
101        logical :: forcing_armcu   = .false.
102        logical :: forcing_rico    = .false.
103        logical :: forcing_radconv = .false.
104        logical :: forcing_toga    = .false.
105        logical :: forcing_twpice  = .false.
106        logical :: forcing_amma    = .false.
107        logical :: forcing_GCM2SCM = .false.
108        logical :: forcing_GCSSold = .false.
109        logical :: forcing_sandu   = .false.
110        logical :: forcing_astex   = .false.
111        logical :: forcing_fire    = .false.
112        integer :: type_ts_forcing ! 0 = SST constant; 1 = SST read from a file
113!                                                            (cf read_tsurf1d.F)
114
115!vertical advection computation
116!       real d_t_z(llm), d_q_z(llm)
117!       real d_t_dyn_z(llm), d_q_dyn_z(llm)
118!       real zz(llm)
119!       real zfact
120
121!flag forcings
122        logical :: nudge_wind=.true.
123        logical :: nudge_thermo=.false.
124        logical :: cptadvw=.true.
125!=====================================================================
126! DECLARATIONS FOR EACH CASE
127!=====================================================================
128!
129#include "1D_decl_cases.h"
130!
131!---------------------------------------------------------------------
132!  Declarations related to vertical discretization:
133!---------------------------------------------------------------------
134      real :: pzero=1.e5
135      real :: play (llm),zlay (llm),sig_s(llm),plev(llm+1)
136      real :: playd(llm),zlayd(llm),ap_amma(llm+1),bp_amma(llm+1),poub
137
138!---------------------------------------------------------------------
139!  Declarations related to variables
140!---------------------------------------------------------------------
141
142      real :: phi(llm)
143      real :: teta(llm),tetal(llm),temp(llm),u(llm),v(llm),w(llm)
144      real :: rlat_rad(1),rlon_rad(1)
145      real :: omega(llm+1),omega2(llm),rho(llm+1)
146      real :: ug(llm),vg(llm),fcoriolis
147      real :: sfdt, cfdt
148      real :: du_phys(llm),dv_phys(llm),dt_phys(llm)
149      real :: dt_dyn(llm)
150      real :: dt_cooling(llm),d_th_adv(llm)
151      real :: alpha
152      real :: ttt
153
154      REAL, ALLOCATABLE, DIMENSION(:,:):: q
155      REAL, ALLOCATABLE, DIMENSION(:,:):: dq
156      REAL, ALLOCATABLE, DIMENSION(:,:):: dq_dyn
157      REAL, ALLOCATABLE, DIMENSION(:,:):: d_q_adv
158!     REAL, ALLOCATABLE, DIMENSION(:):: d_th_adv
159
160!---------------------------------------------------------------------
161!  Initialization of surface variables
162!---------------------------------------------------------------------
163      real :: run_off_lic_0(1)
164      real :: fder(1),snsrf(1,nbsrf),qsurfsrf(1,nbsrf)
165      real :: evap(1,nbsrf),frugs(1,nbsrf)
166      real :: tsoil(1,nsoilmx,nbsrf)
167      real :: agesno(1,nbsrf)
168
169!---------------------------------------------------------------------
170!  Call to phyredem
171!---------------------------------------------------------------------
172      logical :: ok_writedem =.true.
173     
174!---------------------------------------------------------------------
175!  Call to physiq
176!---------------------------------------------------------------------
177      integer, parameter :: longcles=20
178      logical :: firstcall=.true.
179      logical :: lastcall=.false.
180      real :: phis    = 0.0
181      real :: clesphy0(longcles) = 0.0
182      real :: dpsrf
183
184!---------------------------------------------------------------------
185!  Initializations of boundary conditions
186!---------------------------------------------------------------------
187      integer, parameter :: yd = 360
188      real :: phy_nat (yd) = 0.0 ! 0=ocean libre,1=land,2=glacier,3=banquise
189      real :: phy_alb (yd)  ! Albedo land only (old value condsurf_jyg=0.3)
190      real :: phy_sst (yd)  ! SST (will not be used; cf read_tsurf1d.F)
191      real :: phy_bil (yd) = 1.0 ! Ne sert que pour les slab_ocean
192      real :: phy_rug (yd) ! Longueur rugosite utilisee sur land only
193      real :: phy_ice (yd) = 0.0 ! Fraction de glace
194      real :: phy_fter(yd) = 0.0 ! Fraction de terre
195      real :: phy_foce(yd) = 0.0 ! Fraction de ocean
196      real :: phy_fsic(yd) = 0.0 ! Fraction de glace
197      real :: phy_flic(yd) = 0.0 ! Fraction de glace
198
199!---------------------------------------------------------------------
200!  Fichiers et d'autres variables
201!---------------------------------------------------------------------
202      integer :: k,l,i,it=1,mxcalc
203      integer jjmp1
204      parameter (jjmp1=jjm+1-1/jjm)
205      INTEGER nbteta
206      PARAMETER(nbteta=3)
207      REAL dudyn(iim+1,jjmp1,llm)
208      REAL PVteta(1,nbteta)
209      INTEGER read_climoz
210!Al1
211      integer ecrit_slab_oc !1=ecrit,-1=lit,0=no file
212      data ecrit_slab_oc/-1/
213
214!=====================================================================
215! INITIALIZATIONS
216!=====================================================================
217! Initialization of Common turb_forcing
218       dtime_frcg = 0.
219       Turb_fcg_gcssold=.false.
220       hthturb_gcssold = 0.
221       hqturb_gcssold = 0.
222
223!---------------------------------------------------------------------
224! OPTIONS OF THE 1D SIMULATION (lmdz1d.def => unicol.def)
225!---------------------------------------------------------------------
226cAl1
227        call conf_unicol
228cAl1 moves this gcssold var from common fcg_gcssold to
229        Turb_fcg_gcssold = xTurb_fcg_gcssold
230c --------------------------------------------------------------------
231        close(1)
232cAl1
233        write(*,*) 'lmdz1d.def lu => unicol.def'
234
235! forcing_type defines the way the SCM is forced:
236!forcing_type = 0 ==> forcing_les = .true.
237!             initial profiles from file prof.inp.001
238!             no forcing by LS convergence ;
239!             surface temperature imposed ;
240!             radiative cooling may be imposed (iflag_radia=0 in physiq.def)
241!forcing_type = 1 ==> forcing_radconv = .true.
242!             idem forcing_type = 0, but the imposed radiative cooling
243!             is set to 0 (hence, if iflag_radia=0 in physiq.def,
244!             then there is no radiative cooling at all)
245!forcing_type = 2 ==> forcing_toga = .true.
246!             initial profiles from TOGA-COARE IFA files
247!             LS convergence and SST imposed from TOGA-COARE IFA files
248!forcing_type = 3 ==> forcing_GCM2SCM = .true.
249!             initial profiles from the GCM output
250!             LS convergence imposed from the GCM output
251!forcing_type = 4 ==> forcing_twpice = .true.
252!             initial profiles from TWP-ICE cdf file
253!             LS convergence, omega and SST imposed from TWP-ICE files
254!forcing_type = 5 ==> forcing_rico = .true.
255!             initial profiles from RICO files
256!             LS convergence imposed from RICO files
257!forcing_type = 6 ==> forcing_amma = .true.
258!             initial profiles from AMMA nc file
259!             LS convergence, omega and surface fluxes imposed from AMMA file 
260!forcing_type = 40 ==> forcing_GCSSold = .true.
261!             initial profile from GCSS file
262!             LS convergence imposed from GCSS file
263!forcing_type = 50 ==> forcing_fire = .true.
264!             forcing from fire.nc
265!forcing_type = 59 ==> forcing_sandu = .true.
266!             initial profiles from sanduref file: see prof.inp.001
267!             SST varying with time and divergence constante: see ifa_sanduref.txt file
268!             Radiation has to be computed interactively
269!forcing_type = 60 ==> forcing_astex = .true.
270!             initial profiles from file: see prof.inp.001
271!             SST,divergence,ug,vg,ufa,vfa varying with time : see ifa_astex.txt file
272!             Radiation has to be computed interactively
273!forcing_type = 61 ==> forcing_armcu = .true.
274!             initial profiles from file: see prof.inp.001
275!             sensible and latent heat flux imposed: see ifa_arm_cu_1.txt
276!             large scale advective forcing & radiative tendencies applied below 1000m: see ifa_arm_cu_2.txt
277!             use geostrophic wind ug=10m/s vg=0m/s. Duration of the case 53100s
278!             Radiation to be switched off
279!
280      if (forcing_type .eq.0) THEN
281       forcing_les = .true.
282      elseif (forcing_type .eq.1) THEN
283       forcing_radconv = .true.
284      elseif (forcing_type .eq.2) THEN
285       forcing_toga    = .true.
286      elseif (forcing_type .eq.3) THEN
287       forcing_GCM2SCM = .true.
288      elseif (forcing_type .eq.4) THEN
289       forcing_twpice = .true.
290      elseif (forcing_type .eq.5) THEN
291       forcing_rico = .true.
292      elseif (forcing_type .eq.6) THEN
293       forcing_amma = .true.
294      elseif (forcing_type .eq.40) THEN
295       forcing_GCSSold = .true.
296      elseif (forcing_type .eq.50) THEN
297       forcing_fire = .true.
298      elseif (forcing_type .eq.59) THEN
299       forcing_sandu   = .true.
300      elseif (forcing_type .eq.60) THEN
301       forcing_astex   = .true.
302      elseif (forcing_type .eq.61) THEN
303       forcing_armcu = .true.
304       IF(llm.NE.19.AND.llm.NE.40) stop 'Erreur nombre de niveaux !!'
305      else
306       write (*,*) 'ERROR : unknown forcing_type ', forcing_type
307       stop 'Forcing_type should be 0,1,2,3,4,5,6 or 40,59,60,61'
308      ENDIF
309      print*,"forcing type=",forcing_type
310
311! if type_ts_forcing=0, the surface temp of 1D simulation is constant in time
312! (specified by tsurf in lmdz1d.def); if type_ts_forcing=1, the surface temperature
313! varies in time according to a forcing (e.g. forcing_toga) and is passed to read_tsurf1d.F
314! through the common sst_forcing.
315
316        type_ts_forcing = 0
317        if (forcing_toga .or. forcing_sandu .or. forcing_astex)
318     :    type_ts_forcing = 1
319
320!---------------------------------------------------------------------
321!  Definition of the run
322!---------------------------------------------------------------------
323
324      call conf_gcm( 99, .TRUE. , clesphy0 )
325c-----------------------------------------------------------------------
326c   Choix du calendrier
327c   -------------------
328
329c      calend = 'earth_365d'
330      if (calend == 'earth_360d') then
331        call ioconf_calendar('360d')
332        write(*,*)'CALENDRIER CHOISI: Terrestre a 360 jours/an'
333      else if (calend == 'earth_365d') then
334        call ioconf_calendar('noleap')
335        write(*,*)'CALENDRIER CHOISI: Terrestre a 365 jours/an'
336      else if (calend == 'earth_366d') then
337        call ioconf_calendar('all_leap')
338        write(*,*)'CALENDRIER CHOISI: Terrestre bissextile'
339      else if (calend == 'gregorian') then
340        call ioconf_calendar('gregorian') ! not to be used by normal users
341        write(*,*)'CALENDRIER CHOISI: Gregorien'
342      else
343        write (*,*) 'ERROR : unknown calendar ', calend
344        stop 'calend should be 360d,earth_365d,earth_366d,gregorian'
345      endif
346c-----------------------------------------------------------------------
347c
348cc Date :
349c      La date est supposee donnee sous la forme [annee, numero du jour dans
350c      l annee] ; l heure est donnee dans time_ini, lu dans lmdz1d.def.
351c      On appelle ymds2ju pour convertir [annee, jour] en [jour Julien].
352c      Le numero du jour est dans "day". L heure est traitee separement.
353c      La date complete est dans "daytime" (l'unite est le jour).
354      if (nday>0) then
355         fnday=nday
356      else
357         fnday=-nday/float(day_step)
358      endif
359
360c Special case for arm_cu which lasts less than one day : 53100s !! (MPL 20111026)
361      IF(forcing_type .EQ. 61) fnday=53100./86400.
362c Special case for amma which lasts less than one day : 64800s !! (MPL 20120216)
363      IF(forcing_type .EQ. 6) fnday=64800./86400.
364      annee_ref = anneeref
365      mois = 1
366      day_ref = dayref
367      heure = 0.
368      itau_dyn = 0
369      itau_phy = 0
370      call ymds2ju(annee_ref,mois,day_ref,heure,day)
371      day_ini = int(day)
372      day_end = day_ini + fnday
373
374      IF (forcing_type .eq.2) THEN
375! Convert the initial date of Toga-Coare to Julian day
376      call ymds2ju
377     $ (year_ini_toga,mth_ini_toga,day_ini_toga,heure,day_ju_ini_toga)
378
379      ELSEIF (forcing_type .eq.4) THEN
380! Convert the initial date of TWPICE to Julian day
381      call ymds2ju
382     $ (year_ini_twpi,mth_ini_twpi,day_ini_twpi,heure_ini_twpi
383     $ ,day_ju_ini_twpi)
384      ELSEIF (forcing_type .eq.6) THEN
385! Convert the initial date of AMMA to Julian day
386      call ymds2ju
387     $ (year_ini_amma,mth_ini_amma,day_ini_amma,heure_ini_amma
388     $ ,day_ju_ini_amma)
389
390      ELSEIF (forcing_type .eq.59) THEN
391! Convert the initial date of Sandu case to Julian day
392      call ymds2ju
393     $   (year_ini_sandu,mth_ini_sandu,day_ini_sandu,
394     $    time_ini*3600.,day_ju_ini_sandu)
395
396      ELSEIF (forcing_type .eq.60) THEN
397! Convert the initial date of Astex case to Julian day
398      call ymds2ju
399     $   (year_ini_astex,mth_ini_astex,day_ini_astex,
400     $    time_ini*3600.,day_ju_ini_astex)
401
402      ELSEIF (forcing_type .eq.61) THEN
403
404! Convert the initial date of Arm_cu case to Julian day
405      call ymds2ju
406     $ (year_ini_armcu,mth_ini_armcu,day_ini_armcu,heure_ini_armcu
407     $ ,day_ju_ini_armcu)
408      ENDIF
409
410      daytime = day + time_ini/24. ! 1st day and initial time of the simulation
411! Print out the actual date of the beginning of the simulation :
412      call ju2ymds(daytime,year_print, month_print,day_print,sec_print)
413      print *,' Time of beginning : ',
414     $        year_print, month_print, day_print, sec_print
415
416!---------------------------------------------------------------------
417! Initialization of dimensions, geometry and initial state
418!---------------------------------------------------------------------
419      call init_phys_lmdz(1,1,llm,1,(/1/))
420      call suphel
421      call initcomgeomphy
422      call infotrac_init
423
424      if (nqtot>nqmx) STOP'Augmenter nqmx dans lmdz1d.F'
425      allocate(q(llm,nqtot)) ; q(:,:)=0.
426      allocate(dq(llm,nqtot))
427      allocate(dq_dyn(llm,nqtot))
428      allocate(d_q_adv(llm,nqtot))
429!     allocate(d_th_adv(llm))
430
431c
432c   No ozone climatology need be read in this pre-initialization
433c          (phys_state_var_init is called again in physiq)
434      read_climoz = 0
435c
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      rlat=xlat
458      rlon=xlon
459      day1= day_ini
460      time=daytime-day
461      ts_toga(1)=tsurf ! needed by read_tsurf1d.F
462      rho(1)=psurf/(rd*tsurf*(1.+(rv/rd-1.)*qsurf))
463
464!
465!! mpl et jyg le 22/08/2012 :
466!!  pour que les cas a flux de surface imposes marchent
467      IF(.NOT.ok_flux_surf.or.max(abs(wtsurf),abs(wqsurf))>0.) THEN
468       fsens=-wtsurf*rcpd*rho(1)
469       flat=-wqsurf*rlvtt*rho(1)
470       print *,'Flux: ok_flux wtsurf wqsurf',ok_flux_surf,wtsurf,wqsurf
471      ENDIF
472      print*,'Flux sol ',fsens,flat
473!!      ok_flux_surf=.false.
474!!      fsens=-wtsurf*rcpd*rho(1)
475!!      flat=-wqsurf*rlvtt*rho(1)
476!!!!
477
478! Vertical discretization and pressure levels at half and mid levels:
479
480      pa   = 5e4
481!!      preff= 1.01325e5
482      preff = psurf
483      IF (ok_old_disvert) THEN
484        call disvert0(pa,preff,ap,bp,dpres,presnivs,nivsigs,nivsig)
485        print *,'On utilise disvert0'
486      ELSE
487        call disvert(pa,preff,ap,bp,dpres,presnivs,nivsigs,nivsig,
488     :                 scaleheight)
489        print *,'On utilise disvert'
490c       Nouvelle version disvert permettant d imposer ap,bp (modif L.Guez) MPL 18092012
491c       Dans ce cas, on lit ap,bp dans le fichier hybrid.txt
492      ENDIF
493      sig_s=presnivs/preff
494      plev =ap+bp*psurf
495      play = 0.5*(plev(1:llm)+plev(2:llm+1))
496ccc      zlay=-rd*300.*log(play/psurf)/rg ! moved after reading profiles
497
498      IF (forcing_type .eq. 59) THEN
499! pour forcing_sandu, on cherche l'indice le plus proche de 700hpa#3000m
500      write(*,*) '***********************'
501      do l = 1, llm
502       write(*,*) 'l,play(l),presnivs(l): ',l,play(l),presnivs(l)
503       if (trouve_700 .and. play(l).le.70000) then
504         llm700=l
505         print *,'llm700,play=',llm700,play(l)/100.
506         trouve_700= .false.
507       endif
508      enddo
509      write(*,*) '***********************'
510      ENDIF
511
512c
513!=====================================================================
514! EVENTUALLY, READ FORCING DATA :
515!=====================================================================
516
517#include "1D_read_forc_cases.h"
518
519      if (forcing_GCM2SCM) then
520        write (*,*) 'forcing_GCM2SCM not yet implemented'
521        stop 'in initialization'
522      endif ! forcing_GCM2SCM
523
524      print*,'mxcalc=',mxcalc
525      print*,'zlay=',zlay(mxcalc)
526      print*,'play=',play(mxcalc)
527
528cAl1 pour SST forced, appellé depuis ocean_forced_noice
529      ts_cur = tsurf ! SST used in read_tsurf1d
530!=====================================================================
531! Initialisation de la physique :
532!=====================================================================
533
534!  Rq: conf_phys.F90 lit tous les flags de physiq.def; conf_phys appele depuis physiq.F
535!
536! day_step, iphysiq lus dans gcm.def ci-dessus
537! timestep: calcule ci-dessous from rday et day_step
538! ngrid=1
539! llm: defini dans .../modipsl/modeles/LMDZ4/libf/grid/dimension
540! rday: defini dans suphel.F (86400.)
541! day_ini: lu dans run.def (dayref)
542! rlat_rad,rlon-rad: transformes en radian de rlat,rlon lus dans lmdz1d.def (en degres)
543! airefi,zcufi,zcvfi initialises au debut de ce programme
544! rday,ra,rg,rd,rcpd declares dans YOMCST.h et calcules dans suphel.F
545      day_step = float(nsplit_phys)*day_step/float(iphysiq)
546      write (*,*) 'Time step divided by nsplit_phys (=',nsplit_phys,')'
547      timestep =rday/day_step
548      dtime_frcg = timestep
549!
550      zcufi=airefi
551      zcvfi=airefi
552!
553      rlat_rad(:)=rlat(:)*rpi/180.
554      rlon_rad(:)=rlon(:)*rpi/180.
555
556      call iniphysiq(ngrid,llm,rday,day_ini,timestep,
557     .     rlat_rad,rlon_rad,airefi,zcufi,zcvfi,ra,rg,rd,rcpd,(/1/))
558      print*,'apres iniphysiq'
559
560! 2 PARAMETRES QUI DEVRAIENT ETRE LUS DANS run.def MAIS NE LE SONT PAS ICI:
561      co2_ppm= 330.0
562      solaire=1370.0
563
564! Ecriture du startphy avant le premier appel a la physique.
565! On le met juste avant pour avoir acces a tous les champs
566! NB: les clesphy0 seront remplies dans phyredem d'apres les flags lus dans gcm.def
567
568      if (ok_writedem) then
569
570!--------------------------------------------------------------------------
571! pbl_surface_init (called here) and pbl_surface_final (called by phyredem)
572! need : qsol fder snow qsurf evap rugos agesno ftsoil
573!--------------------------------------------------------------------------
574
575        type_ocean = "force"
576        run_off_lic_0(1) = restart_runoff
577        call fonte_neige_init(run_off_lic_0)
578
579        fder=0.
580        snsrf(1,:)=0.        ! couverture de neige des sous surface
581        qsurfsrf(1,:)=qsurf ! humidite de l'air des sous surface
582        evap=0.
583        frugs(1,:)=rugos     ! couverture de neige des sous surface
584        agesno  = xagesno
585        tsoil(:,:,:)=tsurf
586!------ AMMA 2e run avec modele sol et rayonnement actif (MPL 23052012)
587!       tsoil(1,1,1)=299.18
588!       tsoil(1,2,1)=300.08
589!       tsoil(1,3,1)=301.88
590!       tsoil(1,4,1)=305.48
591!       tsoil(1,5,1)=308.00
592!       tsoil(1,6,1)=308.00
593!       tsoil(1,7,1)=308.00
594!       tsoil(1,8,1)=308.00
595!       tsoil(1,9,1)=308.00
596!       tsoil(1,10,1)=308.00
597!       tsoil(1,11,1)=308.00
598!-----------------------------------------------------------------------
599        call pbl_surface_init(qsol, fder, snsrf, qsurfsrf,
600     &                                    evap, frugs, agesno, tsoil)
601
602!------------------ prepare limit conditions for limit.nc -----------------
603!--   Ocean force
604
605        print*,'avant phyredem'
606        pctsrf(1,:)=0.
607        if (nat_surf.eq.0.) then
608          pctsrf(1,is_oce)=1.
609          pctsrf(1,is_ter)=0.
610        else
611          pctsrf(1,is_oce)=0.
612          pctsrf(1,is_ter)=1.
613        end if
614
615        print*,'nat_surf,pctsrf(1,is_oce),pctsrf(1,is_ter)',nat_surf
616     $        ,pctsrf(1,is_oce),pctsrf(1,is_ter)
617
618        zmasq=pctsrf(1,is_ter)+pctsrf(1,is_lic)
619        zpic = zpicinp
620        ftsol=tsurf
621        falb1 = albedo                           
622        falb2 = albedo                           
623        rugoro=rugos
624        t_ancien(1,:)=temp(:)
625        q_ancien(1,:)=q(:,1)
626        pbl_tke=1.e-8
627
628        rain_fall=0.
629        snow_fall=0.
630        solsw=0.
631        sollw=0.
632        radsol=0.
633        rnebcon=0.
634        ratqs=0.
635        clwcon=0.
636        zmea=0.
637        zstd=0.
638        zsig=0.
639        zgam=0.
640        zval=0.
641        zthe=0.
642        sig1=0.
643        w01=0.
644        u_ancien(1,:)=u(:)
645        v_ancien(1,:)=v(:)
646 
647!------------------------------------------------------------------------
648! Make file containing restart for the physics (startphy.nc)
649!
650! NB: List of the variables to be written by phyredem (via put_field):
651! rlon,rlat,zmasq,pctsrf(:,is_ter),pctsrf(:,is_lic),pctsrf(:,is_oce)
652! pctsrf(:,is_sic),ftsol(:,nsrf),tsoil(:,isoil,nsrf),qsurf(:,nsrf)
653! qsol,falb1(:,nsrf),falb2(:,nsrf),evap(:,nsrf),snow(:,nsrf)
654! radsol,solsw,sollw,fder,rain_fall,snow_fall,frugs(:,nsrf)
655! agesno(:,nsrf),zmea,zstd,zsig,zgam,zthe,zpic,zval,rugoro
656! t_ancien,q_ancien,frugs(:,is_oce),clwcon(:,1),rnebcon(:,1),ratqs(:,1)
657! run_off_lic_0,pbl_tke(:,1:klev,nsrf),zmax0,f0,sig1,w01
658! wake_deltat,wake_deltaq,wake_s,wake_cstar,wake_fip
659!------------------------------------------------------------------------
660CAl1 =============== restart option ==========================
661        if (.not.restart) then
662          call phyredem ("startphy.nc")
663        else
664c (desallocations)
665        print*,'callin surf final'
666          call pbl_surface_final(qsol, fder, snsrf, qsurfsrf,
667     &                                    evap, frugs, agesno, tsoil)
668        print*,'after surf final'
669          CALL fonte_neige_final(run_off_lic_0)
670        endif
671
672        ok_writedem=.false.
673        print*,'apres phyredem'
674
675      endif ! ok_writedem
676     
677!------------------------------------------------------------------------
678! Make file containing boundary conditions (limit.nc) **Al1->restartdyn***
679! --------------------------------------------------
680! NB: List of the variables to be written in limit.nc
681!     (by writelim.F, subroutine of 1DUTILS.h):
682!        phy_nat,phy_alb,phy_sst,phy_bil,phy_rug,phy_ice,
683!        phy_fter,phy_foce,phy_flic,phy_fsic)
684!------------------------------------------------------------------------
685      do i=1,yd
686        phy_nat(i)  = nat_surf
687        phy_alb(i)  = albedo
688        phy_sst(i)  = tsurf ! read_tsurf1d will be used instead
689        phy_rug(i)  = rugos
690        phy_fter(i) = pctsrf(1,is_ter)
691        phy_foce(i) = pctsrf(1,is_oce)
692        phy_fsic(i) = pctsrf(1,is_sic)
693        phy_flic(i) = pctsrf(1,is_lic)
694      enddo
695
696C fabrication de limit.nc
697      call writelim (1,phy_nat,phy_alb,phy_sst,phy_bil,phy_rug,
698     &               phy_ice,phy_fter,phy_foce,phy_flic,phy_fsic)
699
700
701      call phys_state_var_end
702cAl1
703      if (restart) then
704        print*,'call to restart dyn 1d'
705        Call dyn1deta0("start1dyn.nc",
706     &              plev,play,phi,phis,presnivs,
707     &              u,v,temp,q,omega2)
708
709       print*,'fnday,annee_ref,day_ref,day_ini',
710     &     fnday,annee_ref,day_ref,day_ini
711c**      call ymds2ju(annee_ref,mois,day_ini,heure,day)
712       day = day_ini
713       day_end = day_ini + nday
714       daytime = day + time_ini/24. ! 1st day and initial time of the simulation
715
716! Print out the actual date of the beginning of the simulation :
717       call ju2ymds(daytime, an, mois, jour, heure)
718       print *,' Time of beginning : y m d h',an, mois,jour,heure/3600.
719
720       day = int(daytime)
721       time=daytime-day
722 
723       print*,'****** intialised fields from restart1dyn *******'
724       print*,'plev,play,phi,phis,presnivs,u,v,temp,q,omega2'
725       print*,'temp(1),q(1,1),u(1),v(1),plev(1),phis :'
726       print*,temp(1),q(1,1),u(1),v(1),plev(1),phis
727c raz for safety
728       do l=1,llm
729         dq_dyn(l,1) = 0.
730       enddo
731      endif
732!Al1 ================  end restart =================================
733      IF (ecrit_slab_oc.eq.1) then
734         open(97,file='div_slab.dat',STATUS='UNKNOWN')
735       elseif (ecrit_slab_oc.eq.0) then
736         open(97,file='div_slab.dat',STATUS='OLD')
737       endif
738!=====================================================================
739! START OF THE TEMPORAL LOOP :
740!=====================================================================
741           
742      do while(it.le.nint(fnday*day_step))
743
744       if (prt_level.ge.1) then
745         print*,'XXXXXXXXXXXXXXXXXXX ITAP,day,time=',
746     .                                it,day,time,nint(fnday*day_step)
747         print*,'PAS DE TEMPS ',timestep
748       endif
749!Al1 demande de restartphy.nc
750       if (it.eq.nint(fnday*day_step)) lastcall=.True.
751
752!---------------------------------------------------------------------
753! Interpolation of forcings in time and onto model levels
754!---------------------------------------------------------------------
755
756#include "1D_interp_cases.h"
757
758      if (forcing_GCM2SCM) then
759        write (*,*) 'forcing_GCM2SCM not yet implemented'
760        stop 'in time loop'
761      endif ! forcing_GCM2SCM
762
763!---------------------------------------------------------------------
764!  Geopotential :
765!---------------------------------------------------------------------
766
767        phi(1)=RD*temp(1)*(plev(1)-play(1))/(.5*(plev(1)+play(1)))
768        do l = 1, llm-1
769          phi(l+1)=phi(l)+RD*(temp(l)+temp(l+1))*
770     .    (play(l)-play(l+1))/(play(l)+play(l+1))
771        enddo
772
773!---------------------------------------------------------------------
774! Listing output for debug prt_level>=1
775!---------------------------------------------------------------------
776       if (prt_level>=1) then
777         print *,' avant physiq : -------- day time ',day,time
778         write(*,*) 'firstcall,lastcall,phis',
779     :               firstcall,lastcall,phis
780         write(*,'(a10,2a4,4a13)') 'BEFOR1 IT=','it','l',
781     :        'presniv','plev','play','phi'
782         write(*,'(a10,2i4,4f13.2)') ('BEFOR1 IT= ',it,l,
783     :         presnivs(l),plev(l),play(l),phi(l),l=1,llm)
784         write(*,'(a11,2a4,a11,6a8)') 'BEFOR2','it','l',
785     :         'presniv','u','v','temp','q1','q2','omega2'
786         write(*,'(a11,2i4,f11.2,5f8.2,e10.2)') ('BEFOR2 IT= ',it,l,
787     :   presnivs(l),u(l),v(l),temp(l),q(l,1),q(l,2),omega2(l),l=1,llm)
788       endif
789
790!---------------------------------------------------------------------
791!   Call physiq :
792!---------------------------------------------------------------------
793
794        call physiq(ngrid,llm,
795     :              firstcall,lastcall,
796     :              day,time,timestep,
797     :              plev,play,phi,phis,presnivs,clesphy0,
798     :              u,v,temp,q,omega2,
799     :              du_phys,dv_phys,dt_phys,dq,dpsrf,
800     :              dudyn,PVteta)
801        firstcall=.false.
802
803!---------------------------------------------------------------------
804! Listing output for debug prt_level>=1
805!---------------------------------------------------------------------
806        if (prt_level>=1) then
807          write(*,'(a11,2a4,4a13)') 'AFTER1 IT=','it','l',
808     :        'presniv','plev','play','phi'
809          write(*,'(a11,2i4,4f13.2)') ('AFTER1 it= ',it,l,
810     :    presnivs(l),plev(l),play(l),phi(l),l=1,llm)
811          write(*,'(a11,2a4,a11,6a8)') 'AFTER2','it','l',
812     :         'presniv','u','v','temp','q1','q2','omega2'
813          write(*,'(a11,2i4,f11.2,5f8.2,e10.2)') ('AFTER2 it= ',it,l,
814     :    presnivs(l),u(l),v(l),temp(l),q(l,1),q(l,2),omega2(l),l=1,llm)
815          write(*,'(a11,2a4,a11,5a8)') 'AFTER3','it','l',
816     :         'presniv','du_phys','dv_phys','dt_phys','dq1','dq2'
817           write(*,'(a11,2i4,f11.2,5f8.2)') ('AFTER3 it= ',it,l,
818     :      presnivs(l),86400*du_phys(l),86400*dv_phys(l),
819     :       86400*dt_phys(l),86400*dq(l,1),dq(l,2),l=1,llm)
820          write(*,*) 'dpsrf',dpsrf
821        endif
822!---------------------------------------------------------------------
823!   Add physical tendencies :
824!---------------------------------------------------------------------
825
826       fcoriolis=2.*sin(rpi*xlat/180.)*romega
827       if (forcing_radconv .or. forcing_fire) then
828         fcoriolis=0.0
829         dt_cooling=0.0
830         d_th_adv=0.0
831         d_q_adv=0.0
832       endif
833
834       if (forcing_toga .or. forcing_GCSSold .or. forcing_twpice
835     :    .or.forcing_amma) then
836         fcoriolis=0.0 ; ug=0. ; vg=0.
837       endif
838         if(forcing_rico) then
839          dt_cooling=0.
840        endif
841
842      print*, 'fcoriolis ', fcoriolis, xlat
843
844        du_age(1:mxcalc)=
845     : fcoriolis*(v(1:mxcalc)-vg(1:mxcalc))
846       dv_age(1:mxcalc)=
847     : -fcoriolis*(u(1:mxcalc)-ug(1:mxcalc))
848
849!!!!!!!!!!!!!!!!!!!!!!!!
850! Geostrophic wind
851!!!!!!!!!!!!!!!!!!!!!!!!
852       sfdt = sin(0.5*fcoriolis*timestep)
853       cfdt = cos(0.5*fcoriolis*timestep)
854!
855        du_age(1:mxcalc)= -2.*sfdt/timestep*
856     s          (sfdt*(u(1:mxcalc)-ug(1:mxcalc)) -
857     s           cfdt*(v(1:mxcalc)-vg(1:mxcalc))  )
858!!     : fcoriolis*(v(1:mxcalc)-vg(1:mxcalc))
859!
860       dv_age(1:mxcalc)= -2.*sfdt/timestep*
861     s          (cfdt*(u(1:mxcalc)-ug(1:mxcalc)) +
862     s           sfdt*(v(1:mxcalc)-vg(1:mxcalc))  )
863!!     : -fcoriolis*(u(1:mxcalc)-ug(1:mxcalc))
864!
865!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
866!         call  writefield_phy('dv_age' ,dv_age,llm)
867!         call  writefield_phy('du_age' ,du_age,llm)
868!         call  writefield_phy('du_phys' ,du_phys,llm)
869!         call  writefield_phy('u_tend' ,u,llm)
870!         call  writefield_phy('u_g' ,ug,llm)
871!
872!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
873!! Increment state variables
874!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
875! pour les cas sandu et astex, on reclacule u,v,q,temp et teta dans 1D_nudge_sandu_astex.h
876! au dessus de 700hpa, on relaxe vers les profils initiaux
877      if (forcing_sandu .OR. forcing_astex) then
878#include "1D_nudge_sandu_astex.h"
879      else
880        u(1:mxcalc)=u(1:mxcalc) + timestep*(
881     s              du_phys(1:mxcalc)
882     s             +du_age(1:mxcalc) )           
883        v(1:mxcalc)=v(1:mxcalc) + timestep*(
884     s              dv_phys(1:mxcalc)
885     s             +dv_age(1:mxcalc) )
886        q(1:mxcalc,:)=q(1:mxcalc,:)+timestep*(
887     s                dq(1:mxcalc,:)
888     s               +d_q_adv(1:mxcalc,:) )
889
890        if (prt_level.ge.1) then
891          print *,
892     :    'physiq-> temp(1),dt_phys(1),d_th_adv(1),dt_cooling(1) ',
893     :              temp(1),dt_phys(1),d_th_adv(1),dt_cooling(1)
894           print*,du_phys
895           print*, v
896           print*, vg
897        endif
898
899        temp(1:mxcalc)=temp(1:mxcalc)+timestep*(
900     .              dt_phys(1:mxcalc)
901     .             +d_th_adv(1:mxcalc)
902     .             +dt_cooling(1:mxcalc))  ! Taux de chauffage ou refroid.
903
904      endif  ! forcing_sandu or forcing_astex
905
906        teta=temp*(pzero/play)**rkappa
907!
908!---------------------------------------------------------------------
909!   Nudge soil temperature if requested
910!---------------------------------------------------------------------
911
912      IF (nudge_tsoil) THEN
913       ftsoil(1,isoil_nudge,:) = ftsoil(1,isoil_nudge,:)
914     .  -timestep/tau_soil_nudge*(ftsoil(1,isoil_nudge,:)-Tsoil_nudge)
915      ENDIF
916
917!---------------------------------------------------------------------
918!   Add large-scale tendencies (advection, etc) :
919!---------------------------------------------------------------------
920
921ccc nrlmd
922ccc        tmpvar=teta
923ccc        call advect_vert(llm,omega,timestep,tmpvar,plev)
924ccc
925ccc        teta(1:mxcalc)=tmpvar(1:mxcalc)
926ccc        tmpvar(:)=q(:,1)
927ccc        call advect_vert(llm,omega,timestep,tmpvar,plev)
928ccc        q(1:mxcalc,1)=tmpvar(1:mxcalc)
929ccc        tmpvar(:)=q(:,2)
930ccc        call advect_vert(llm,omega,timestep,tmpvar,plev)
931ccc        q(1:mxcalc,2)=tmpvar(1:mxcalc)
932
933!---------------------------------------------------------------------
934!   Air temperature :
935!---------------------------------------------------------------------       
936        if (lastcall) then
937          print*,'Pas de temps final ',it
938          call ju2ymds(daytime, an, mois, jour, heure)
939          print*,'a la date : a m j h',an, mois, jour ,heure/3600.
940        endif
941
942c  incremente day time
943c        print*,'daytime bef',daytime,1./day_step
944        daytime = daytime+1./day_step
945!Al1dbg
946        day = int(daytime+0.1/day_step)
947c        time = max(daytime-day,0.0)
948cAl1&jyg: correction de bug
949ccc        time = real(mod(it,day_step))/day_step
950        time = time_ini/24.+real(mod(it,day_step))/day_step
951c        print*,'daytime nxt time',daytime,time
952        it=it+1
953
954      enddo
955
956!Al1
957      if (ecrit_slab_oc.ne.-1) close(97)
958
959!Al1 Call to 1D equivalent of dynredem (an,mois,jour,heure ?)
960! -------------------------------------
961       call dyn1dredem("restart1dyn.nc",
962     :              plev,play,phi,phis,presnivs,
963     :              u,v,temp,q,omega2)
964
965        CALL abort_gcm ('lmdz1d   ','The End  ',0)
966
967      end
968
969#include "1DUTILS.h"
970#include "1Dconv.h"
971
972
Note: See TracBrowser for help on using the repository browser.