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

Last change on this file since 2186 was 2181, checked in by jyg, 10 years ago

1/ Update of the loop on sub-surfaces in pbl_surface_mod.F90 : all initializations
are removed from the loop and put before the loop.
2/ The possibility of nudging RH, T, u, and v is added in lmdz1d.F90 (and in
1DUTILS.h).

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