source: trunk/LMDZ.MARS/libf/phymars/dyn1d/init_testphys1d_mod.F90 @ 4063

Last change on this file since 4063 was 4059, checked in by jbclement, 7 days ago

Mars PCM:
Deletion of 'qsurf' in the information provided by "start1D.txt" since it is unused and unnecessary until now. Addition of an error-leading check to prevent an unintended situation (which should not happen) regarding the presence of starting files.
JBC

File size: 34.0 KB
Line 
1MODULE init_testphys1d_mod
2
3implicit none
4
5contains
6
7SUBROUTINE init_testphys1d(start1Dname,startfiname,therestart1D,therestartfi,ngrid,nlayer,odpref,               &
8                           nq,q,time,psurf,u,v,temp,ndt,ptif,pks,dttestphys,zqsat,dq,dqdyn,day0,day,gru,grv,w,  &
9                           play,plev,latitude,longitude,cell_area,ctrl_h2ovap,relaxtime_h2ovap,qref_h2ovap,      &
10                           ctrl_h2oice,relaxtime_h2oice,qref_h2oice)
11
12use ioipsl_getincom,          only: getin ! To use 'getin'
13use comcstfi_h,               only: pi, rad, omeg, g, mugaz, rcp, r, cpp
14use time_phylmdz_mod,         only: daysec, steps_per_sol, outputs_per_sol, iphysiq
15use planete_h,                only: year_day, periheli, aphelie, peri_day, obliquit, emin_turb, lmixmin
16use surfdat_h,                only: albedodat, z0_default, z0, emissiv, emisice, albedice, iceradius, dtemisice,      &
17                                    zmea, zstd, zsig, zgam, zthe, hmons, summit, base, phisfi, watercaptag, watercap, &
18                                    tsurf, emis, qsurf, perennial_co2ice, ini_surfdat_h_slope_var, end_surfdat_h_slope_var
19use infotrac,                 only: nqtot, tname, nqperes, nqfils
20use read_profile_mod,         only: read_profile
21use profile_temp_mod,         only: profile_temp
22use iostart,                  only: open_startphy, get_var, close_startphy
23use physics_distribution_mod, only: init_physics_distribution
24use comsoil_h,                only: volcapa, nsoilmx, inertiesoil, inertiedat, layer, mlayer, flux_geo, tsoil, &
25                                    adsorption_soil, qsoil, igcm_h2o_ice_soil, porosity_reg,                   &
26                                    ini_comsoil_h_slope_var, end_comsoil_h_slope_var, ads_massive_ice
27use comvert_mod,              only: ap, bp, aps, bps, pa, preff, presnivs, pseudoalt, scaleheight
28use dimradmars_mod,           only: tauvis, totcloudfrac, albedo, ini_dimradmars_mod_slope_var, end_dimradmars_mod_slope_var
29use regular_lonlat_mod,       only: init_regular_lonlat
30use mod_interface_dyn_phys,   only: init_interface_dyn_phys
31use geometry_mod,             only: init_geometry, init_geometry_cell_area_for_outputs
32use dimphy,                   only: init_dimphy
33use comgeomfi_h,              only: sinlat, ini_fillgeom
34use slope_mod,                only: theta_sl, psi_sl
35use comslope_mod,             only: def_slope, subslope_dist
36use dust_param_mod,           only: tauscaling
37use tracer_mod,               only: igcm_co2
38use logic_mod,                only: hybrid
39use vertical_layers_mod,      only: init_vertical_layers
40use inichim_newstart_mod,     only: inichim_newstart
41use mod_grid_phy_lmdz,        only: regular_lonlat
42use phys_state_var_init_mod,  only: phys_state_var_init
43use turb_mod,                 only: q2
44use nonoro_gwd_ran_mod,       only: du_nonoro_gwd, dv_nonoro_gwd
45use conf_phys_mod,            only: conf_phys
46use paleoclimate_mod,         only: paleoclimate, ini_paleoclimate_h, end_paleoclimate_h, h2o_ice_depth
47use comslope_mod,             only: nslope, subslope_dist, ini_comslope_h, end_comslope_h
48use co2condens_mod,           only: CO2cond_ps
49use callkeys_mod,             only: water, photochem, callthermos
50use call_dayperi_mod,         only: call_dayperi
51use initracer_mod,            only: initracer
52! Mostly for XIOS outputs:
53use mod_const_mpi,            only: COMM_LMDZ
54
55implicit none
56
57include "dimensions.h"
58
59!=======================================================================
60! Arguments
61!=======================================================================
62integer,      intent(in) :: ngrid, nlayer
63real,         intent(in) :: odpref                     ! DOD reference pressure (Pa)
64character(*), intent(in) :: start1Dname, startfiname   ! Name of starting files for 1D
65logical,      intent(in) :: therestart1D, therestartfi ! Use of starting files for 1D
66
67integer, intent(inout) :: nq
68
69real, dimension(:,:,:), allocatable, intent(out) :: q          ! tracer mixing ratio (e.g. kg/kg)
70real,                                intent(out) :: time       ! time (0<time<1; time=0.5 at noon)
71real,                                intent(out) :: psurf      ! Surface pressure
72real, dimension(nlayer),             intent(out) :: u, v       ! zonal, meridional wind
73real, dimension(nlayer),             intent(out) :: temp       ! temperature at the middle of the layers
74integer,                             intent(out) :: ndt
75real,                                intent(out) :: ptif, pks
76real,                                intent(out) :: dttestphys ! testphys1d timestep
77real, dimension(:),     allocatable, intent(out) :: zqsat
78real, dimension(:,:,:), allocatable, intent(out) :: dq, dqdyn  ! Physical and dynamical tandencies
79integer,                             intent(out) :: day0       ! initial (sol ; =0 at Ls=0) and final date
80real,                                intent(out) :: day        ! date during the run
81real,                                intent(out) :: gru, grv   ! prescribed "geostrophic" background wind
82real, dimension(nlayer),             intent(out) :: w          ! "Dummy wind" in 1D
83real, dimension(nlayer),             intent(out) :: play       ! Pressure at the middle of the layers (Pa)
84real, dimension(nlayer + 1),         intent(out) :: plev       ! intermediate pressure levels (pa)
85real, dimension(1),                  intent(out) :: latitude, longitude, cell_area
86logical,                             intent(out) :: ctrl_h2ovap, ctrl_h2oice           ! Flags to prescribe atmospheric water vapour/ice profiles
87real,                                intent(out) :: relaxtime_h2ovap, relaxtime_h2oice ! Relaxation time towards the atmospheric water vapour/ice profiles
88real, dimension(nlayer),             intent(out) :: qref_h2ovap, qref_h2oice           ! User-defined atmospheric water vapour/ice profiles
89
90!=======================================================================
91! Local variables
92!=======================================================================
93integer                   :: ierr, iq, j, ilayer, isoil, nlevel, nsoil, flagthermo, flagh2o, nslope_tmp
94integer                   :: dayn ! Final date
95real, dimension(nlayer)   :: zlay ! altitude estimee dans les couches (km)
96real, dimension(0:nlayer) :: zkm, teta
97
98! Dummy variables along "dynamics scalar grid"
99real, dimension(:,:,:,:), allocatable :: qdyn
100real, dimension(:,:),     allocatable :: psdyn
101
102! Use of starting files for 1D
103logical                :: found
104character(30)          :: header
105real, dimension(100)   :: tab_cntrl
106
107! New flag to compute paleo orbital configurations + few variables JN
108logical :: paleomars
109real    :: halfaxe, eccentric, Lsperi
110
111! MVals: isotopes as in the dynamics (CRisi)
112integer                                  :: ipere, ierr0
113character(30), dimension(:), allocatable :: tnom_transp ! transporting fluid short name
114character(80)                            :: line        ! to store a line of text
115logical                                  :: continu
116
117! LL: Possibility to add subsurface ice
118real    :: ice_depth          ! depth of the ice table, ice_depth < 0. means no ice
119real    :: inertieice         ! ice thermal inertia
120integer :: iref
121
122! LL: Subsurface water ice
123real :: rho_H2O_ice           ! Density (kg/m^3) of water ice
124
125!=======================================================================
126! Code
127!=======================================================================
128!------------------------------------------------------
129! Prescribed constants to be set here
130!------------------------------------------------------
131! Mars planetary constants
132! ------------------------
133rad = 3397200.            ! mars radius (m) ~3397200 m
134daysec = 88775.           ! length of a sol (s) ~88775 s
135omeg = 4.*asin(1.)/daysec ! rotation rate (rad.s-1)
136g = 3.72                  ! gravity (m.s-2) ~3.72
137mugaz = 43.49             ! atmosphere molar mass (g.mol-1) ~43.49
138rcp = .256793             ! = r/cp ~0.256793
139r = 8.314511*1000./mugaz
140cpp = r/rcp
141year_day = 669            ! length of year (sols) ~668.6
142periheli = 206.66         ! minimum sun-mars distance (Mkm) ~206.66
143aphelie = 249.22          ! maximum sun-mars distance (Mkm) ~249.22
144halfaxe = 227.94          ! demi-grand axe de l'ellipse
145peri_day = 485.           ! perihelion date (sols since N. Spring)
146obliquit = 25.2           ! Obliquity (deg) ~25.2
147eccentric = 0.0934        ! Eccentricity (0.0934)
148Lsperi = 251.             ! Solar longitude of perihelion
149
150! Planetary Boundary Layer and Turbulence parameters
151! --------------------------------------------------
152z0_default = 1.e-2 ! surface roughness (m) ~0.01
153emin_turb = 1.e-6  ! minimal turbulent energy ~1.e-8
154lmixmin = 30       ! mixing length ~100
155
156! Cap properties and surface emissivities
157! ---------------------------------------
158emissiv = 0.95         ! Bare ground emissivity ~.95
159emisice(1) = 0.95      ! Northern cap emissivity
160emisice(2) = 0.95      ! Southern cap emisssivity
161albedice(1) = 0.5      ! Northern cap albedo
162albedice(2) = 0.5      ! Southern cap albedo
163iceradius(1) = 100.e-6 ! mean scat radius of CO2 snow (north)
164iceradius(2) = 100.e-6 ! mean scat radius of CO2 snow (south)
165dtemisice(1) = 2.      ! time scale for snow metamorphism (north)
166dtemisice(2) = 2.      ! time scale for snow metamorphism (south)
167
168! Water ice properties
169! ---------------------------------------
170inertieice = 2100.
171rho_H2O_ice = 920.
172
173! Mesh surface (not a very usefull quantity in 1D)
174! ------------------------------------------------
175cell_area(1) = 1.
176
177! Check if there is a 'traceur.def' file and process it
178! Load tracer names from file 'traceur.def'
179! -----------------------------------------------------
180open(90,file = 'traceur.def',status = 'old',form = 'formatted',iostat = ierr)
181if (ierr /= 0) then
182    write(*,*) 'Cannot find required file "traceur.def"'
183    write(*,*) ' If you want to run with tracers, I need it'
184    write(*,*) ' ... might as well stop here ...'
185    error stop
186else
187    write(*,*) "init_testphys1d: Reading file traceur.def"
188    ! read number of tracers:
189    read(90,*,iostat = ierr) nq
190    nqtot = nq ! set value of nqtot (in infotrac module) as nq
191    if (ierr /= 0) then
192        write(*,*) "init_testphys1d: error reading number of tracers"
193        write(*,*) "   (first line of traceur.def) "
194        error stop
195    endif
196    if (nq < 1) then
197        write(*,*) "init_testphys1d: error number of tracers"
198        write(*,*) "is nq=",nq," but must be >=1!"
199        error stop
200    endif
201endif
202
203! Allocate arrays
204allocate(tname(nq),q(1,nlayer,nq),zqsat(nlayer))
205allocate(dq(1,nlayer,nq),dqdyn(1,nlayer,nq),tnom_transp(nq))
206
207! Read tracer names from file traceur.def
208do iq = 1,nq
209    read(90,'(80a)',iostat = ierr) line ! store the line from traceur.def
210    if (ierr /= 0) then
211        error stop 'init_testphys1d: error reading tracer names...'
212    endif
213    ! if format is tnom_0, tnom_transp (isotopes)
214    read(line,*,iostat = ierr0) tname(iq),tnom_transp(iq)
215    if (ierr0 /= 0) then
216        read(line,*) tname(iq)
217        tnom_transp(iq) = 'air'
218    endif
219enddo
220close(90)
221
222! Isotopes: as in the 3D case we have to determine father/son relations for isotopes and carrying fluid
223allocate(nqfils(nqtot))
224nqperes = 0
225nqfils = 0
226do iq = 1,nqtot
227    if (tnom_transp(iq) == 'air') then
228        ! ceci est un traceur père
229        write(*,*) 'Le traceur',iq,', appele ',trim(tname(iq)),', est un pere'
230        nqperes = nqperes + 1
231    else !if (tnom_transp(iq) == 'air') then
232        ! ceci est un fils. Qui est son père?
233        write(*,*) 'Le traceur',iq,', appele ',trim(tname(iq)),', est un fils'
234        continu = .true.
235        ipere = 1
236        do while (continu)
237            if (tnom_transp(iq) == tname(ipere)) then
238                ! Son père est ipere
239                write(*,*) 'Le traceur',iq,'appele ',trim(tname(iq)),' est le fils de ',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 > nqtot) then
245                    write(*,*) 'Le traceur',iq,'appele ',trim(tname(iq)),', est orpelin.'
246                    call abort_gcm('infotrac_init','Un traceur est orphelin',1)
247                endif !if (ipere > nqtot) then
248            endif !if (tnom_transp(iq) == tnom_0(ipere)) then
249        enddo !do while (continu)
250    endif !if (tnom_transp(iq) == 'air') then
251enddo !do iq=1,nqtot
252write(*,*) 'nqperes=',nqperes
253write(*,*) 'nqfils=',nqfils
254
255#ifdef CPP_XIOS
256    call init_physics_distribution(regular_lonlat,4,1,1,1,nlayer,COMM_LMDZ)
257#else
258    call init_physics_distribution(regular_lonlat,4,1,1,1,nlayer,1)
259#endif
260
261! Date and local time at beginning of run
262! ---------------------------------------
263if (.not. therestartfi) then
264    ! Date (in sols since spring solstice) at beginning of run
265    day0 = 0 ! default value for day0
266    write(*,*) 'Initial date (in martian sols; =0 at Ls=0)?'
267    call getin("day0",day0)
268    day = float(day0)
269    write(*,*) " day0 = ",day0
270    ! Local time at beginning of run
271    time = 0 ! default value for time
272    write(*,*)'Initial local time (in hours, between 0 and 24)?'
273    call getin("time",time)
274    write(*,*)" time = ",time
275    time = time/24. ! convert time (hours) to fraction of sol
276else
277    call open_startphy(trim(startfiname))
278    call get_var("controle",tab_cntrl,found)
279    if (.not. found) then
280        call abort_physic("open_startphy","tabfi: Failed reading <controle> array",1)
281    else
282        write(*,*)'tabfi: tab_cntrl',tab_cntrl
283    endif
284    day0 = tab_cntrl(3)
285    day = float(day0)
286
287    call get_var("Time",time,found)
288    call close_startphy
289endif !(.not. therestartfi)
290
291! Discretization (Definition of grid and time steps)
292! --------------------------------------------------
293nlevel = nlayer + 1
294nsoil = nsoilmx
295
296steps_per_sol = 48 ! Default value for day_step (old name for steps_per_sol)
297write(*,*)'Number of time steps per sol?'
298call getin("day_step",steps_per_sol)
299write(*,*) " steps_per_sol (aka day_step) = ",steps_per_sol
300
301outputs_per_sol = steps_per_sol ! Default value for outputs_per_sol
302
303ndt = 10 ! Default value for ndt
304write(*,*)'Number of sols to run?'
305call getin("ndt",ndt)
306write(*,*) " ndt = ",ndt
307
308dayn = day0 + ndt
309ndt = ndt*steps_per_sol
310dttestphys = daysec/steps_per_sol
311
312! Imposed surface pressure
313! ------------------------
314psurf = 610. ! Default value for psurf
315pa = 20.     ! Default value for transition pressure (for hybrid coord.)
316preff = 610. ! Default value for reference surface pressure
317write(*,*) 'Surface pressure (Pa)?'
318if (.not. therestart1D) then
319    call getin("psurf",psurf)
320else
321    open(3,file = trim(start1Dname),status = "old",action = "read")
322    read(3,*) header, psurf, pa, preff
323endif
324write(*,*) " psurf = ", psurf
325write(*,*) " pa    = ", pa
326write(*,*) " preff = ", preff
327
328! Compute pressures and altitudes of atmospheric levels
329! -----------------------------------------------------
330! Vertical Coordinates
331! """"""""""""""""""""
332hybrid = .true.
333write(*,*)'Hybrid coordinates?'
334call getin("hybrid",hybrid)
335write(*,*) " hybrid = ", hybrid
336
337call disvert_noterre
338! Now that disvert has been called, initialize module vertical_layers_mod
339call init_vertical_layers(nlayer,preff,scaleheight,ap,bp,aps,bps,presnivs,pseudoalt)
340
341plev = ap + psurf*bp
342play = aps + psurf*bps
343zlay = -200.*r*log(play/plev(1))/g
344
345! Coefficient to control the surface pressure change
346! To be defined in "callphys.def" so that the PEM can read it asa well
347! If CO2cond_ps = 0, surface pressure is kept constant; If CO2cond_ps = 1, surface pressure varies normally
348! CO2cond_ps = 300*400/144.37e6 = 8.31e-4 (ratio of polar cap surface over planetary surface) is a typical value for tests
349! --------------------------------------------------------------------
350CO2cond_ps = 1. ! Default value
351write(*,*) 'Coefficient to control the surface pressure change?'
352call getin("CO2cond_ps",CO2cond_ps)
353if (CO2cond_ps < 0. .or. CO2cond_ps > 1.) then
354    write(*,*) 'Value for ''CO2cond_ps'' is not between 0 and 1.'
355    error stop 'Please, specify a correct value!'
356endif
357
358! Aerosol properties
359! ------------------
360tauvis = 0.2 ! Default value for tauvis (dust opacity)
361write(*,'("Reference dust opacity at ",f4.0," Pa?")') odpref
362call getin("tauvis",tauvis)
363write(*,*) " tauvis =",tauvis
364
365! Orbital parameters
366! ------------------
367if (.not. therestartfi) then
368    paleomars = .false. ! Default: no water ice reservoir
369    call getin("paleomars",paleomars)
370    write(*,*) "paleomars=", paleomars
371    if (paleomars) then
372        write(*,*) "Orbital parameters taken from callphys.def"
373
374        write(*,*)'Obliquity?'
375        call getin("obliquit",obliquit)
376        write(*,*) " obliquit =",obliquit
377
378        write(*,*) 'Eccentricity (0<e<1)?'
379        call getin('eccentric',eccentric)
380        write(*,*)"eccentric =",eccentric
381        if (eccentric < 0. .or. eccentric > 1.) then
382            write(*,*) 'Value for ''eccentric'' is not between 0 and 1.'
383            error stop 'Please, specify a correct value!'
384        endif
385
386        call getin('Lsperi',Lsperi)
387        write(*,*) 'Solar longitude of perihelion (0<Ls<360)?'
388        write(*,*)"Lsperi =",Lsperi
389        if (eccentric < 0. .or. eccentric > 360.) then
390            write(*,*) 'Value for ''Lsperi'' is not between 0 and 360.'
391            error stop 'Please, specify a correct value!'
392        endif
393       
394        Lsperi = Lsperi*pi/180. ! Convert to rad for peri_day
395        periheli = halfaxe*(1. - eccentric)
396        aphelie = halfaxe*(1. + eccentric)
397        call call_dayperi(Lsperi,eccentric,peri_day,year_day)
398        write(*,*) "Corresponding orbital params for PCM"
399        write(*,*) "   periheli =",periheli
400        write(*,*) "   aphelie  =",aphelie
401        write(*,*) "   peri_day =",peri_day
402    else
403        write(*,*) "Default present-day orbital parameters"
404     endif
405else
406    write(*,*) "Orbital parameters taken from startfi.nc"
407endif !(.not. therestartfi)
408
409! Latitude/longitude
410! ------------------
411latitude = 0. ! default value for latitude
412write(*,*) 'Latitude (in degrees)?'
413call getin("latitude",latitude(1))
414write(*,*) " latitude =",latitude
415latitude = latitude*pi/180.
416longitude = 0.
417longitude = longitude*pi/180.
418
419! Some initializations (some of which have already been done above!)
420! and loads parameters set in callphys.def and allocates some arrays
421! Mars possible matter with dttestphys in input and include!!!
422! Initializations below should mimick what is done in iniphysiq for 3D GCM
423call init_interface_dyn_phys
424call init_regular_lonlat(1,1,longitude,latitude,(/0.,0./),(/0.,0./))
425call init_geometry(1,longitude,latitude,(/0.,0.,0.,0./),(/0.,0.,0.,0./),cell_area,[1])
426call init_geometry_cell_area_for_outputs(1,cell_area)
427! Ehouarn: init_vertial_layers called later (because disvert not called yet)
428! call init_vertical_layers(nlayer,preff,scaleheight,ap,bp,aps,bps,presnivs,pseudoalt)
429call init_dimphy(1,nlayer) ! Initialize dimphy module
430call phys_state_var_init(1,llm,nq,tname,day0,dayn,time,daysec,dttestphys,rad,g,r,cpp,nqperes,nqfils) ! MVals: variables isotopes
431call ini_fillgeom(1,latitude,longitude,(/1.0/))
432call conf_phys(1,llm,nq)
433call initracer(ngrid,nq)
434
435! In 1D model physics are called every time step
436! ovverride iphysiq value that has been set by conf_phys
437if (iphysiq /= 1) then
438    write(*,*) "init_testphys1d: setting iphysiq=1"
439    iphysiq = 1
440endif
441
442! Initialize local slope parameters
443! ---------------------------------
444nslope_tmp = nslope
445write(*,*)'Number of slopes?' ! Number of subgrid scale slopes
446call getin("nslope",nslope_tmp)
447write(*,*) " nslope = ", nslope_tmp
448if (.not. therestartfi) then
449    nslope = nslope_tmp
450    call end_comslope_h
451    call ini_comslope_h(1,nslope)
452    select case (nslope)
453        case(1)
454            ! Sub-slopes parameters (minimum/maximun angles)
455            def_slope(1) = -90.
456            def_slope(2) = 90.
457            ! Fraction of subslopes in mesh
458            subslope_dist(1,1) = 1.
459        case(3)
460            ! Sub-slopes parameters (minimum/maximun angles)
461            def_slope(1) = -50.
462            def_slope(2) = -3.
463            def_slope(3) = 3.
464            def_slope(4) = 50.
465            ! Fraction of subslopes in mesh
466            subslope_dist(1,1) = 0.1
467            subslope_dist(1,2) = 0.8
468            subslope_dist(1,3) = 0.1
469        case(5)
470            ! Sub-slopes parameters (minimum/maximun angles)
471            def_slope(1) = -43.
472            def_slope(2) = -9.
473            def_slope(3) = -3.
474            def_slope(4) = 3.
475            def_slope(5) = 9.
476            def_slope(6) = 43.
477            ! Fraction of subslopes in mesh
478            subslope_dist(1,1) = 0.025
479            subslope_dist(1,2) = 0.075
480            subslope_dist(1,3) = 0.8
481            subslope_dist(1,4) = 0.075
482            subslope_dist(1,5) = 0.025
483        case(7)
484            ! Sub-slopes parameters (minimum/maximun angles)
485            def_slope(1) = -43.
486            def_slope(2) = -19.
487            def_slope(3) = -9.
488            def_slope(4) = -3.
489            def_slope(5) = 3.
490            def_slope(6) = 9.
491            def_slope(7) = 19.
492            def_slope(8) = 43.
493            ! Fraction of subslopes in mesh
494            subslope_dist(1,1) = 0.01
495            subslope_dist(1,2) = 0.02
496            subslope_dist(1,3) = 0.06
497            subslope_dist(1,4) = 0.8
498            subslope_dist(1,5) = 0.06
499            subslope_dist(1,6) = 0.02
500            subslope_dist(1,7) = 0.01
501        case default
502            error stop 'Number of slopes not possible: nslope should be 1, 3, 5 or 7!'
503    end select
504
505    call end_surfdat_h_slope_var
506    call ini_surfdat_h_slope_var(1,nq,nslope)
507    call end_comsoil_h_slope_var
508    call ini_comsoil_h_slope_var(1,nslope)
509    call end_dimradmars_mod_slope_var
510    call ini_dimradmars_mod_slope_var(1,nslope)
511    call end_paleoclimate_h
512    call ini_paleoclimate_h(1,nslope)
513else
514    if (nslope /= nslope_tmp) then
515        error stop 'Number of slopes not compatible between "startfi.nc" and "run.def".'
516    else
517        write(*,*) 'Be careful: subslopes distribution and definition are taken from "startfi.nc" and not initialized by default!'
518    endif
519endif
520
521! Slope inclination angle (deg) 0: horizontal, 90: vertical
522theta_sl = 0. ! default: no inclination
523call getin("slope_inclination",theta_sl)
524
525! Slope orientation (deg)
526! 0 == Northward, 90 == Eastward, 180 == Southward, 270 == Westward
527psi_sl = 0. ! default value
528call getin("slope_orientation",psi_sl)
529
530! Initialize tracers (surface and atmosphere)
531!--------------------------------------------
532write(*,*) "init_testphys1d: initializing tracers"
533if (.not. therestart1D) then
534    call read_profile(nq,nlayer,qsurf(1,:,1),q,play)
535    do iq = 1,nq
536        qsurf(1,iq,:) = qsurf(1,iq,1)
537    enddo
538else
539    do iq = 1,nq
540        read(3,*,iostat = ierr) header, (q(1,ilayer,iq), ilayer = 1,nlayer)
541        if (ierr /= 0) error stop 'Not enough atmospheric layers defined in the file "'//trim(start1Dname)//'" for the tracer "'//trim(header)//'"!'
542        if (trim(tname(iq)) /= trim(header)) then
543            write(*,*) 'Tracer names between "traceur.def" and "'//trim(start1Dname)//'" do not match!'
544            write(*,*) 'Please, write the tracer names in the same order for both files.'
545            error stop
546        endif
547    enddo
548endif
549
550! Initialize albedo / soil thermal inertia
551! ----------------------------------------
552if (.not. therestartfi) then
553    albedodat(1) = 0.2 ! default value for albedodat
554    write(*,*)'Albedo of bare ground?'
555    call getin("albedo",albedodat(1))
556    write(*,*) " albedo = ",albedodat(1)
557    albedo(1,:,:) = albedodat(1)
558
559    inertiedat(1,1) = 400 ! default value for inertiedat
560    write(*,*)'Soil thermal inertia (SI)?'
561    call getin("inertia",inertiedat(1,1))
562    write(*,*) " inertia = ",inertiedat(1,1)
563
564    ice_depth = -1 ! default value: no ice
565    call getin("subsurface_ice_depth",ice_depth)
566    write(*,*) " subsurface_ice_depth = ",ice_depth
567    h2o_ice_depth(:,:) = ice_depth
568    z0(1) = z0_default ! default value for roughness
569    write(*,*) 'Surface roughness length z0 (m)?'
570    call getin("z0",z0(1))
571    write(*,*) " z0 = ",z0(1)
572endif !(.not. therestartfi)
573
574! For the gravity wave scheme
575! ---------------------------
576zmea(1) = 0.
577zstd(1) = 0.
578zsig(1) = 0.
579zgam(1) = 0.
580zthe(1) = 0.
581
582! For the non-orographic gravity waves scheme
583du_nonoro_gwd(1,:) = 0
584dv_nonoro_gwd(1,:) = 0
585
586! For the slope wind scheme
587! -------------------------
588hmons(1) = 0.
589write(*,*)'hmons is initialized to ',hmons(1)
590summit(1) = 0.
591write(*,*)'summit is initialized to ',summit(1)
592base(1) = 0.
593
594! Default values initializing the coefficients calculated later
595! -------------------------------------------------------------
596tauscaling(1) = 1.   ! calculated in aeropacity_mod.F
597totcloudfrac(1) = 1. ! calculated in watercloud_mod.F
598
599! Specific initializations for "physiq"
600! -------------------------------------
601! surface geopotential is not used (or useful) since in 1D
602! everything is controled by surface pressure
603phisfi(1) = 0.
604
605! Initialization to take into account prescribed winds
606! ----------------------------------------------------
607ptif = 2.*omeg*sinlat(1)
608
609! Geostrophic wind
610gru = 10. ! default value for gru
611write(*,*)'zonal eastward component of the geostrophic wind (m/s)?'
612call getin("u",gru)
613write(*,*) " u = ",gru
614grv = 0. !default value for grv
615write(*,*)'meridional northward component of the geostrophic wind (m/s)?'
616call getin("v",grv)
617write(*,*) " v = ",grv
618
619! Initialize winds for first time step
620!-------------------------------------
621if (.not. therestart1D) then
622    u = gru
623    v = grv
624else
625    read(3,*) header, (u(ilayer), ilayer = 1,nlayer)
626    read(3,*) header, (v(ilayer), ilayer = 1,nlayer)
627endif
628w = 0. ! default: no vertical wind
629
630! Initialize turbulent kinetic energy
631!------------------------------------
632q2 = 0.
633
634! CO2 ice on the surface
635! ----------------------
636if (.not. therestartfi) then
637    qsurf(1,igcm_co2,:) = 0. ! default value for co2ice
638    write(*,*)'Initial CO2 ice on the surface (kg.m-2)'
639    call getin("co2ice",qsurf(1,igcm_co2,1))
640    qsurf(1,igcm_co2,:) = qsurf(1,igcm_co2,1)
641    write(*,*) " co2ice = ",qsurf(1,igcm_co2,1)
642endif !(.not. therestartfi)
643
644! Emissivity
645! ----------
646if (.not. therestartfi) then
647    emis = emissiv
648    if (qsurf(1,igcm_co2,1) == 1.) then
649        emis = emisice(1) ! northern hemisphere
650        if (latitude(1) < 0) emis = emisice(2) ! southern hemisphere
651    endif
652endif !(.not. therestartfi)
653
654! Initialize temperature profile
655! ------------------------------
656pks = psurf**rcp
657
658if (.not. therestart1D) then
659    ! Altitude in km in "profile_temp": divide zlay by 1000
660    zkm(0) = 0.
661    zkm(1:) = zlay/1000.
662
663    call profile_temp(nlayer + 1,zkm,teta)
664    tsurf = teta(0)
665    temp = teta(1:)
666else
667    read(3,*) header, (tsurf(1,j), j = 1,size(tsurf,2)), (temp(ilayer), ilayer = 1,nlayer)
668    close(3)
669endif
670
671! Initialize soil properties, temperature and content
672! ------------------------------------------
673volcapa = 1.e6 ! volumetric heat capacity
674
675if (.not. therestartfi) qsoil = 0.
676
677if (.not. therestartfi) then
678    ! Initialize depths
679    ! -----------------
680    do isoil = 0,nsoil - 1
681       mlayer(isoil) = 2.e-4*(1 + isoil**2.9*(1 - exp(-real(isoil)/20.))) ! mid layer depth
682    enddo
683    do isoil = 1,nsoil - 1
684       layer(isoil) = (mlayer(isoil) + mlayer(isoil - 1))/2
685    enddo
686    layer(nsoil) = 2*mlayer(nsoil - 1) - mlayer(nsoil - 2)
687
688    ! Creating the new soil inertia table if there is subsurface ice:
689    if (ice_depth > 0) then
690        iref = 1 ! ice/regolith boundary index
691        if (ice_depth < layer(1)) then
692            inertiedat(1,1) = sqrt(layer(1)/((ice_depth/inertiedat(1,1)**2) + ((layer(1) - ice_depth)/inertieice**2)))
693            inertiedat(1,2:) = inertieice
694            if(adsorption_soil) then
695            ! add subsurface ice in qsoil if one runs with adsorption
696               if(ads_massive_ice) then
697                   qsoil(1,1,igcm_h2o_ice_soil,:) = (ice_depth-layer(2))/(layer(1) - layer(2))*rho_H2O_ice ! assume no porosity in the ice
698                   qsoil(1,2:,igcm_h2o_ice_soil,:) = rho_H2O_ice
699               else
700                   qsoil(1,1,igcm_h2o_ice_soil,:) = (ice_depth-layer(2))/(layer(1) - layer(2))*porosity_reg*rho_H2O_ice
701                   qsoil(1,2:,igcm_h2o_ice_soil,:) = porosity_reg*rho_H2O_ice
702               endif           
703            endif
704        else ! searching for the ice/regolith boundary:
705            do isoil = 1,nsoil
706                if ((ice_depth >= layer(isoil)) .and. (ice_depth < layer(isoil + 1))) then
707                    iref = isoil + 1
708                    exit
709                endif
710            enddo
711            ! We then change the soil inertia table:
712            inertiedat(1,:iref - 1) = inertiedat(1,1)
713            ! We compute the transition in layer(iref)
714            inertiedat(1,iref) =  sqrt((layer(iref) - layer(iref - 1))/(((ice_depth - layer(iref - 1))/inertiedat(1,1)**2) + ((layer(iref) - ice_depth)/inertieice**2)))
715            ! Finally, we compute the underlying ice:
716            inertiedat(1,iref + 1:) = inertieice
717           
718            if(adsorption_soil) then
719            ! add subsurface ice in qsoil if one runs with adsorption
720                qsoil(1,:iref - 1,igcm_h2o_ice_soil,:) = 0.
721                if(ads_massive_ice) then
722                   qsoil(1,iref+1:,igcm_h2o_ice_soil,:) = rho_H2O_ice
723                   qsoil(1,iref,igcm_h2o_ice_soil,:) = (ice_depth-layer(iref))/(layer(iref-1) - layer(iref))*rho_H2O_ice
724                else
725                   qsoil(1,iref+1:,igcm_h2o_ice_soil,:) = porosity_reg*rho_H2O_ice
726                   qsoil(1,iref,igcm_h2o_ice_soil,:) = (ice_depth-layer(iref))/(layer(iref-1) - layer(iref))*porosity_reg*rho_H2O_ice
727                endif
728            endif
729        endif ! (ice_depth < layer(1))
730    else ! ice_depth < 0 all is set to surface thermal inertia
731        inertiedat(1,:) = inertiedat(1,1) ! soil thermal inertia
732    endif ! ice_depth > 0
733   
734    do isoil = 1,nsoil
735        inertiesoil(1,isoil,:) = inertiedat(1,isoil)
736        tsoil(1,isoil,:) = tsurf(1,:) ! soil temperature
737    enddo
738endif !(.not. therestartfi)
739
740! Initialize subsurface geothermal flux
741! -------------------------------------
742flux_geo = 0.
743call getin("flux_geo",flux_geo(1,1))
744write(*,*) " flux_geo = ",flux_geo(1,1)
745flux_geo = flux_geo(1,1)
746
747! Initialize traceurs
748! -------------------
749if (photochem .or. callthermos) then
750    write(*,*) 'Initializing chemical species'
751    ! flagthermo=0: initialize over all atmospheric layers
752    flagthermo = 0
753    ! check if "h2o_vap" has already been initialized
754    ! (it has been if there is a "profile_h2o_vap" file around)
755    inquire(file = "profile_h2o_vap",exist = found)
756    if (found) then
757        flagh2o = 0 ! 0: do not initialize h2o_vap
758    else
759        flagh2o = 1 ! 1: initialize h2o_vap in inichim_newstart
760    endif
761
762    ! hack to accomodate that inichim_newstart assumes that
763    ! q & psurf arrays are on the dynamics scalar grid
764    allocate(qdyn(2,1,llm,nq),psdyn(2,1))
765    qdyn(1,1,1:llm,1:nq) = q(1,1:llm,1:nq)
766    psdyn(1:2,1) = psurf
767    call inichim_newstart(ngrid,nq,qdyn,qsurf(1,:,1),psdyn,flagh2o,flagthermo)
768    q(1,1:llm,1:nq) = qdyn(1,1,1:llm,1:nq)
769endif
770
771! Check if the surface is a water ice reservoir
772! ---------------------------------------------
773if (.not. therestartfi) watercap(1,:) = 0 ! Initialize watercap
774watercaptag(1) = .false. ! Default: no water ice reservoir
775write(*,*)'Water ice cap on ground?'
776call getin("watercaptag",watercaptag)
777write(*,*) " watercaptag = ",watercaptag
778
779! Initialize perennial_co2ice
780! ---------------------------
781if (.not. therestartfi .and. paleoclimate) perennial_co2ice = 0.
782
783! Prescription of atmospheric water vapour/ice profiles
784! ----------------------------------------------
785ctrl_h2ovap = .false. ; ctrl_h2oice = .false. ! Default: free atmospheric water vapour/ice profiles
786relaxtime_h2ovap = -1. ; relaxtime_h2oice = -1. ! Default: no time relaxation
787qref_h2ovap = 0. ; qref_h2oice = 0. ! Default
788if (water) then
789    write(*,*)'Force atmospheric water vapour profile (mixing ratio in kg/kg)?'
790    call getin('ctrl_h2ovap',ctrl_h2ovap)
791    write(*,*) 'ctrl_h2ovap = ', ctrl_h2ovap
792    if (ctrl_h2ovap) then
793        write(*,*) 'Atmospheric water vapour profile is prescribed'
794        open(10,file = 'profile_ref_h2o_vap',status = 'old',form = 'formatted',action = 'read',iostat = ierr)
795        if (ierr == 0) then
796            do ilayer = 1,nlayer
797                read(10,*,iostat=ierr) qref_h2ovap(ilayer)
798                if (ierr /= 0) error stop 'Not enough atmospheric layers defined in the file "profile_ref_h2o_vap"!'
799            enddo
800        else
801            error stop 'File "profile_ref_h2o_vap" was not found!'
802        endif
803        close(10)
804
805        write(*,*) 'Relax atmospheric water vapour profile with time constant (s)?'
806       
807        call getin('relaxtime_h2ovap',relaxtime_h2ovap)
808        write(*,*) 'relaxtime_h2ovap = ', relaxtime_h2ovap
809        if (relaxtime_h2ovap < 0.) then
810            write(*,*) 'Atmospheric water vapour profile is not relaxed (forcing).'
811        else
812            write(*,*) 'Atmospheric water vapour profile is relaxed towards the profile in "profile_ref_h2o_vap"'
813        endif
814    else
815        write(*,*) 'Free atmospheric water vapour profile'
816    endif
817
818    write(*,*)'Force atmospheric water ice profile (mixing ratio in kg/kg)?'
819    call getin('ctrl_h2oice',ctrl_h2oice)
820    write(*,*) 'ctrl_h2oice = ', ctrl_h2oice
821    if (ctrl_h2oice) then
822        write(*,*) 'Atmospheric water ice profile is prescribed'
823        open(10,file = 'profile_ref_h2o_ice',status = 'old',form = 'formatted',action = 'read',iostat = ierr)
824        if (ierr == 0) then
825            do ilayer = 1,nlayer
826                read(10,*,iostat=ierr) qref_h2oice(ilayer)
827                if (ierr /= 0) error stop 'Not enough atmospheric layers defined in the file "profile_ref_h2o_ice"!'
828            enddo
829        else
830            error stop 'File "profile_ref_h2o_ice" was not found!'
831        endif
832        close(10)
833
834        write(*,*) 'Relax atmospheric water ice profile with time constant (s)?'
835       
836        call getin('relaxtime_h2oice',relaxtime_h2oice)
837        write(*,*) 'relaxtime_h2oice = ', relaxtime_h2oice
838        if (relaxtime_h2oice < 0.) then
839            write(*,*) 'Atmospheric water ice profile is not relaxed (forcing).'
840        else
841            write(*,*) 'Atmospheric water ice profile is relaxed towards the profile in "profile_ref_h2o_ice"'
842        endif
843    else
844        write(*,*) 'Free atmospheric water ice profile'
845    endif
846endif
847
848END SUBROUTINE init_testphys1d
849
850END MODULE init_testphys1d_mod
Note: See TracBrowser for help on using the repository browser.