67 REAL,
PARAMETER :: tsim = 5.0e-3
71 REAL,
PARAMETER :: re = 1.0e+4
72 REAL,
PARAMETER :: ma = 1.0e+2
73 REAL,
PARAMETER :: mdisk = 1.0e-2
75 REAL,
PARAMETER :: rhomin = 1.0e-8
77 INTEGER,
PARAMETER :: vistype = beta
79 REAL,
PARAMETER :: pl_exp = 0.0
80 REAL,
PARAMETER :: tvis = 4./3.*re
81 REAL,
PARAMETER :: tinit = 1.0e-3
85 INTEGER,
PARAMETER :: mgeo = spherical
86 INTEGER,
PARAMETER :: xres = 200
87 INTEGER,
PARAMETER :: yres = 16
88 INTEGER,
PARAMETER :: zres = 6
89 REAL,
PARAMETER :: rmin = 0.1
90 REAL,
PARAMETER :: rmax = 2.0
91 REAL,
PARAMETER :: gpar = 1.0
93 INTEGER,
PARAMETER :: onum = 100
94 CHARACTER(LEN=256),
PARAMETER &
96 CHARACTER(LEN=256),
PARAMETER &
97 :: ofname =
'pringle3d'
99 REAL,
PARAMETER :: csiso = 1./ma
101 CLASS(
fosite),
ALLOCATABLE :: sim
105 CALL sim%InitFosite()
108 CALL initdata(sim,sim%Mesh, sim%Physics, sim%Timedisc)
124 TYPE(
dict_typ),
POINTER :: mesh, physics, boundary, datafile, &
125 grav, vis, timedisc, fluxes, sources
126 REAL :: x1,x2,y1,y2,z1,z2,cvis
142 bc(west) = no_gradients
143 bc(east) = no_gradients
152 bc(bottom) = reflecting
156 cvis = (xres/1207.0)**1.9
164 bc(west) = no_gradients
165 bc(east) = no_gradients
168 bc(bottom) = no_gradients
169 bc(top) = no_gradients
186 bc(south) = absorbing
187 bc(north) = absorbing
188 bc(bottom) = periodic
192 cvis = (xres/1207.0)**1.9
194 CALL sim%Error(
"pringle::MakeConfig",
"mesh geometry not supported for Pringle disk")
199 "meshtype" / midpoint, &
216 "western" / bc(west), &
217 "eastern" / bc(east), &
218 "southern" / bc(south), &
219 "northern" / bc(north), &
220 "bottomer" / bc(bottom), &
226 "problem" / euler_isotherm, &
227 "units" / geometrical, &
234 "variables" / primitive, &
235 "limiter" / vanleer, &
240 "stype" / viscosity, &
241 "vismodel" / vistype, &
242 "dynconst" / (1./re), &
243 "exponent" / pl_exp, &
249 "pmass/gtype"/ pointmass, &
250 "pmass/mass" / 1.0, &
251 "pmass/outbound" / 0)
261 "stoptime" / (tsim * tvis), &
263 "dtlimit" / (1e-10 * tvis), &
265 "maxiter" / 100000000, &
267 "tol_abs" / (/1e-40,1e-5,1e-8,1e-5/))
278 "filename" / (trim(odir) // trim(ofname)), &
283 "physics" / physics, &
284 "boundary" / boundary, &
286 "sources" / sources, &
287 "timedisc" / timedisc, &
288 "datafile" / datafile)
297 CLASS(
fosite),
INTENT(INOUT) :: Sim
298 CLASS(physics_base),
INTENT(IN) :: Physics
299 CLASS(mesh_base),
INTENT(IN) :: Mesh
300 CLASS(timedisc_base),
INTENT(INOUT) :: Timedisc
305 CLASS(geometry_base),
ALLOCATABLE :: geo_cyl
306 TYPE(
dict_typ),
POINTER :: geo_config
311 REAL,
DIMENSION(Mesh%IGMIN:Mesh%IGMAX,Mesh%JGMIN:Mesh%JGMAX,Mesh%KGMIN:Mesh%KGMAX,3) &
313 REAL :: r,z,vr,vphi,sden,mu,height,mass
314 CHARACTER(LEN=64) :: value
316 geo_config =>
dict(
"geometry" / cylindrical)
321 CALL geo_cyl%Convert2Curvilinear(mesh%bccart,bccyl)
324 SELECT TYPE(bwest => timedisc%Boundary%boundary(west)%p)
325 CLASS IS (boundary_custom)
326 CALL bwest%SetCustomBoundaries(mesh,physics, &
327 (/custom_nograd,custom_outflow,custom_kepler,custom_nograd/))
329 SELECT TYPE(beast => timedisc%Boundary%boundary(east)%p)
330 CLASS IS (boundary_custom)
331 CALL beast%SetCustomBoundaries(mesh,physics, &
332 (/custom_nograd,custom_outflow,custom_kepler,custom_nograd/))
336 sp => sim%Sources%GetSourcesPointer(gravity)
337 IF (.NOT.
ASSOCIATED(sp))
CALL physics%Error(
"pringle::InitData",
"no gravity term initialized")
345 CALL timedisc%Error(
"pringle::InitData",
"only pringle and beta viscosity possible")
349 SELECT TYPE(pvar => timedisc%pvar)
350 TYPE IS(statevector_eulerisotherm)
351 DO k=mesh%KGMIN,mesh%KGMAX
352 DO j=mesh%JGMIN,mesh%JGMAX
353 DO i=mesh%IMIN,mesh%IGMAX
359 vphi = r/(r**2+z**2)**0.75
361 height = csiso * sqrt(r**3)
370 pvar%density%data3d(i,j,k) = rhomin + sden / (sqrt(2*pi) * height) &
371 * exp(-0.5*(z/height)**2)
379 mass = sum(mesh%volume%data1d(:)*pvar%density%data1d(:),mask=mesh%without_ghost_zones%mask1d)
381 CALL mpi_allreduce(mpi_in_place,mass,1,default_mpi_real,mpi_sum, &
382 mesh%comm_cart,ierror)
385 pvar%density%data1d(:) = pvar%density%data1d(:) * mdisk / mass
386 pvar%velocity%data1d(:) = 0.0
387 CALL gp%UpdateGravity(mesh,sim%Physics,sim%Fluxes,pvar,0.0,0.0)
388 pvar%velocity%data4d(:,:,:,1:physics%VDIM) = &
389 timedisc%GetCentrifugalVelocity(mesh,sim%Physics,sim%Fluxes,sim%Sources,(/0.,0.,1./), &
392 CALL timedisc%Error(
"pringle::InitData",
"only isothermal physics possible")
395 IF (
ASSOCIATED(timedisc%solution))
THEN
397 timedisc%solution = timedisc%pvar
401 CALL physics%Convert2Conservative(timedisc%pvar,timedisc%cvar)
403 CALL mesh%Info(
" DATA-----> initial condition: " //
"pringle disk")
404 WRITE(
value,
"(ES11.3)") tvis
405 CALL mesh%Info(
" " //
"viscous timescale: " //trim(
value))
406 WRITE(
value,
"(ES11.3)") re
407 CALL mesh%Info(
" " //
"Reynolds number: " //trim(
value))
408 WRITE(
value,
"(ES11.3)") ma
409 CALL mesh%Info(
" " //
"Mach number: " //trim(
value))
410 WRITE(
value,
"(ES11.3)") mdisk
411 CALL mesh%Info(
" " //
"Disk mass: " //trim(
value))
413 CALL geo_cyl%Finalize()
421 REAL,
INTENT(IN) :: mu,tau,r
422 REAL,
INTENT(OUT) :: Sigma,vr
424 REAL :: lambda,x,r4l,invr4l
427 r4l = r**(0.25/lambda)
428 x = 2*lambda**2 / tau * r**(0.25/lambda)
430 sigma = mdisk/sqrt((4*pi)**3 * tau) * r4l**(1.5-9*lambda) &
431 * exp(-lambda**2/tau*(1.0-r4l)**2)
432 vr = -1.5/re * ((x*tau)/(2*lambda*lambda))**(4*lambda-2.) &
433 * (1.0 - 0.5*x/lambda * (x*tau/(2*lambda*lambda) - 1.0))
440 REAL,
INTENT(IN) :: mu,tau,r
441 REAL,
INTENT(OUT) :: Sigma,vr
446 r4l = r**(0.25/lambda)
447 x = 2*lambda**2 / tau * r4l
448 sigma = mdisk/(4*pi) * (lambda/tau) * r4l**2 * r**(-2.25) &
449 * exp(-lambda**2/tau*(1.0+r4l**2)) &
451 vr = -1.5/re * ((x*tau)/(2*lambda*lambda))**(4*lambda-2.) &
452 * (1.0 - 0.5*x/lambda * (x*tau/(2*lambda*lambda) &
460 REAL,
INTENT(IN) :: nu,x
465 IF (aint(nu).EQ.nu)
THEN
468 inu = fgsl_sf_bessel_icn(n,x)
472 inu = (0.5*x)**nu / fgsl_sf_gamma(1.0+nu) * fgsl_sf_hyperg_0f1(1.0+nu,(0.5*x)**2)
subroutine initdata(Mesh, Physics, Fluxes, Timedisc)
subroutine makeconfig(Sim, config)
Dictionary for generic data types.
recursive subroutine, public deletedict(root)
Delete the dictionary 'root' and all subnodes.
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...
elemental real function, public asinh(x)
inverse hyperbolic sine function
constructor for geometry class
subroutine new_geometry(Geometry, config)
generic source terms module providing functionaly common to all source terms
generic gravity terms module providing functionaly common to all gravity terms
subroutine viscousring_analytic(mu, tau, r, Sigma, vr)
real function modified_bessel_inu(nu, x)
subroutine viscousring_approx(mu, tau, r, Sigma, vr)