62 REAL,
PARAMETER :: tyear = 4.05 *
day
63 REAL,
PARAMETER :: rplanet = 0.788*
rearth
64 REAL,
PARAMETER :: theta0 = 0.0
65 REAL,
PARAMETER :: phi0 = 0.0
66 REAL,
PARAMETER ::
omega = 2*
pi/tyear/8
67 REAL,
PARAMETER :: fsun =
omega/(2*
pi) - &
69 REAL,
PARAMETER :: mass = 0.297*
mearth
70 REAL,
PARAMETER :: gacc =
gn*mass/rplanet**2
72 REAL,
PARAMETER :: sm_axis = 0.02219*
au
73 REAL,
PARAMETER :: excent = 0.0
75 REAL,
PARAMETER :: rstar = 0.121*
rsun
76 REAL,
PARAMETER :: tstar = 2511
77 REAL,
PARAMETER :: fstar = (rstar/
au)**2 &
80 REAL,
PARAMETER :: tsim = 10 * tyear
81 INTEGER,
PARAMETER :: xres = 1
82 INTEGER,
PARAMETER :: yres = 16
83 INTEGER,
PARAMETER :: zres = 32
85 REAL,
PARAMETER :: gamma = 1.4
86 REAL,
PARAMETER :: p0 = 1.0e7
87 REAL,
PARAMETER :: mu = 2.8586e-2
88 REAL,
PARAMETER :: t0 = 258.1
90 REAL,
PARAMETER :: albedo = 0.306
92 INTEGER,
PARAMETER :: onum = 10
93 CHARACTER(LEN=256),
PARAMETER &
95 CHARACTER(LEN=256),
PARAMETER &
96 :: ofname =
'planet2d'
98 CLASS(
fosite),
ALLOCATABLE :: sim
108 CALL sim%InitFosite()
112 IF (
ALLOCATED(sim%Sources))
THEN
113 sp => sim%Sources%GetSourcesPointer(planet_cooling)
114 IF (
ASSOCIATED(sp))
THEN
123 IF (.NOT.
ASSOCIATED(spcool)) &
124 CALL sim%Error(
"planet2d",
"planetary cooling sources not initialized")
126 CALL initdata(sim%Mesh, sim%Physics, sim%Timedisc%pvar, sim%Timedisc%cvar)
131 tmean = sum(spcool%T_s%data1d(:)*sim%Mesh%volume%data1d(:),mask=sim%Mesh%without_ghost_zones%mask1d(:)) &
132 / sum(sim%Mesh%volume%data1d(:),mask=sim%Mesh%without_ghost_zones%mask1d(:))
136 tap_check_close(t0,tmean,t0*0.01,
"deviation of final mean surface temperature from initial T0 < 1%")
148 TYPE(
dict_typ),
POINTER :: mesh, physics, boundary, datafile, &
149 timedisc, fluxes, heating, sources, &
150 grav, pmass, cooling, rotframe
154 "meshtype" / midpoint, &
155 "geometry" / spherical_planet, &
156 "decomposition" / (/ 1, -1, 1/), &
162 "xmax" / (
pi-0.05), &
167 "gparam" / rplanet, &
168 "output/rotation" / 0, &
169 "output/volume" / 1, &
174 "western" / reflecting, &
175 "eastern" / reflecting, &
176 "southern" / periodic, &
177 "northern" / periodic, &
178 "bottomer" / reflecting, &
179 "topper" / reflecting)
192 "variables" / primitive, &
203 "stype" / rotating_frame, &
204 "disable_centaccel" / 1)
208 "stype" / planet_cooling, &
209 "output/Qcool" / 1, &
211 "output/RHO_s" / 1, &
215 "mean_distance" / (sm_axis*(1. + 0.5*excent**2)), &
222 "stype" / planet_heating, &
223 "output/Qstar" / 1, &
230 "semimajoraxis" / sm_axis, &
231 "excentricity" / excent, &
235 "rotframe" / rotframe, &
236 "cooling" / cooling, &
244 "tol_rel" / 0.0095, &
245 "dtlimit" / 1.0e-15, &
246 "maxiter" / 1000000000)
249 "fileformat" / vtk, &
250 "filename" / (trim(odir) // trim(ofname)), &
255 "physics" / physics, &
256 "boundary" / boundary, &
258 "timedisc" / timedisc, &
259 "sources" / sources, &
260 "datafile" / datafile)
267 CLASS(physics_base),
INTENT(IN) :: Physics
268 CLASS(mesh_base),
INTENT(IN) :: Mesh
269 CLASS(marray_compound),
POINTER,
INTENT(INOUT) :: pvar,cvar
273 SELECT TYPE(p => pvar)
274 CLASS IS(statevector_euler)
278 p%density%data1d(:) = p0/gacc
283 p%pressure%data1d(:) = gamma*physics%Constants%RG / mu * t0 * p%density%data1d(:)
286 p%velocity%data1d(:) = 0.0
288 CALL mesh%Error(
"planet2d::InitData",
"state vector not supported")
291 CALL physics%Convert2Conservative(pvar,cvar)
292 CALL mesh%Info(
" DATA-----> initial condition: 2D planetary atmosphere")
subroutine initdata(Mesh, Physics, Fluxes, Timedisc)
subroutine makeconfig(Sim, config)
Dictionary for generic data types.
type(dict_typ) function, pointer, public dict(n1, n2, n3, n4, n5, n6, n7, n8, n9, n10, n11, n12, n13, n14, n15, n16, n17, n18, n19, n20)
Construct a new dictionary from several key/value pairs. Together with the Assign subroutine and over...
module for SI units and physical constants
real, parameter, public mjupiter
jupiter mass [kg]
real, parameter, public mearth
earth mass [kg]
real, parameter, public gn
gravitational constant [m^3/kg/s^2]
real, parameter, public rsun
sun radius [m]
real, parameter, public ke
electron scattering opacity [m^2/kg]
real, parameter, public c
vacuum speed of light [m/s]
real, parameter, public day
length of a day [s]
real, parameter, public sb
Stefan-Boltzmann constant [W/m^2/K^4].
real, parameter, public au
astronomical unit [m]
real, parameter, public na
Avogadro constant [1/mol].
real, parameter, public kb
Boltzmann constant [J/K].
real, parameter, public pi
real, parameter, public rjupiter
jupiter radius [m]
real, parameter, public rearth
earth radius [m]
real, parameter, public msun
sun mass [kg]
generic source terms module providing functionaly common to all source terms
gray cooling of planetary atmospheres
pure real function, dimension(size(ephir%data2d, dim=1)) omega(ephir, velocity)