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

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

Some cleanup: remove (unused) clesph0 from dynamics.
EM

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