source: LMDZ5/trunk/libf/phylmd/lmdz1d.F90 @ 2023

Last change on this file since 2023 was 2023, checked in by Ehouarn Millour, 10 years ago

Quick fix to be once more able to compile gcm with makelmdz_fcm.
EM

File size: 37.1 KB
Line 
1#ifdef CPP_1D
2#include "../dyn3d/mod_const_mpi.F90"
3#include "../dyn3d_common/control_mod.F90"
4#include "../dyn3d_common/infotrac.F90"
5#include "../dyn3d_common/disvert.F90"
6
7
8      PROGRAM lmdz1d
9
10      USE ioipsl, only: ju2ymds, ymds2ju, ioconf_calendar
11      use phys_state_var_mod
12      use comgeomphy
13      use dimphy
14      use surface_data, only : type_ocean,ok_veget
15      use pbl_surface_mod, only : ftsoil, pbl_surface_init,                     &
16     &                            pbl_surface_final
17      use fonte_neige_mod, only : fonte_neige_init, fonte_neige_final
18
19      use infotrac ! new
20      use control_mod
21      USE indice_sol_mod
22      USE phyaqua_mod
23
24      implicit none
25#include "dimensions.h"
26#include "YOMCST.h"
27#include "temps.h"
28!!#include "control.h"
29#include "iniprint.h"
30#include "clesphys.h"
31#include "dimsoil.h"
32!#include "indicesol.h"
33
34#include "comvert.h"
35#include "compar1d.h"
36#include "flux_arp.h"
37#include "tsoilnudge.h"
38#include "fcg_gcssold.h"
39!!!#include "fbforcing.h"
40
41!=====================================================================
42! DECLARATIONS
43!=====================================================================
44
45!---------------------------------------------------------------------
46!  Externals
47!---------------------------------------------------------------------
48      external fq_sat
49      real fq_sat
50
51!---------------------------------------------------------------------
52!  Arguments d' initialisations de la physique (USER DEFINE)
53!---------------------------------------------------------------------
54
55      integer, parameter :: ngrid=1
56      real :: zcufi    = 1.
57      real :: zcvfi    = 1.
58
59!-      real :: nat_surf
60!-      logical :: ok_flux_surf
61!-      real :: fsens
62!-      real :: flat
63!-      real :: tsurf
64!-      real :: rugos
65!-      real :: qsol(1:2)
66!-      real :: qsurf
67!-      real :: psurf
68!-      real :: zsurf
69!-      real :: albedo
70!-
71!-      real :: time     = 0.
72!-      real :: time_ini
73!-      real :: xlat
74!-      real :: xlon
75!-      real :: wtsurf
76!-      real :: wqsurf
77!-      real :: restart_runoff
78!-      real :: xagesno
79!-      real :: qsolinp
80!-      real :: zpicinp
81!-
82      real :: fnday
83      real :: day, daytime
84      real :: day1
85      real :: heure
86      integer :: jour
87      integer :: mois
88      integer :: an
89 
90!---------------------------------------------------------------------
91!  Declarations related to forcing and initial profiles
92!---------------------------------------------------------------------
93
94        integer :: kmax = llm
95        integer llm700,nq1,nq2
96        INTEGER, PARAMETER :: nlev_max=1000, nqmx=1000
97        real timestep, frac
98        real height(nlev_max),tttprof(nlev_max),qtprof(nlev_max)
99        real  uprof(nlev_max),vprof(nlev_max),e12prof(nlev_max)
100        real  ugprof(nlev_max),vgprof(nlev_max),wfls(nlev_max)
101        real  dqtdxls(nlev_max),dqtdyls(nlev_max)
102        real  dqtdtls(nlev_max),thlpcar(nlev_max)
103        real  qprof(nlev_max,nqmx)
104
105!        integer :: forcing_type
106        logical :: forcing_les     = .false.
107        logical :: forcing_armcu   = .false.
108        logical :: forcing_rico    = .false.
109        logical :: forcing_radconv = .false.
110        logical :: forcing_toga    = .false.
111        logical :: forcing_twpice  = .false.
112        logical :: forcing_amma    = .false.
113        logical :: forcing_GCM2SCM = .false.
114        logical :: forcing_GCSSold = .false.
115        logical :: forcing_sandu   = .false.
116        logical :: forcing_astex   = .false.
117        logical :: forcing_fire    = .false.
118        integer :: type_ts_forcing ! 0 = SST constant; 1 = SST read from a file
119!                                                            (cf read_tsurf1d.F)
120
121!vertical advection computation
122!       real d_t_z(llm), d_q_z(llm)
123!       real d_t_dyn_z(llm), d_q_dyn_z(llm)
124!       real zz(llm)
125!       real zfact
126
127!flag forcings
128        logical :: nudge_wind=.true.
129        logical :: nudge_thermo=.false.
130        logical :: cptadvw=.true.
131!=====================================================================
132! DECLARATIONS FOR EACH CASE
133!=====================================================================
134!
135#include "1D_decl_cases.h"
136!
137!---------------------------------------------------------------------
138!  Declarations related to vertical discretization:
139!---------------------------------------------------------------------
140      real :: pzero=1.e5
141      real :: play (llm),zlay (llm),sig_s(llm),plev(llm+1)
142      real :: playd(llm),zlayd(llm),ap_amma(llm+1),bp_amma(llm+1),poub
143
144!---------------------------------------------------------------------
145!  Declarations related to variables
146!---------------------------------------------------------------------
147
148      real :: phi(llm)
149      real :: teta(llm),tetal(llm),temp(llm),u(llm),v(llm),w(llm)
150      real :: rlat_rad(1),rlon_rad(1)
151      real :: omega(llm+1),omega2(llm),rho(llm+1)
152      real :: ug(llm),vg(llm),fcoriolis
153      real :: sfdt, cfdt
154      real :: du_phys(llm),dv_phys(llm),dt_phys(llm)
155      real :: dt_dyn(llm)
156      real :: dt_cooling(llm),d_th_adv(llm)
157      real :: alpha
158      real :: ttt
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!     REAL, ALLOCATABLE, DIMENSION(:):: d_th_adv
165
166!---------------------------------------------------------------------
167!  Initialization of surface variables
168!---------------------------------------------------------------------
169      real :: run_off_lic_0(1)
170      real :: fder(1),snsrf(1,nbsrf),qsurfsrf(1,nbsrf)
171      real :: evap(1,nbsrf),frugs(1,nbsrf)
172      real :: tsoil(1,nsoilmx,nbsrf)
173      real :: agesno(1,nbsrf)
174
175!---------------------------------------------------------------------
176!  Call to phyredem
177!---------------------------------------------------------------------
178      logical :: ok_writedem =.true.
179     
180!---------------------------------------------------------------------
181!  Call to physiq
182!---------------------------------------------------------------------
183      integer, parameter :: longcles=20
184      logical :: firstcall=.true.
185      logical :: lastcall=.false.
186      real :: phis    = 0.0
187      real :: clesphy0(longcles) = 0.0
188      real :: dpsrf
189
190!---------------------------------------------------------------------
191!  Initializations of boundary conditions
192!---------------------------------------------------------------------
193      integer, parameter :: yd = 360
194      real :: phy_nat (yd) = 0.0 ! 0=ocean libre,1=land,2=glacier,3=banquise
195      real :: phy_alb (yd)  ! Albedo land only (old value condsurf_jyg=0.3)
196      real :: phy_sst (yd)  ! SST (will not be used; cf read_tsurf1d.F)
197      real :: phy_bil (yd) = 1.0 ! Ne sert que pour les slab_ocean
198      real :: phy_rug (yd) ! Longueur rugosite utilisee sur land only
199      real :: phy_ice (yd) = 0.0 ! Fraction de glace
200      real :: phy_fter(yd) = 0.0 ! Fraction de terre
201      real :: phy_foce(yd) = 0.0 ! Fraction de ocean
202      real :: phy_fsic(yd) = 0.0 ! Fraction de glace
203      real :: phy_flic(yd) = 0.0 ! Fraction de glace
204
205!---------------------------------------------------------------------
206!  Fichiers et d'autres variables
207!---------------------------------------------------------------------
208      integer :: 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!---------------------------------------------------------------------
232!Al1
233        call conf_unicol
234!Al1 moves this gcssold var from common fcg_gcssold to
235        Turb_fcg_gcssold = xTurb_fcg_gcssold
236! --------------------------------------------------------------------
237        close(1)
238!Al1
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 )
331!-----------------------------------------------------------------------
332!   Choix du calendrier
333!   -------------------
334
335!      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
352!-----------------------------------------------------------------------
353!
354!c Date :
355!      La date est supposee donnee sous la forme [annee, numero du jour dans
356!      l annee] ; l heure est donnee dans time_ini, lu dans lmdz1d.def.
357!      On appelle ymds2ju pour convertir [annee, jour] en [jour Julien].
358!      Le numero du jour est dans "day". L heure est traitee separement.
359!      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
366! Special case for arm_cu which lasts less than one day : 53100s !! (MPL 20111026)
367      IF(forcing_type .EQ. 61) fnday=53100./86400.
368! 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 = int(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!     allocate(d_th_adv(llm))
436
437!
438!   No ozone climatology need be read in this pre-initialization
439!          (phys_state_var_init is called again in physiq)
440      read_climoz = 0
441!
442      call phys_state_var_init(read_climoz)
443
444      if (ngrid.ne.klon) then
445         print*,'stop in inifis'
446         print*,'Probleme de dimensions :'
447         print*,'ngrid = ',ngrid
448         print*,'klon  = ',klon
449         stop
450      endif
451!!!=====================================================================
452!!! Feedback forcing values for Gateaux differentiation (al1)
453!!!=====================================================================
454!!! Surface Planck forcing bracketing call radiation
455!!      surf_Planck = 0.
456!!      surf_Conv   = 0.
457!!      write(*,*) 'Gateaux-dif Planck,Conv:',surf_Planck,surf_Conv
458!!! a mettre dans le lmdz1d.def ou autre
459!!
460!!
461      qsol = qsolinp
462      qsurf = fq_sat(tsurf,psurf/100.)
463      rlat=xlat
464      rlon=xlon
465      day1= day_ini
466      time=daytime-day
467      ts_toga(1)=tsurf ! needed by read_tsurf1d.F
468      rho(1)=psurf/(rd*tsurf*(1.+(rv/rd-1.)*qsurf))
469
470!
471!! mpl et jyg le 22/08/2012 :
472!!  pour que les cas a flux de surface imposes marchent
473      IF(.NOT.ok_flux_surf.or.max(abs(wtsurf),abs(wqsurf))>0.) THEN
474       fsens=-wtsurf*rcpd*rho(1)
475       flat=-wqsurf*rlvtt*rho(1)
476       print *,'Flux: ok_flux wtsurf wqsurf',ok_flux_surf,wtsurf,wqsurf
477      ENDIF
478      print*,'Flux sol ',fsens,flat
479!!      ok_flux_surf=.false.
480!!      fsens=-wtsurf*rcpd*rho(1)
481!!      flat=-wqsurf*rlvtt*rho(1)
482!!!!
483
484! Vertical discretization and pressure levels at half and mid levels:
485
486      pa   = 5e4
487!!      preff= 1.01325e5
488      preff = psurf
489      IF (ok_old_disvert) THEN
490        call disvert0(pa,preff,ap,bp,dpres,presnivs,nivsigs,nivsig)
491        print *,'On utilise disvert0'
492      ELSE
493        call disvert(pa,preff,ap,bp,dpres,presnivs,nivsigs,nivsig,          &
494     &                 scaleheight)
495        print *,'On utilise disvert'
496!       Nouvelle version disvert permettant d imposer ap,bp (modif L.Guez) MPL 18092012
497!       Dans ce cas, on lit ap,bp dans le fichier hybrid.txt
498      ENDIF
499      sig_s=presnivs/preff
500      plev =ap+bp*psurf
501      play = 0.5*(plev(1:llm)+plev(2:llm+1))
502!cc      zlay=-rd*300.*log(play/psurf)/rg ! moved after reading profiles
503
504      IF (forcing_type .eq. 59) THEN
505! pour forcing_sandu, on cherche l'indice le plus proche de 700hpa#3000m
506      write(*,*) '***********************'
507      do l = 1, llm
508       write(*,*) 'l,play(l),presnivs(l): ',l,play(l),presnivs(l)
509       if (trouve_700 .and. play(l).le.70000) then
510         llm700=l
511         print *,'llm700,play=',llm700,play(l)/100.
512         trouve_700= .false.
513       endif
514      enddo
515      write(*,*) '***********************'
516      ENDIF
517
518!
519!=====================================================================
520! EVENTUALLY, READ FORCING DATA :
521!=====================================================================
522
523#include "1D_read_forc_cases.h"
524
525      if (forcing_GCM2SCM) then
526        write (*,*) 'forcing_GCM2SCM not yet implemented'
527        stop 'in initialization'
528      endif ! forcing_GCM2SCM
529
530      print*,'mxcalc=',mxcalc
531      print*,'zlay=',zlay(mxcalc)
532      print*,'play=',play(mxcalc)
533
534!Al1 pour SST forced, appellé depuis ocean_forced_noice
535      ts_cur = tsurf ! SST used in read_tsurf1d
536!=====================================================================
537! Initialisation de la physique :
538!=====================================================================
539
540!  Rq: conf_phys.F90 lit tous les flags de physiq.def; conf_phys appele depuis physiq.F
541!
542! day_step, iphysiq lus dans gcm.def ci-dessus
543! timestep: calcule ci-dessous from rday et day_step
544! ngrid=1
545! llm: defini dans .../modipsl/modeles/LMDZ4/libf/grid/dimension
546! rday: defini dans suphel.F (86400.)
547! day_ini: lu dans run.def (dayref)
548! rlat_rad,rlon-rad: transformes en radian de rlat,rlon lus dans lmdz1d.def (en degres)
549! airefi,zcufi,zcvfi initialises au debut de ce programme
550! rday,ra,rg,rd,rcpd declares dans YOMCST.h et calcules dans suphel.F
551      day_step = float(nsplit_phys)*day_step/float(iphysiq)
552      write (*,*) 'Time step divided by nsplit_phys (=',nsplit_phys,')'
553      timestep =rday/day_step
554      dtime_frcg = timestep
555!
556      zcufi=airefi
557      zcvfi=airefi
558!
559      rlat_rad(:)=rlat(:)*rpi/180.
560      rlon_rad(:)=rlon(:)*rpi/180.
561
562      call iniphysiq(ngrid,llm,rday,day_ini,timestep,                        &
563     &     rlat_rad,rlon_rad,airefi,zcufi,zcvfi,ra,rg,rd,rcpd,(/1/))
564      print*,'apres iniphysiq'
565
566! 2 PARAMETRES QUI DEVRAIENT ETRE LUS DANS run.def MAIS NE LE SONT PAS ICI:
567      co2_ppm= 330.0
568      solaire=1370.0
569
570! Ecriture du startphy avant le premier appel a la physique.
571! On le met juste avant pour avoir acces a tous les champs
572! NB: les clesphy0 seront remplies dans phyredem d'apres les flags lus dans gcm.def
573
574      if (ok_writedem) then
575
576!--------------------------------------------------------------------------
577! pbl_surface_init (called here) and pbl_surface_final (called by phyredem)
578! need : qsol fder snow qsurf evap rugos agesno ftsoil
579!--------------------------------------------------------------------------
580
581        type_ocean = "force"
582        run_off_lic_0(1) = restart_runoff
583        call fonte_neige_init(run_off_lic_0)
584
585        fder=0.
586        snsrf(1,:)=0.        ! couverture de neige des sous surface
587        qsurfsrf(1,:)=qsurf ! humidite de l'air des sous surface
588        evap=0.
589        frugs(1,:)=rugos     ! couverture de neige des sous surface
590        agesno  = xagesno
591        tsoil(:,:,:)=tsurf
592!------ AMMA 2e run avec modele sol et rayonnement actif (MPL 23052012)
593!       tsoil(1,1,1)=299.18
594!       tsoil(1,2,1)=300.08
595!       tsoil(1,3,1)=301.88
596!       tsoil(1,4,1)=305.48
597!       tsoil(1,5,1)=308.00
598!       tsoil(1,6,1)=308.00
599!       tsoil(1,7,1)=308.00
600!       tsoil(1,8,1)=308.00
601!       tsoil(1,9,1)=308.00
602!       tsoil(1,10,1)=308.00
603!       tsoil(1,11,1)=308.00
604!-----------------------------------------------------------------------
605        call pbl_surface_init(qsol, fder, snsrf, qsurfsrf,                  &
606     &                                    evap, frugs, agesno, tsoil)
607
608!------------------ prepare limit conditions for limit.nc -----------------
609!--   Ocean force
610
611        print*,'avant phyredem'
612        pctsrf(1,:)=0.
613        if (nat_surf.eq.0.) then
614          pctsrf(1,is_oce)=1.
615          pctsrf(1,is_ter)=0.
616        else
617          pctsrf(1,is_oce)=0.
618          pctsrf(1,is_ter)=1.
619        end if
620
621        print*,'nat_surf,pctsrf(1,is_oce),pctsrf(1,is_ter)',nat_surf         &
622     &        ,pctsrf(1,is_oce),pctsrf(1,is_ter)
623
624        zmasq=pctsrf(1,is_ter)+pctsrf(1,is_lic)
625        zpic = zpicinp
626        ftsol=tsurf
627        falb1 = albedo                           
628        falb2 = albedo                           
629        rugoro=rugos
630        t_ancien(1,:)=temp(:)
631        q_ancien(1,:)=q(:,1)
632        pbl_tke=1.e-8
633
634        rain_fall=0.
635        snow_fall=0.
636        solsw=0.
637        sollw=0.
638        radsol=0.
639        rnebcon=0.
640        ratqs=0.
641        clwcon=0.
642        zmea=0.
643        zstd=0.
644        zsig=0.
645        zgam=0.
646        zval=0.
647        zthe=0.
648        sig1=0.
649        w01=0.
650        u_ancien(1,:)=u(:)
651        v_ancien(1,:)=v(:)
652 
653!------------------------------------------------------------------------
654! Make file containing restart for the physics (startphy.nc)
655!
656! NB: List of the variables to be written by phyredem (via put_field):
657! rlon,rlat,zmasq,pctsrf(:,is_ter),pctsrf(:,is_lic),pctsrf(:,is_oce)
658! pctsrf(:,is_sic),ftsol(:,nsrf),tsoil(:,isoil,nsrf),qsurf(:,nsrf)
659! qsol,falb1(:,nsrf),falb2(:,nsrf),evap(:,nsrf),snow(:,nsrf)
660! radsol,solsw,sollw,fder,rain_fall,snow_fall,frugs(:,nsrf)
661! agesno(:,nsrf),zmea,zstd,zsig,zgam,zthe,zpic,zval,rugoro
662! t_ancien,q_ancien,frugs(:,is_oce),clwcon(:,1),rnebcon(:,1),ratqs(:,1)
663! run_off_lic_0,pbl_tke(:,1:klev,nsrf),zmax0,f0,sig1,w01
664! wake_deltat,wake_deltaq,wake_s,wake_cstar,wake_fip
665!------------------------------------------------------------------------
666!Al1 =============== restart option ==========================
667        if (.not.restart) then
668          call phyredem ("startphy.nc")
669        else
670! (desallocations)
671        print*,'callin surf final'
672          call pbl_surface_final(qsol, fder, snsrf, qsurfsrf,               &
673     &                                    evap, frugs, agesno, tsoil)
674        print*,'after surf final'
675          CALL fonte_neige_final(run_off_lic_0)
676        endif
677
678        ok_writedem=.false.
679        print*,'apres phyredem'
680
681      endif ! ok_writedem
682     
683!------------------------------------------------------------------------
684! Make file containing boundary conditions (limit.nc) **Al1->restartdyn***
685! --------------------------------------------------
686! NB: List of the variables to be written in limit.nc
687!     (by writelim.F, subroutine of 1DUTILS.h):
688!        phy_nat,phy_alb,phy_sst,phy_bil,phy_rug,phy_ice,
689!        phy_fter,phy_foce,phy_flic,phy_fsic)
690!------------------------------------------------------------------------
691      do i=1,yd
692        phy_nat(i)  = nat_surf
693        phy_alb(i)  = albedo
694        phy_sst(i)  = tsurf ! read_tsurf1d will be used instead
695        phy_rug(i)  = rugos
696        phy_fter(i) = pctsrf(1,is_ter)
697        phy_foce(i) = pctsrf(1,is_oce)
698        phy_fsic(i) = pctsrf(1,is_sic)
699        phy_flic(i) = pctsrf(1,is_lic)
700      enddo
701
702! fabrication de limit.nc
703      call writelim (1,phy_nat,phy_alb,phy_sst,phy_bil,phy_rug,             &
704     &               phy_ice,phy_fter,phy_foce,phy_flic,phy_fsic)
705
706
707      call phys_state_var_end
708!Al1
709      if (restart) then
710        print*,'call to restart dyn 1d'
711        Call dyn1deta0("start1dyn.nc",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
716!**      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
732! 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)=fcoriolis*(v(1:mxcalc)-vg(1:mxcalc))
850       dv_age(1:mxcalc)=-fcoriolis*(u(1:mxcalc)-ug(1:mxcalc))
851
852!!!!!!!!!!!!!!!!!!!!!!!!
853! Geostrophic wind
854!!!!!!!!!!!!!!!!!!!!!!!!
855       sfdt = sin(0.5*fcoriolis*timestep)
856       cfdt = cos(0.5*fcoriolis*timestep)
857!
858        du_age(1:mxcalc)= -2.*sfdt/timestep*                                &
859     &          (sfdt*(u(1:mxcalc)-ug(1:mxcalc)) -                          &
860     &           cfdt*(v(1:mxcalc)-vg(1:mxcalc))  )
861!!     : fcoriolis*(v(1:mxcalc)-vg(1:mxcalc))
862!
863       dv_age(1:mxcalc)= -2.*sfdt/timestep*                                 &
864     &          (cfdt*(u(1:mxcalc)-ug(1:mxcalc)) +                           &
865     &           sfdt*(v(1:mxcalc)-vg(1:mxcalc))  )
866!!     : -fcoriolis*(u(1:mxcalc)-ug(1:mxcalc))
867!
868!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
869!         call  writefield_phy('dv_age' ,dv_age,llm)
870!         call  writefield_phy('du_age' ,du_age,llm)
871!         call  writefield_phy('du_phys' ,du_phys,llm)
872!         call  writefield_phy('u_tend' ,u,llm)
873!         call  writefield_phy('u_g' ,ug,llm)
874!
875!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
876!! Increment state variables
877!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
878! pour les cas sandu et astex, on reclacule u,v,q,temp et teta dans 1D_nudge_sandu_astex.h
879! au dessus de 700hpa, on relaxe vers les profils initiaux
880      if (forcing_sandu .OR. forcing_astex) then
881#include "1D_nudge_sandu_astex.h"
882      else
883        u(1:mxcalc)=u(1:mxcalc) + timestep*(                                &
884     &              du_phys(1:mxcalc)                                       &
885     &             +du_age(1:mxcalc) )           
886        v(1:mxcalc)=v(1:mxcalc) + timestep*(                                 &
887     &              dv_phys(1:mxcalc)                                       &
888     &             +dv_age(1:mxcalc) )
889        q(1:mxcalc,:)=q(1:mxcalc,:)+timestep*(                              &
890     &                dq(1:mxcalc,:)                                        &
891     &               +d_q_adv(1:mxcalc,:) )
892
893        if (prt_level.ge.1) then
894          print *,                                                          &
895     &    'physiq-> temp(1),dt_phys(1),d_th_adv(1),dt_cooling(1) ',         &
896     &              temp(1),dt_phys(1),d_th_adv(1),dt_cooling(1)
897           print*,du_phys
898           print*, v
899           print*, vg
900        endif
901
902        temp(1:mxcalc)=temp(1:mxcalc)+timestep*(                            &
903     &              dt_phys(1:mxcalc)                                       &
904     &             +d_th_adv(1:mxcalc)                                      &
905     &             +dt_cooling(1:mxcalc))  ! Taux de chauffage ou refroid.
906
907      endif  ! forcing_sandu or forcing_astex
908
909        teta=temp*(pzero/play)**rkappa
910!
911!---------------------------------------------------------------------
912!   Nudge soil temperature if requested
913!---------------------------------------------------------------------
914
915      IF (nudge_tsoil) THEN
916       ftsoil(1,isoil_nudge,:) = ftsoil(1,isoil_nudge,:)                     &
917     &  -timestep/tau_soil_nudge*(ftsoil(1,isoil_nudge,:)-Tsoil_nudge)
918      ENDIF
919
920!---------------------------------------------------------------------
921!   Add large-scale tendencies (advection, etc) :
922!---------------------------------------------------------------------
923
924!cc nrlmd
925!cc        tmpvar=teta
926!cc        call advect_vert(llm,omega,timestep,tmpvar,plev)
927!cc
928!cc        teta(1:mxcalc)=tmpvar(1:mxcalc)
929!cc        tmpvar(:)=q(:,1)
930!cc        call advect_vert(llm,omega,timestep,tmpvar,plev)
931!cc        q(1:mxcalc,1)=tmpvar(1:mxcalc)
932!cc        tmpvar(:)=q(:,2)
933!cc        call advect_vert(llm,omega,timestep,tmpvar,plev)
934!cc        q(1:mxcalc,2)=tmpvar(1:mxcalc)
935
936!---------------------------------------------------------------------
937!   Air temperature :
938!---------------------------------------------------------------------       
939        if (lastcall) then
940          print*,'Pas de temps final ',it
941          call ju2ymds(daytime, an, mois, jour, heure)
942          print*,'a la date : a m j h',an, mois, jour ,heure/3600.
943        endif
944
945!  incremente day time
946!        print*,'daytime bef',daytime,1./day_step
947        daytime = daytime+1./day_step
948!Al1dbg
949        day = int(daytime+0.1/day_step)
950!        time = max(daytime-day,0.0)
951!Al1&jyg: correction de bug
952!cc        time = real(mod(it,day_step))/day_step
953        time = time_ini/24.+real(mod(it,day_step))/day_step
954!        print*,'daytime nxt time',daytime,time
955        it=it+1
956
957      enddo
958
959!Al1
960      if (ecrit_slab_oc.ne.-1) close(97)
961
962!Al1 Call to 1D equivalent of dynredem (an,mois,jour,heure ?)
963! -------------------------------------
964       call dyn1dredem("restart1dyn.nc",                                    &
965     &              plev,play,phi,phis,presnivs,                            &
966     &              u,v,temp,q,omega2)
967
968        CALL abort_gcm ('lmdz1d   ','The End  ',0)
969
970      end
971
972#include "1DUTILS.h"
973#include "1Dconv.h"
974
975#endif
976
Note: See TracBrowser for help on using the repository browser.