55 REAL,
PARAMETER :: tsim = 2.0e-2
59 REAL,
PARAMETER :: re = 1.0e+4
60 REAL,
PARAMETER :: ma = 1.0e+2
61 REAL,
PARAMETER :: mdisk = 1.0e-2
63 REAL,
PARAMETER :: rhomin = 1.0e-20
66 INTEGER,
PARAMETER :: vistype = powerlaw
67 REAL,
PARAMETER :: pl_exp = 0.0
68 REAL,
PARAMETER :: tvis = 4./3.*re
69 REAL,
PARAMETER :: tinit = 1.0e-3
80 INTEGER,
PARAMETER :: mgeo = cylindrical
83 INTEGER,
PARAMETER :: xres = 400
84 INTEGER,
PARAMETER :: yres = 1
85 INTEGER,
PARAMETER :: zres = 1
86 REAL,
PARAMETER :: rmin = 0.05
87 REAL,
PARAMETER :: rmax = 2.0
88 REAL,
PARAMETER :: gpar = 1.0
90 INTEGER,
PARAMETER :: onum = 10
91 CHARACTER(LEN=256),
PARAMETER &
93 CHARACTER(LEN=256),
PARAMETER &
96 REAL,
PARAMETER :: csiso = 1./ma
101 CLASS(
fosite),
ALLOCATABLE :: sim
102 REAL,
DIMENSION(:),
ALLOCATABLE :: sigma
103 REAL :: sum_numer, sum_denom
104 INTEGER :: n,den,vx,vy
105 REAL :: next_output_time
111 CALL sim%InitFosite()
114 IF (sim%GetRank().EQ.0)
THEN
127 CALL initdata(sim,sim%Mesh, sim%Physics, sim%Timedisc, sim%Sources)
128 next_output_time = sim%Datafile%time
129 IF (
ASSOCIATED(sim%Timedisc%solution))
THEN
131 DO WHILE((sim%Timedisc%maxiter.LE.0).OR.(sim%iter.LE.sim%Timedisc%maxiter))
134 IF (next_output_time.LT.sim%Datafile%time)
THEN
136 next_output_time = sim%Datafile%time
137 SELECT TYPE(pvar => sim%Timedisc%solution)
138 TYPE IS (statevector_eulerisotherm)
140 DO j=sim%Mesh%JMIN,sim%Mesh%JMAX
141 DO i=sim%Mesh%IMIN,sim%Mesh%IMAX
142 IF (vistype.EQ.beta)
THEN
144 sim%Mesh%radius%bcenter(i,j,k), &
145 pvar%density%data3d(i,j,k),pvar%velocity%data4d(i,j,k,1))
148 sim%Mesh%radius%bcenter(i,j,k), &
149 pvar%density%data3d(i,j,k),pvar%velocity%data4d(i,j,k,1))
157 CALL sim%Error(
"pringle",
"computation of analytical solution requires GSL support")
162 ok = .NOT.sim%aborted
163 ALLOCATE(sigma(sim%Physics%VNUM))
165 IF (
ASSOCIATED(sim%Timedisc%solution))
THEN
166 DO n=1,sim%Physics%VNUM
169 sum_numer = sum(abs(sim%Timedisc%pvar%data4d(sim%Mesh%IMIN:sim%Mesh%IMAX,&
170 sim%Mesh%JMIN:sim%Mesh%JMAX,sim%Mesh%KMIN:sim%Mesh%KMAX,n) &
171 -sim%Timedisc%solution%data4d(sim%Mesh%IMIN:sim%Mesh%IMAX,&
172 sim%Mesh%JMIN:sim%Mesh%JMAX,sim%Mesh%KMIN:sim%Mesh%KMAX,n)))
173 sum_denom = sum(abs(sim%Timedisc%solution%data4d(sim%Mesh%IMIN:sim%Mesh%IMAX,&
174 sim%Mesh%JMIN:sim%Mesh%JMAX,sim%Mesh%KMIN:sim%Mesh%KMAX,n)))
176 IF (sim%GetRank().GT.0)
THEN
177 CALL mpi_reduce(sum_numer,sum_numer,1,default_mpi_real,mpi_sum,0,mpi_comm_world,sim%ierror)
178 CALL mpi_reduce(sum_denom,sum_denom,1,default_mpi_real,mpi_sum,0,mpi_comm_world,sim%ierror)
180 CALL mpi_reduce(mpi_in_place,sum_numer,1,default_mpi_real,mpi_sum,0,mpi_comm_world,sim%ierror)
181 CALL mpi_reduce(mpi_in_place,sum_denom,1,default_mpi_real,mpi_sum,0,mpi_comm_world,sim%ierror)
183 sigma(n) = sum_numer / sum_denom
192 den = sim%Physics%DENSITY
193 vx = sim%Physics%XVELOCITY
194 vy = sim%Physics%YVELOCITY
197 IF (sim%GetRank().EQ.0)
THEN
199tap_check(ok,
"stoptime reached")
202tap_check_small(sigma(den),5.0e-02,
"density deviation < 5%")
203tap_check_small(sigma(vx),9.0e-02,
"radial velocity deviation < 9%")
204tap_check_small(sigma(vy),2.0e-04,
"azimuthal velocity deviation < 0.02%")
211 IF (
ALLOCATED(sigma))
DEALLOCATE(sigma)
226 TYPE(
dict_typ),
POINTER :: mesh, physics, boundary, datafile, &
227 grav, vis, timedisc, fluxes, sources
228 REAL :: x1,x2,y1,y2,z1,z2,cvis
249 bc(bottom) = no_gradients
250 bc(top) = no_gradients
259 bc(west) = no_gradients
260 bc(east) = no_gradients
263 bc(bottom) = no_gradients
264 bc(top) = no_gradients
277 bc(south) = no_gradients
278 bc(north) = no_gradients
279 bc(bottom) = periodic
283 CALL sim%Error(
"pringle::MakeConfig",
"mesh geometry not supported for Pringle disk")
288 "meshtype" / midpoint, &
299 "fargo/method" / 0, &
304 "western" / bc(west), &
305 "eastern" / bc(east), &
306 "southern" / bc(south), &
307 "northern" / bc(north), &
308 "bottomer" / bc(bottom), &
314 "problem" / euler_isotherm, &
315 "units" / geometrical, &
322 "variables" / primitive, &
323 "limiter" / vanleer, &
328 "stype" / viscosity, &
329 "vismodel" / vistype, &
330 "dynconst" / (1./re), &
331 "exponent" / pl_exp, &
337 "pmass/gtype"/ pointmass, &
338 "pmass/mass" / 1.0, &
341 "pmass/outbound" / 0)
351 "stoptime" / (tsim * tvis), &
353 "dtlimit" / (1e-10 * tvis), &
355 "maxiter" / 100000000, &
357 "tol_abs" / (/1e-16,1e-2,1e-6/))
361 IF (mgeo.EQ.cylindrical) &
362 CALL setattr(timedisc,
"output/solution",1)
368 "filename" / (trim(odir) // trim(ofname)), &
373 "physics" / physics, &
374 "boundary" / boundary, &
376 "sources" / sources, &
377 "timedisc" / timedisc, &
378 "datafile" / datafile)
381 SUBROUTINE initdata(Sim,Mesh,Physics,Timedisc,Sources)
386 CLASS(
fosite),
INTENT(INOUT) :: Sim
387 CLASS(physics_base),
INTENT(IN) :: Physics
388 CLASS(mesh_base),
INTENT(IN) :: Mesh
389 CLASS(timedisc_base),
INTENT(INOUT) :: Timedisc
390 CLASS(sources_list),
ALLOCATABLE,
INTENT(IN) :: Sources
396 REAL,
DIMENSION(Mesh%IGMIN:Mesh%IGMAX,Mesh%JGMIN:Mesh%JGMAX,Mesh%KGMIN:Mesh%KGMAX,3) &
398 REAL,
DIMENSION(Mesh%IGMIN:Mesh%IGMAX,Mesh%JGMIN:Mesh%JGMAX,Mesh%KGMIN:Mesh%KGMAX) :: radius
399 REAL :: r,vr,vphi,sden,mu
400 CHARACTER(LEN=64) :: value
404 IF (
ALLOCATED(sources))
THEN
405 sp => sources%GetSourcesPointer(gravity)
406 IF (
ASSOCIATED(sp))
THEN
415 IF (.NOT.
ASSOCIATED(gp)) &
416 CALL physics%Error(
"pringle::InitData",
"no gravity term initialized")
418 IF (abs(x0).LE.tiny(x0).AND.abs(y0).LE.tiny(y0))
THEN
420 radius(:,:,:) = mesh%radius%bcenter(:,:,:)
421 posvec(:,:,:,:) = mesh%posvec%bcenter(:,:,:,:)
428 CALL mesh%Geometry%Convert2Curvilinear(mesh%bcenter,posvec,posvec)
432 posvec(:,:,:,:) = mesh%posvec%bcenter(:,:,:,:) - posvec(:,:,:,:)
434 radius(:,:,:) = sqrt(posvec(:,:,:,1)**2+posvec(:,:,:,2)**2+posvec(:,:,:,3)**2)
442 ephi(:,:,:,1) = -posvec(:,:,:,2)/radius(:,:,:)
443 ephi(:,:,:,2) = posvec(:,:,:,1)/radius(:,:,:)
447 SELECT TYPE(bwest => timedisc%Boundary%boundary(west)%p)
448 CLASS IS (boundary_custom)
449 CALL bwest%SetCustomBoundaries(mesh,physics, &
450 (/custom_nograd,custom_outflow,custom_kepler/))
452 SELECT TYPE(beast => timedisc%Boundary%boundary(east)%p)
453 CLASS IS (boundary_custom)
454 CALL beast%SetCustomBoundaries(mesh,physics, &
455 (/custom_logexpol,custom_outflow,custom_kepler/))
464 CALL timedisc%Error(
"pringle::InitData",
"only pringle and beta viscosity possible")
468 SELECT TYPE(pvar => timedisc%pvar)
469 TYPE IS(statevector_eulerisotherm)
470 DO k=mesh%KGMIN,mesh%KGMAX
471 DO j=mesh%JGMIN,mesh%JGMAX
472 DO i=mesh%IMIN,mesh%IGMAX
485 pvar%density%data3d(i,j,k) = rhomin + sden
487 pvar%velocity%data4d(i,j,k,1:2) = &
488 vr*posvec(i,j,k,1:2)/r + vphi*ephi(i,j,k,1:2)
499 CALL timedisc%Error(
"pringle::InitData",
"only isothermal physics possible")
502 IF (
ASSOCIATED(timedisc%solution))
THEN
504 timedisc%solution = timedisc%pvar
508 SELECT CASE(mesh%fargo%GetType())
510 SELECT CASE(mesh%fargo%GetDirection())
512 timedisc%w(:,:) = sqrt(1.0/radius(mesh%IMIN,:,:))
514 timedisc%w(:,:) = sqrt(1.0/radius(:,mesh%IMIN,:))
516 timedisc%w(:,:) = sqrt(1.0/radius(:,:,mesh%KMIN))
521 CALL physics%Convert2Conservative(timedisc%pvar,timedisc%cvar)
523 CALL mesh%Info(
" DATA-----> initial condition: " //
"viscous spreading ring")
524 WRITE(
value,
"(ES11.3)") tvis
525 CALL mesh%Info(
" " //
"viscous timescale: " //trim(
value))
526 WRITE(
value,
"(ES11.3)") re
527 CALL mesh%Info(
" " //
"Reynolds number: " //trim(
value))
528 WRITE(
value,
"(ES11.3)") ma
529 CALL mesh%Info(
" " //
"Mach number: " //trim(
value))
537 REAL,
INTENT(IN) :: mu,tau,r
538 REAL,
INTENT(OUT) :: Sigma,vr
540 REAL :: lambda,x,r4l,invr4l
543 r4l = r**(0.25/lambda)
544 x = 2*lambda**2 / tau * r**(0.25/lambda)
546 sigma = mdisk/sqrt((4*pi)**3 * tau) * r4l**(1.5-9*lambda) &
547 * exp(-lambda**2/tau*(1.0-r4l)**2)
548 vr = -1.5/re * ((x*tau)/(2*lambda*lambda))**(4*lambda-2.) &
549 * (1.0 - 0.5*x/lambda * (x*tau/(2*lambda*lambda) - 1.0))
556 REAL,
INTENT(IN) :: mu,tau,r
557 REAL,
INTENT(OUT) :: Sigma,vr
562 r4l = r**(0.25/lambda)
563 x = 2*lambda**2 / tau * r4l
564 sigma = mdisk/(4*pi) * (lambda/tau) * r4l**2 * r**(-2.25) &
565 * exp(-lambda**2/tau*(1.0+r4l**2)) &
567 vr = -1.5/re * ((x*tau)/(2*lambda*lambda))**(4*lambda-2.) &
568 * (1.0 - 0.5*x/lambda * (x*tau/(2*lambda*lambda) &
576 REAL,
INTENT(IN) :: nu,x
581 IF (aint(nu).EQ.nu)
THEN
584 inu = fgsl_sf_bessel_icn(n,x)
588 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.
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
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)