source: trunk/LMDZ.COMMON/libf/evolution/init_phys_1d_mod.F90 @ 2998

Last change on this file since 2998 was 2980, checked in by romain.vande, 3 years ago

Mars PEM :

Adapt PEM to 1d runs.
Cleaning of names and unused variables.
Correct minor errors.
Adapt and correct reshape_xios_output utilitary for 1d diagfi output.

RV

File size: 20.4 KB
Line 
1module init_phys_1d_mod
2
3implicit none
4
5contains
6
7#ifdef CPP_1D
8
9subroutine init_phys_1d(llm_in,nqtot_in,v,u,temp,q,psurf,time,ap_out,bp_out)
10
11!-------------- init_phys_1d -------------
12!
13! This function is a copy of the test_phys_1d.F90 initialisation part.
14! It is meant to initialize all the variables in modules and for the ouput
15! that can't be initialzed via the startfi file.
16!
17! Basically it reads run.def and save variables values
18!
19!--------------              -------------
20
21      USE ioipsl_getincom, only: getin
22      use comcstfi_h, only: pi, rad, omeg, g, mugaz, rcp, r, cpp
23      use time_phylmdz_mod, only: daysec, dtphys, day_step, &
24                                 ecritphy, iphysiq
25      use planete_h, only: year_day, periheli, aphelie, peri_day,&
26                          obliquit, emin_turb, lmixmin
27      use surfdat_h, only: albedodat, z0_default, emissiv, emisice,&
28                          albedice, iceradius, dtemisice, z0,&
29                          zmea, zstd, zsig, zgam, zthe, phisfi,&
30                          watercaptag, watercap, hmons, summit, base,&
31                           tsurf, emis,qsurf       
32      use infotrac, only: nqtot, tname, nqperes,nqfils
33      use regular_lonlat_mod, only: init_regular_lonlat
34      use mod_grid_phy_lmdz, only : regular_lonlat
35      USE comvert_mod, ONLY: ap,bp,aps,bps,pa,preff,sig,&
36                            presnivs,pseudoalt,scaleheight
37      use dimradmars_mod, only: tauvis,totcloudfrac
38      use mod_interface_dyn_phys, only: init_interface_dyn_phys
39      use geometry_mod, only: init_geometry
40      use dimphy, only : init_dimphy
41      USE phys_state_var_init_mod, ONLY: phys_state_var_init
42      use comgeomfi_h, only: sinlat, ini_fillgeom
43      use slope_mod, only: theta_sl, psi_sl
44      use comslope_mod, only: def_slope,subslope_dist,def_slope_mean
45      use dimradmars_mod, only: tauvis,totcloudfrac
46      use dust_param_mod, only: tauscaling
47      use tracer_mod, only: igcm_h2o_vap,igcm_h2o_ice,igcm_co2,noms
48      USE logic_mod, ONLY: hybrid
49      USE vertical_layers_mod, ONLY: init_vertical_layers
50      use comsoil_h, only: volcapa, layer, mlayer, inertiedat, &
51                          inertiesoil,nsoilmx,flux_geo
52      USE read_profile_mod, ONLY: read_profile
53      use inichim_newstart_mod, only: inichim_newstart
54      use physics_distribution_mod, only: init_physics_distribution
55      use iostart, only: open_startphy,get_var, close_startphy
56
57      include "dimensions.h"
58      integer, parameter :: ngrid = 1 !(2+(jjm-1)*iim - 1/jjm)
59      integer, parameter :: nlayer = llm
60      include "callkeys.h"
61      include "netcdf.inc"
62
63!-------------- INPUT VARIABLES -------------
64
65      integer, intent(in) :: llm_in,nqtot_in
66
67!-------------- OUTPUT VARIABLES -------------
68
69      REAL,INTENT(OUT):: u(nlayer),v(nlayer)  ! zonal, meridional wind
70      REAL,INTENT(OUT):: temp(nlayer)   ! temperature at the middle of the layers
71      REAL,INTENT(OUT) :: q(2,llm_in,nqtot_in) ! tracer mixing ratio (e.g. kg/kg)
72      REAL,INTENT(OUT):: psurf(2)
73      REAL,INTENT(OUT):: time             ! time (0<time<1 ; time=0.5 a midi)
74      REAL,INTENT(OUT) :: ap_out(llm_in+1),bp_out(llm_in+1)
75
76
77!-------------- LOCAL VARIABLES -------------
78
79      REAL halfaxe, excentric, Lsperi
80      integer :: nq=1 ! number of tracers
81      real :: latitude(1), longitude(1), cell_area(1)
82
83!   MVals: isotopes as in the dynamics (CRisi)
84      INTEGER :: ifils,ipere,generation
85      CHARACTER(len=30), ALLOCATABLE, DIMENSION(:) :: tnom_transp ! transporting fluid short name
86      CHARACTER(len=80) :: line ! to store a line of text     
87      INTEGER ierr0
88      LOGICAL :: continu
89
90      logical :: present,found
91      INTEGER nlevel,nsoil,ndt
92      REAL gru,grv   ! prescribed "geostrophic" background wind
93      REAL pks, ptif, w(nlayer)
94      REAL q2(nlayer+1)   ! Turbulent Kinetic Energy
95      REAL play(nlayer)   ! Pressure at the middle of the layers (Pa)
96      REAL plev(nlayer+1) ! intermediate pressure levels (pa)
97      REAL zlay(nlayer)   ! altitude estimee dans les couches (km)
98      REAL tmp1(0:nlayer),tmp2(0:nlayer)
99      INTEGER flagthermo,flagh2o
100      REAL atm_wat_profile
101      REAL,ALLOCATABLE :: zqsat(:) ! useful for (atm_wat_profile=2)
102
103      INTEGER :: ierr,iq,ilayer,ilevel,isoil
104      INTEGER day0,dayn          ! date initial (sol ; =0 a Ls=0) and final
105      REAL day           ! date durant le run
106      REAL tab_cntrl(100)
107      real,parameter :: odpref=610. ! DOD reference pressure (Pa)
108
109!   LL: Subsurface geothermal flux
110      real :: flux_geo_tmp
111
112      real,allocatable :: qdyn(:,:,:,:),psdyn(:,:)
113
114
115! check if 'run.def' file is around (otherwise reading parameters
116! from callphys.def via getin() routine won't work.
117      open(99,file='run.def',status='old',form='formatted',&
118          iostat=ierr)
119      if (ierr.ne.0) then
120        write(*,*) 'Cannot find required file "run.def"'
121        write(*,*) '  (which should contain some input parameters'
122        write(*,*) '   along with the following line:'
123        write(*,*) '   INCLUDEDEF=callphys.def'
124        write(*,*) '   )'
125        write(*,*) ' ... might as well stop here ...'
126        stop
127      else
128        close(99)
129      endif
130
131! ------------------------------------------------------
132!  Prescribed constants to be set here
133! ------------------------------------------------------
134
135!      if(.not. startfile_1D ) then
136
137      pi=2.E+0*asin(1.E+0)
138
139!     Mars planetary constants
140!     ----------------------------
141      rad=3397200.               ! mars radius (m)  ~3397200 m
142      daysec=88775.              ! length of a sol (s)  ~88775 s
143      omeg=4.*asin(1.)/(daysec)  ! rotation rate (rad.s-1)
144      g=3.72                     ! gravity (m.s-2) ~3.72 
145      mugaz=43.49                ! atmosphere mola mass (g.mol-1) ~43.49
146      rcp=.256793                ! = r/cp  ~0.256793
147      r= 8.314511E+0 *1000.E+0/mugaz
148      cpp= r/rcp
149      year_day = 669             ! lenght of year (sols) ~668.6
150      periheli = 206.66          ! minimum sun-mars distance (Mkm) ~206.66
151      aphelie = 249.22           ! maximum sun-mars distance (Mkm) ~249.22
152      halfaxe = 227.94           ! demi-grand axe de l'ellipse
153      peri_day =  485.           ! perihelion date (sols since N. Spring)
154      obliquit = 25.2            ! Obliquity (deg) ~25.2         
155      excentric = 0.0934         ! Eccentricity (0.0934)         
156 
157!     Planetary Boundary Layer and Turbulence parameters
158!     --------------------------------------------------
159      z0_default =  1.e-2        ! surface roughness (m) ~0.01
160      emin_turb = 1.e-6          ! minimal turbulent energy ~1.e-8
161      lmixmin = 30               ! mixing length ~100
162 
163!     cap properties and surface emissivities
164!     ----------------------------------------------------
165      emissiv= 0.95              ! Bare ground emissivity ~.95
166      emisice(1)=0.95            ! Northern cap emissivity
167      emisice(2)=0.95            ! Southern cap emisssivity
168      albedice(1)=0.5            ! Northern cap albedo
169      albedice(2)=0.5            ! Southern cap albedo
170      iceradius(1) = 100.e-6     ! mean scat radius of CO2 snow (north)
171      iceradius(2) = 100.e-6     ! mean scat radius of CO2 snow (south)
172      dtemisice(1) = 2.          ! time scale for snow metamorphism (north)
173      dtemisice(2) = 2.          ! time scale for snow metamorphism (south)
174
175!      endif ! .not. startfile_1D
176
177!     mesh surface (not a very usefull quantity in 1D)
178!     ----------------------------------------------------
179      cell_area(1)=1.E+0
180
181
182! check if there is a 'traceur.def' file
183! and process it.
184      ! load tracer names from file 'traceur.def'
185        open(90,file='traceur.def',status='old',form='formatted',&
186            iostat=ierr)
187        if (ierr.ne.0) then
188          write(*,*) 'Cannot find required file "traceur.def"'
189          write(*,*) ' If you want to run with tracers, I need it'
190          write(*,*) ' ... might as well stop here ...'
191          stop
192        else
193          nqtot=nqtot_in ! set value of nqtot (in infotrac module) as nq
194          nq=nqtot_in
195        endif
196        ! allocate arrays:
197        allocate(tname(nq))
198        allocate(qsurf(1,1,nq))
199        allocate(tnom_transp(nq))
200       
201        ! read tracer names from file traceur.def
202        do iq=1,nq
203          read(90,'(80a)',iostat=ierr) line ! store the line from traceur.def
204          if (ierr.ne.0) then
205            write(*,*) 'testphys1d: error reading tracer names...'
206            stop
207          endif
208          ! if format is tnom_0, tnom_transp (isotopes)
209          read(line,*,iostat=ierr0) tname(iq),tnom_transp(iq)
210          if (ierr0.ne.0) then
211            read(line,*) tname(iq)
212            tnom_transp(iq)='air'
213          endif
214
215        enddo
216        close(90)
217
218       ! Isotopes: as in the 3D case we have to determine father/son relations for isotopes and carrying fluid
219       ALLOCATE(nqfils(nqtot))
220       nqperes=0
221       nqfils(:)=0 
222       DO iq=1,nqtot
223       if (tnom_transp(iq) == 'air') then
224         ! ceci est un traceur père
225         WRITE(*,*) 'Le traceur',iq,', appele ',&
226               trim(tname(iq)),', est un pere'
227         nqperes=nqperes+1
228       else !if (tnom_transp(iq) == 'air') then
229         ! ceci est un fils. Qui est son père?
230         WRITE(*,*) 'Le traceur',iq,', appele ',&
231                     trim(tname(iq)),', est un fils'
232         continu=.true.
233         ipere=1
234         do while (continu)           
235           if (tnom_transp(iq) .eq. tname(ipere)) then
236             ! Son père est ipere
237             WRITE(*,*) 'Le traceur',iq,'appele ',&
238        trim(tname(iq)),' est le fils de ',&
239        ipere,'appele ',trim(tname(ipere))
240             nqfils(ipere)=nqfils(ipere)+1         
241             continu=.false.
242           else !if (tnom_transp(iq) == tnom_0(ipere)) then
243             ipere=ipere+1
244             if (ipere.gt.nqtot) then
245                 WRITE(*,*) 'Le traceur',iq,'appele ',&
246                trim(tname(iq)),', est orpelin.'
247                 CALL abort_gcm('infotrac_init',&
248                       'Un traceur est orphelin',1)
249             endif !if (ipere.gt.nqtot) then
250           endif !if (tnom_transp(iq) == tnom_0(ipere)) then
251         enddo !do while (continu)
252       endif !if (tnom_transp(iq) == 'air') then
253       enddo !DO iq=1,nqtot
254       WRITE(*,*) 'nqperes=',nqperes   
255       WRITE(*,*) 'nqfils=',nqfils
256
257        ! initialize tracers here:
258        write(*,*) "testphys1d: initializing tracers"
259        call read_profile(nq, nlayer, qsurf, q(1,:,:))
260        q(2,:,:)=q(1,:,:)
261
262      call init_physics_distribution(regular_lonlat,4,&
263                                    1,1,1,nlayer,1)
264
265      call open_startphy("startfi_evol.nc")
266      call get_var("controle",tab_cntrl,found)
267       if (.not.found) then
268         call abort_physic("open_startphy", &
269             "tabfi: Failed reading <controle> array",1)
270       else
271         write(*,*)'tabfi: tab_cntrl',tab_cntrl
272       endif
273       day0 = tab_cntrl(3)
274       day=float(day0)
275
276       call get_var("Time",time,found)
277
278       call close_startphy
279
280
281!  Discretization (Definition of grid and time steps)
282!  --------------
283!
284      nlevel=nlayer+1
285      nsoil=nsoilmx
286
287      day_step=48 ! default value for day_step
288      PRINT *,'Number of time steps per sol ?'
289      call getin("day_step",day_step)
290      write(*,*) " day_step = ",day_step
291
292      ecritphy=day_step ! default value for ecritphy, output every time step
293
294      ndt=10 ! default value for ndt
295      PRINT *,'Number of sols to run ?'
296      call getin("ndt",ndt)
297      write(*,*) " ndt = ",ndt
298
299      dayn=day0+ndt
300      ndt=ndt*day_step     
301      dtphys=daysec/day_step 
302
303! Imposed surface pressure
304! ------------------------------------
305!
306
307      psurf(1)=610. ! default value for psurf
308      PRINT *,'Surface pressure (Pa) ?'
309      call getin("psurf",psurf(1))
310      write(*,*) " psurf = ",psurf(1)
311      psurf(2)=psurf(1)
312
313! Reference pressures
314      pa=20.   ! transition pressure (for hybrid coord.)
315      preff=610.      ! reference surface pressure
316 
317! Aerosol properties
318! --------------------------------
319      tauvis=0.2 ! default value for tauvis (dust opacity)
320      write(*,'("Reference dust opacity at ",f4.0," Pa ?")')odpref
321      call getin("tauvis",tauvis)
322      write(*,*) " tauvis = ",tauvis
323
324! Orbital parameters
325! ------------------
326
327!  latitude/longitude
328!  ------------------
329      latitude(1)=0 ! default value for latitude
330      PRINT *,'latitude (in degrees) ?'
331      call getin("latitude",latitude(1))
332      write(*,*) " latitude = ",latitude
333      latitude=latitude*pi/180.E+0
334      longitude=0.E+0
335      longitude=longitude*pi/180.E+0
336
337!  some initializations (some of which have already been
338!  done above!) and loads parameters set in callphys.def
339!  and allocates some arrays
340!Mars possible matter with dtphys in input and include!!!
341! Initializations below should mimick what is done in iniphysiq for 3D GCM
342      call init_interface_dyn_phys
343      call init_regular_lonlat(1,1,longitude,latitude,&
344                        (/0.,0./),(/0.,0./))
345      call init_geometry(1,longitude,latitude,&
346                        (/0.,0.,0.,0./),(/0.,0.,0.,0./),&
347                        cell_area)
348! Ehouarn: init_vertial_layers called later (because disvert not called yet)
349!      call init_vertical_layers(nlayer,preff,scaleheight,
350!     &                      ap,bp,aps,bps,presnivs,pseudoalt)
351      call init_dimphy(1,nlayer) ! Initialize dimphy module
352      call phys_state_var_init(1,llm,nq,tname,&
353               day0,dayn,time,&
354               daysec,dtphys,&
355               rad,g,r,cpp,&
356               nqperes,nqfils)! MVals: variables isotopes
357      call ini_fillgeom(1,latitude,longitude,(/1.0/))
358      call conf_phys(1,llm,nq)
359
360      ! in 1D model physics are called every time step
361      ! ovverride iphysiq value that has been set by conf_phys
362      if (iphysiq/=1) then
363        write(*,*) "testphys1d: setting iphysiq=1"
364        iphysiq=1
365      endif
366
367!  Initialize albedo / soil thermal inertia
368!  ----------------------------------------
369!
370
371
372! Initialize local slope parameters (only matters if "callslope"
373! is .true. in callphys.def)
374      ! slope inclination angle (deg) 0: horizontal, 90: vertical
375      theta_sl(1)=0.0 ! default: no inclination
376      call getin("slope_inclination",theta_sl(1))
377      ! slope orientation (deg)
378      ! 0 == Northward, 90 == Eastward, 180 == Southward, 270 == Westward
379      psi_sl(1)=0.0 ! default value
380      call getin("slope_orientation",psi_sl(1))
381     
382      ! sub-slopes parameters (assuming no sub-slopes distribution for now).
383      def_slope(1)=-90 ! minimum slope angle
384      def_slope(2)=90 ! maximum slope angle
385      subslope_dist(1,1)=1 ! fraction of subslopes in mesh
386!
387!  for the gravity wave scheme
388!  ---------------------------------
389!
390      zmea(1)=0.E+0
391      zstd(1)=0.E+0
392      zsig(1)=0.E+0
393      zgam(1)=0.E+0
394      zthe(1)=0.E+0
395!
396!  for the slope wind scheme
397!  ---------------------------------
398
399      hmons(1)=0.E+0
400      PRINT *,'hmons is initialized to ',hmons(1)
401      summit(1)=0.E+0
402      PRINT *,'summit is initialized to ',summit(1)
403      base(1)=0.E+0
404!
405!  Default values initializing the coefficients calculated later
406!  ---------------------------------
407
408      tauscaling(1)=1. ! calculated in aeropacity_mod.F
409      totcloudfrac(1)=1. ! calculated in watercloud_mod.F     
410
411!   Specific initializations for "physiq"
412!   -------------------------------------
413!   surface geopotential is not used (or useful) since in 1D
414!   everything is controled by surface pressure
415      phisfi(1)=0.E+0
416
417!   Initialization to take into account prescribed winds
418!   ------------------------------------------------------
419      ptif=2.E+0*omeg*sinlat(1)
420 
421!    geostrophic wind
422      gru=10. ! default value for gru
423      PRINT *,'zonal eastward component of the geostrophic wind (m/s) ?'
424      call getin("u",gru)
425      write(*,*) " u = ",gru
426      grv=0. !default value for grv
427      PRINT *,'meridional northward component of the geostrophic',&
428     ' wind (m/s) ?'
429      call getin("v",grv)
430      write(*,*) " v = ",grv
431
432!     Initialize winds  for first time step
433      DO ilayer=1,nlayer
434         u(ilayer)=gru
435         v(ilayer)=grv
436         w(ilayer)=0 ! default: no vertical wind
437      ENDDO
438
439!     Initialize turbulente kinetic energy
440      DO ilevel=1,nlevel
441         q2(ilevel)=0.E+0
442      ENDDO
443
444!  CO2 ice on the surface
445!  -------------------
446      ! get the index of co2 tracer (not known at this stage)
447      igcm_co2=0
448      do iq=1,nq
449        if (trim(tname(iq))=="co2") then
450          igcm_co2=iq
451        endif
452      enddo
453      if (igcm_co2==0) then
454        write(*,*) "testphys1d error, missing co2 tracer!"
455        stop
456      endif
457
458!  Compute pressures and altitudes of atmospheric levels
459!  ----------------------------------------------------------------
460
461!    Vertical Coordinates
462!    """"""""""""""""""""
463      hybrid=.true.
464      PRINT *,'Hybrid coordinates ?'
465      call getin("hybrid",hybrid)
466      write(*,*) " hybrid = ", hybrid
467
468      CALL  disvert_noterre
469      ! now that disvert has been called, initialize module vertical_layers_mod
470      call init_vertical_layers(nlayer,preff,scaleheight,&
471                           ap,bp,aps,bps,presnivs,pseudoalt)
472
473      DO ilevel=1,nlevel
474        plev(ilevel)=ap(ilevel)+psurf(1)*bp(ilevel)
475      ENDDO
476
477      DO ilayer=1,nlayer
478        play(ilayer)=aps(ilayer)+psurf(1)*bps(ilayer)
479      ENDDO
480
481      DO ilayer=1,nlayer
482        zlay(ilayer)=-200.E+0 *r*log(play(ilayer)/plev(1))&
483        /g
484      ENDDO
485
486
487!  Initialize temperature profile
488!  --------------------------------------
489      pks=psurf(1)**rcp
490
491! altitude in km in profile: divide zlay by 1000
492      tmp1(0)=0.E+0
493      DO ilayer=1,nlayer
494        tmp1(ilayer)=zlay(ilayer)/1000.E+0
495      ENDDO
496
497      call profile(nlayer+1,tmp1,tmp2)
498
499      tsurf=tmp2(0)
500      DO ilayer=1,nlayer
501        temp(ilayer)=tmp2(ilayer)
502      ENDDO
503     
504
505! Initialize soil properties and temperature
506! ------------------------------------------
507      volcapa=1.e6 ! volumetric heat capacity
508
509      flux_geo_tmp=0.
510      call getin("flux_geo",flux_geo_tmp)
511      flux_geo(:,:) = flux_geo_tmp
512
513! Initialize depths
514! -----------------
515      do isoil=0,nsoil-1
516        mlayer(isoil)=2.e-4*(2.**(isoil-0.5)) ! mid-layer depth
517      enddo
518      do isoil=1,nsoil
519        layer(isoil)=2.e-4*(2.**(isoil-1)) ! layer depth
520      enddo
521
522!    Initialize traceurs
523!    ---------------------------
524
525      if (photochem.or.callthermos) then
526         write(*,*) 'Initializing chemical species'
527         ! flagthermo=0: initialize over all atmospheric layers
528         flagthermo=0
529         ! check if "h2o_vap" has already been initialized
530         ! (it has been if there is a "profile_h2o_vap" file around)
531         inquire(file="profile_h2o_vap",exist=present)
532         if (present) then
533           flagh2o=0 ! 0: do not initialize h2o_vap
534         else
535           flagh2o=1 ! 1: initialize h2o_vap in inichim_newstart
536         endif
537         
538         ! hack to accomodate that inichim_newstart assumes that
539         ! q & psurf arrays are on the dynamics scalar grid
540         allocate(qdyn(2,1,llm,nq),psdyn(2,1))
541         qdyn(1,1,1:llm,1:nq)=q(1,1:llm,1:nq)
542         psdyn(1:2,1)=psurf(1)
543         call inichim_newstart(ngrid,nq,qdyn,qsurf,psdyn,&
544                              flagh2o,flagthermo)
545         q(1,1:llm,1:nq)=qdyn(1,1,1:llm,1:nq)
546      endif
547
548! Check if the surface is a water ice reservoir
549! --------------------------------------------------
550      watercaptag(1)=.false. ! Default: no water ice reservoir
551      print *,'Water ice cap on ground ?'
552      call getin("watercaptag",watercaptag)
553      write(*,*) " watercaptag = ",watercaptag
554     
555! Check if the atmospheric water profile is specified
556! ---------------------------------------------------
557      ! Adding an option to force atmospheric water values JN
558      atm_wat_profile=-1 ! Default: free atm wat profile
559      print *,'Force atmospheric water vapor profile ?'
560      call getin("atm_wat_profile",atm_wat_profile)
561      write(*,*) "atm_wat_profile = ", atm_wat_profile
562      if (atm_wat_profile.EQ.-1) then
563        write(*,*) "Free atmospheric water vapor profile"
564        write(*,*) "Total water is conserved on the column"
565      else if (atm_wat_profile.EQ.0) then
566        write(*,*) "Dry atmospheric water vapor profile"
567      else if ((atm_wat_profile.GT.0).and.(atm_wat_profile.LE.1)) then
568        write(*,*) "Prescribed atmospheric water vapor MMR profile"
569        write(*,*) "Unless it reaches saturation (maximal value)"
570        write(*,*) "MMR chosen : ", atm_wat_profile
571      endif
572
573      ap_out=ap
574      bp_out=bp
575
576end subroutine init_phys_1d
577
578#endif
579
580end module init_phys_1d_mod
Note: See TracBrowser for help on using the repository browser.