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

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

Make iniphysiq a module.
Fix call to iniphysiq in lmdz1d (missing arguments and arrays of wrong sizes).
EM

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