pringle3d.f90
Go to the documentation of this file.
1!#############################################################################
2!# #
3!# fosite - 3D hydrodynamical simulation program #
4!# module: pringle3d.f90 #
5!# #
6!# Copyright (C) 2008-2019 #
7!# Bjoern Sperling <sperling@astrophysik.uni-kiel.de> #
8!# Tobias Illenseer <tillense@astrophysik.uni-kiel.de> #
9!# #
10!# This program is free software; you can redistribute it and/or modify #
11!# it under the terms of the GNU General Public License as published by #
12!# the Free Software Foundation; either version 2 of the License, or (at #
13!# your option) any later version. #
14!# #
15!# This program is distributed in the hope that it will be useful, but #
16!# WITHOUT ANY WARRANTY; without even the implied warranty of #
17!# MERCHANTABILITY OR FITNESS FOR A PARTICULAR PURPOSE, GOOD TITLE or #
18!# NON INFRINGEMENT. See the GNU General Public License for more #
19!# details. #
20!# #
21!# You should have received a copy of the GNU General Public License #
22!# along with this program; if not, write to the Free Software #
23!# Foundation, Inc., 675 Mass Ave, Cambridge, MA 02139, USA. #
24!# #
25!#############################################################################
26
27!----------------------------------------------------------------------------!
60!----------------------------------------------------------------------------!
61PROGRAM pringle3d
62 USE fosite_mod
63 USE common_dict
64 IMPLICIT NONE
65 !--------------------------------------------------------------------------!
66 ! simulation parameters
67 REAL, PARAMETER :: tsim = 5.0e-3 ! simulation time [TVIS] see below
68 ! should be < 1 ~ viscous time scale
69 ! these are the basic parameters; for stable solutions one requires
70 ! MA >> 1 and RE > MA
71 REAL, PARAMETER :: re = 1.0e+4 ! Reynolds number (at R=1)
72 REAL, PARAMETER :: ma = 1.0e+2 ! Mach number (at R=1)
73 REAL, PARAMETER :: mdisk = 1.0e-2 ! disk mass << 1 = central mass
74 ! lower limit for nitial density
75 REAL, PARAMETER :: rhomin = 1.0e-8 ! minimal initial density
76 ! viscosity prescription
77 INTEGER, PARAMETER :: vistype = beta
78! INTEGER, PARAMETER :: VISTYPE = POWERLAW
79 REAL, PARAMETER :: pl_exp = 0.0 ! exponent for power law viscosity
80 REAL, PARAMETER :: tvis = 4./3.*re ! viscous time scale
81 REAL, PARAMETER :: tinit = 1.0e-3 ! time for initial condition [TVIS]
82 ! mesh settings
83! INTEGER, PARAMETER :: MGEO = CYLINDRICAL
84! INTEGER, PARAMETER :: MGEO = LOGCYLINDRICAL
85 INTEGER, PARAMETER :: mgeo = spherical !!! ATTENTION: not applicable in 1D
86 INTEGER, PARAMETER :: xres = 200 ! x-resolution
87 INTEGER, PARAMETER :: yres = 16 ! y-resolution
88 INTEGER, PARAMETER :: zres = 6 ! z-resolution
89 REAL, PARAMETER :: rmin = 0.1 ! min radius of comp. domain
90 REAL, PARAMETER :: rmax = 2.0 ! max radius of comp. domain
91 REAL, PARAMETER :: gpar = 1.0 ! geometry scaling parameter
92 ! output parameters
93 INTEGER, PARAMETER :: onum = 100 ! number of output data sets
94 CHARACTER(LEN=256), PARAMETER & ! output data dir
95 :: odir = './'
96 CHARACTER(LEN=256), PARAMETER & ! output data file name
97 :: ofname = 'pringle3d'
98 ! derived parameters
99 REAL, PARAMETER :: csiso = 1./ma ! isothermal speed of sound
100 !--------------------------------------------------------------------------!
101 CLASS(fosite), ALLOCATABLE :: sim
102 !--------------------------------------------------------------------------!
103
104 ALLOCATE(sim)
105 CALL sim%InitFosite()
106 CALL makeconfig(sim, sim%config)
107 CALL sim%Setup()
108 CALL initdata(sim,sim%Mesh, sim%Physics, sim%Timedisc)
109 CALL sim%Run()
110 CALL sim%Finalize()
111 DEALLOCATE(sim)
112
113CONTAINS
114
115 SUBROUTINE makeconfig(Sim, config)
116 USE functions, ONLY : asinh
117 IMPLICIT NONE
118 !------------------------------------------------------------------------!
119 CLASS(fosite) :: Sim
120 TYPE(dict_typ),POINTER :: config
121 !------------------------------------------------------------------------!
122 ! Local variable declaration
123 INTEGER :: bc(6)
124 TYPE(dict_typ), POINTER :: mesh, physics, boundary, datafile, &
125 grav, vis, timedisc, fluxes, sources
126 REAL :: x1,x2,y1,y2,z1,z2,cvis
127 !------------------------------------------------------------------------!
128 INTENT(INOUT) :: sim
129 !------------------------------------------------------------------------!
130 ! mesh settings and boundary conditions
131 SELECT CASE(mgeo)
132 CASE(cylindrical)
133 x1 = rmin
134 x2 = rmax
135 y1 = -pi
136 y2 = pi
137! z1 = 0.0
138! z2 = 0.0
139 z1 = -0.1
140 z2 = 0.1
141! bc(WEST) = AXIS
142 bc(west) = no_gradients
143 bc(east) = no_gradients
144! bc(WEST) = ABSORBING
145! bc(EAST) = ABSORBING
146! bc(WEST) = CUSTOM
147! bc(EAST) = CUSTOM
148 bc(south) = periodic
149 bc(north) = periodic
150! bc(BOTTOM) = ABSORBING
151! bc(TOP) = ABSORBING
152 bc(bottom) = reflecting
153 bc(top) = reflecting
154 !!! ATTENTION: empirical formula found for standard parameters
155 !!! seems to work for RMIN=0.05 and XRES <= 1000
156 cvis = (xres/1207.0)**1.9
157 CASE(logcylindrical)
158 x1 = log(rmin/gpar)
159 x2 = log(rmax/gpar)
160 y1 = -pi
161 y2 = pi
162 z1 = 0.0
163 z2 = 0.0
164 bc(west) = no_gradients
165 bc(east) = no_gradients
166 bc(south) = periodic
167 bc(north) = periodic
168 bc(bottom) = no_gradients
169 bc(top) = no_gradients
170 cvis = 0.1
171 CASE(spherical)
172 x1 = rmin
173 x2 = rmax
174 y1 = pi/2 - 0.02*pi
175 y2 = pi/2 + 0.02*pi
176 z1 = -pi
177 z2 = pi
178! bc(WEST) = ABSORBING
179 bc(east) = absorbing
180! bc(WEST) = NO_GRADIENTS
181! bc(EAST) = NO_GRADIENTS
182 bc(west) = custom
183! bc(EAST) = CUSTOM
184! bc(SOUTH) = NO_GRADIENTS
185! bc(NORTH) = NO_GRADIENTS
186 bc(south) = absorbing
187 bc(north) = absorbing
188 bc(bottom) = periodic
189 bc(top) = periodic
190 !!! ATTENTION: empirical formula found for standard parameters
191 !!! seems to work for RMIN=0.05 and XRES <= 1000
192 cvis = (xres/1207.0)**1.9
193 CASE DEFAULT
194 CALL sim%Error("pringle::MakeConfig","mesh geometry not supported for Pringle disk")
195 END SELECT
196
197 ! mesh settings
198 mesh => dict( &
199 "meshtype" / midpoint, &
200 "geometry" / mgeo, &
201 "inum" / xres, &
202 "jnum" / yres, &
203 "knum" / zres, &
204 "xmin" / x1, &
205 "xmax" / x2, &
206 "ymin" / y1, &
207 "ymax" / y2, &
208 "zmin" / z1, &
209 "zmax" / z2, &
210! "use_fargo" / 1, &
211! "fargo" / 2, &
212 "gparam" / gpar)
213
214 ! boundary conditions
215 boundary => dict( &
216 "western" / bc(west), &
217 "eastern" / bc(east), &
218 "southern" / bc(south), &
219 "northern" / bc(north), &
220 "bottomer" / bc(bottom), &
221 "topper" / bc(top))
222
223
224 ! physics settings
225 physics => dict( &
226 "problem" / euler_isotherm, &
227 "units" / geometrical, &
228 "cs" / csiso)
229
230 ! flux calculation and reconstruction method
231 fluxes => dict( &
232 "order" / linear, &
233 "fluxtype" / kt, &
234 "variables" / primitive, &
235 "limiter" / vanleer, &
236 "theta" / 1.2)
237
238 ! source term due to viscosity
239 vis => dict( &
240 "stype" / viscosity, &
241 "vismodel" / vistype, &
242 "dynconst" / (1./re), &
243 "exponent" / pl_exp, &
244 "cvis" / cvis)
245
246 ! source term due to a point mass
247 grav => dict( &
248 "stype" / gravity, &
249 "pmass/gtype"/ pointmass, & ! grav. accel. of a point mass
250 "pmass/mass" / 1.0, & ! non-dim. mass of the accreting object
251 "pmass/outbound" / 0) ! disable accretion
252
253 ! combine all source terms
254 sources => dict( &
255 "viscosity" / vis, &
256 "grav" / grav)
257
258 ! time discretization settings
259 timedisc => dict( &
260 "method" / ssprk, &
261 "stoptime" / (tsim * tvis), &
262 "cfl" / 0.4, &
263 "dtlimit" / (1e-10 * tvis), &
264! "rhstype" / 1, &
265 "maxiter" / 100000000, &
266 "tol_rel" / 0.01, &
267 "tol_abs" / (/1e-40,1e-5,1e-8,1e-5/))
268
269! ! enable output of analytical solution
270! #ifdef HAVE_FGSL
271! IF (MGEO.EQ.CYLINDRICAL) &
272! CALL SetAttr(timedisc,"output/solution",1)
273! #endif
274
275 ! initialize data input/output
276 datafile => dict( &
277 "fileformat"/ vtk, &
278 "filename" / (trim(odir) // trim(ofname)), &
279 "count" / onum)
280
281 config => dict( &
282 "mesh" / mesh, &
283 "physics" / physics, &
284 "boundary" / boundary, &
285 "fluxes" / fluxes, &
286 "sources" / sources, &
287 "timedisc" / timedisc, &
288 "datafile" / datafile)
289 END SUBROUTINE makeconfig
290
291 SUBROUTINE initdata(Sim,Mesh,Physics,Timedisc)
295 IMPLICIT NONE
296 !------------------------------------------------------------------------!
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
301 !------------------------------------------------------------------------!
302 ! Local variable declaration
303 CLASS(sources_base), POINTER :: sp
304 CLASS(sources_gravity), POINTER :: gp => null()
305 CLASS(geometry_base), ALLOCATABLE :: geo_cyl
306 TYPE(dict_typ), POINTER :: geo_config
307 INTEGER :: i,j,k
308#ifdef PARALLEL
309 INTEGER :: ierror
310#endif
311 REAL, DIMENSION(Mesh%IGMIN:Mesh%IGMAX,Mesh%JGMIN:Mesh%JGMAX,Mesh%KGMIN:Mesh%KGMAX,3) &
312 :: bccyl
313 REAL :: r,z,vr,vphi,sden,mu,height,mass
314 CHARACTER(LEN=64) :: value
315 !------------------------------------------------------------------------!
316 geo_config => dict("geometry" / cylindrical)
317 CALL new_geometry(geo_cyl,geo_config)
318 CALL deletedict(geo_config)
319
320 ! compute cylindrical coordinates for each cell bary center
321 CALL geo_cyl%Convert2Curvilinear(mesh%bccart,bccyl)
322
323 ! custom boundary conditions if requested
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/))
328 END SELECT
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/))
333 END SELECT
334
335 ! get gravitational acceleration
336 sp => sim%Sources%GetSourcesPointer(gravity)
337 IF (.NOT.ASSOCIATED(sp)) CALL physics%Error("pringle::InitData","no gravity term initialized")
338
339 SELECT CASE(vistype)
340 CASE(beta)
341 mu=1.0
342 CASE(powerlaw)
343 mu=pl_exp
344 CASE DEFAULT
345 CALL timedisc%Error("pringle::InitData","only pringle and beta viscosity possible")
346 END SELECT
347
348 ! initial condition
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
354 ! distance to axis (r-coordinate in cylindrical coordinates)
355 r = bccyl(i,j,k,1)
356 ! z-coordinate
357 z = bccyl(i,j,k,3)
358 ! Keplerian velocity
359 vphi = r/(r**2+z**2)**0.75
360 ! compute pressure scale height
361 height = csiso * sqrt(r**3) ! h = cs / Omega = cs / SQRT(GM/r**3), with GM=1
362#ifdef HAVE_FGSL
363 ! use exact solution at time TINIT
364 CALL viscousring_analytic(mu,tinit,r,sden,vr)
365#else
366 ! use approximation for TINIT << 1
367 ! no need for modified Bessel functions
368 CALL viscousring_approx(mu,tinit,r,sden,vr)
369#endif
370 pvar%density%data3d(i,j,k) = rhomin + sden / (sqrt(2*pi) * height) &
371 * exp(-0.5*(z/height)**2)
372! pvar%velocity%data4d(i,j,k,1) = vr
373! pvar%velocity%data4d(i,j,k,2) = vphi
374! pvar%velocity%data4d(i,j,k,3) = 0.0
375 END DO
376 END DO
377 END DO
378 ! determine disk mass
379 mass = sum(mesh%volume%data1d(:)*pvar%density%data1d(:),mask=mesh%without_ghost_zones%mask1d)
380#ifdef PARALLEL
381 CALL mpi_allreduce(mpi_in_place,mass,1,default_mpi_real,mpi_sum, &
382 mesh%comm_cart,ierror)
383#endif
384 ! rescale disk mass
385 pvar%density%data1d(:) = pvar%density%data1d(:) * mdisk / mass
386 pvar%velocity%data1d(:) = 0.0 ! initialize all velocities with zero
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./), &
390 gp%accel%data4d)
391 CLASS DEFAULT
392 CALL timedisc%Error("pringle::InitData","only isothermal physics possible")
393 END SELECT
394
395 IF (ASSOCIATED(timedisc%solution)) THEN
396 ! store initial condition in exact solution array
397 timedisc%solution = timedisc%pvar
398 END IF
399
400 ! transform to conservative variables
401 CALL physics%Convert2Conservative(timedisc%pvar,timedisc%cvar)
402
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))
412
413 CALL geo_cyl%Finalize()
414 DEALLOCATE(geo_cyl)
415 END SUBROUTINE initdata
416
417
418 SUBROUTINE viscousring_approx(mu,tau,r,Sigma,vr)
419 IMPLICIT NONE
420 !--------------------------------------------------------------------!
421 REAL, INTENT(IN) :: mu,tau,r
422 REAL, INTENT(OUT) :: Sigma,vr
423 !--------------------------------------------------------------------!
424 REAL :: lambda,x,r4l,invr4l
425 !--------------------------------------------------------------------!
426 lambda = 1./(4.0-mu)
427 r4l = r**(0.25/lambda)
428 x = 2*lambda**2 / tau * r**(0.25/lambda)
429 invr4l = 1./r4l
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))
434 END SUBROUTINE viscousring_approx
435
436#ifdef HAVE_FGSL
437 SUBROUTINE viscousring_analytic(mu,tau,r,Sigma,vr)
438 IMPLICIT NONE
439 !--------------------------------------------------------------------!
440 REAL, INTENT(IN) :: mu,tau,r
441 REAL, INTENT(OUT) :: Sigma,vr
442 !--------------------------------------------------------------------!
443 REAL :: lambda,r4l,x
444 !--------------------------------------------------------------------!
445 lambda = 1./(4.0-mu)
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)) &
450 * modified_bessel_inu(lambda,x)
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) &
453 - modified_bessel_inu(1+lambda,x)/modified_bessel_inu(lambda,x)))
454 END SUBROUTINE viscousring_analytic
455
456 FUNCTION modified_bessel_inu(nu,x) RESULT(Inu)
457 USE fgsl
458 IMPLICIT NONE
459 !--------------------------------------------------------------------!
460 REAL, INTENT(IN) :: nu,x
461 REAL :: inu
462 !--------------------------------------------------------------------!
463 INTEGER :: n
464 !--------------------------------------------------------------------!
465 IF (aint(nu).EQ.nu) THEN
466 ! nu is an integer -> use modified bessel functions of integer order
467 n = int(nu)
468 inu = fgsl_sf_bessel_icn(n,x)
469 ELSE
470 ! nu is not an integer -> use representation of modified bessel functions
471 ! through hypergeometric functions 0F1
472 inu = (0.5*x)**nu / fgsl_sf_gamma(1.0+nu) * fgsl_sf_hyperg_0f1(1.0+nu,(0.5*x)**2)
473 END IF
474 END FUNCTION modified_bessel_inu
475#endif
476
477END PROGRAM pringle3d
subroutine initdata(Mesh, Physics, Fluxes, Timedisc)
Definition: bondi2d.f90:274
subroutine makeconfig(Sim, config)
Definition: bondi2d.f90:165
Dictionary for generic data types.
Definition: common_dict.f90:61
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...
Mathematical functions.
Definition: functions.f90:48
real, parameter pi
Definition: functions.f90:52
elemental real function, public asinh(x)
inverse hyperbolic sine function
Definition: functions.f90:100
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
program pringle3d
Definition: pringle3d.f90:61
subroutine viscousring_analytic(mu, tau, r, Sigma, vr)
Definition: pringle.f90:554
real function modified_bessel_inu(nu, x)
Definition: pringle.f90:573
subroutine viscousring_approx(mu, tau, r, Sigma, vr)
Definition: pringle.f90:535
main fosite class
Definition: fosite.f90:71