source: LMDZ5/branches/testing/libf/phy1d/lmdz1d.F @ 1921

Last change on this file since 1921 was 1921, checked in by Laurent Fairhead, 10 years ago

Merged trunk changes r1909:1920 into testing branch

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