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

Last change on this file since 3207 was 3207, checked in by jbclement, 12 months ago

Mars PCM:
Modification of r3200 for the 1D model: the 'CO2cond_ps' coefficient (now defined in "co2condens_mod.F") is applied not only to update the surface pressure but also all the tracers and the winds from now on. It is done in "physiq_mod.F" just after the call to 'co2condens'.
JBC

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