source: LMDZ5/trunk/libf/phylmd/dyn1d/lmdz1d.F90 @ 2351

Last change on this file since 2351 was 2351, checked in by Ehouarn Millour, 9 years ago

More on physics/dynamics separation and cleanup:

  • Set things up so that all physics-related initializations are done via iniphysiq.
  • Created a "geometry_mod.F90" module in phy_common to store information on the loacl grid (i.e. replaces comgeomphy) and moreover give these variables more obvious names (e.g.: rlond => longitude, rlatd => latitude, airephy => cell_area).
  • removed obsolete comgeomphy.h and comgeomphy.F90

EM

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