source: LMDZ5/trunk/libf/phy1d/lmdz1d.F @ 1948

Last change on this file since 1948 was 1948, checked in by idelkadi, 11 years ago
  • Option pour tourner sans thermiques ni ajustement sec si ifalg_thermals<0 (phylmd/physiq.F90)
  • Modifications permettant d'imposer des profils initiaux de traceurs

lus dans un tracer.inp.001 (phylmd/1DUTILS.h, 1D_read_forc_cases.h et lmdz1d.F)

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