59 REAL,
PARAMETER :: gn = 6.6742d-11
60 REAL,
PARAMETER :: cc = 2.99792458d+8
61 REAL,
PARAMETER :: msun = 1.989d+30
62 REAL,
PARAMETER :: au = 1.49597870691e+11
63 REAL,
PARAMETER :: parsec = 0.5/pi*1.296e+6 * au
64 REAL,
PARAMETER :: year = 3.15576e+7
65 REAL,
PARAMETER :: eta_acc = 0.1
67 REAL,
PARAMETER :: tsim = 1.0e+4*year
68 REAL,
PARAMETER :: m_bh = 1.0e+6 * msun
69 REAL,
PARAMETER :: m_disk = 1.0e+0 * m_bh
70 REAL,
PARAMETER :: rs = 2 * gn * m_bh / cc**2
71 REAL,
PARAMETER :: lscale = 1.0e+3 * rs
72 REAL,
PARAMETER :: etasgs = 1.0e-6
73 REAL,
PARAMETER :: gamma = 1.4
75 REAL,
PARAMETER :: rmin = 5.0d+0 * lscale
76 REAL,
PARAMETER :: rmax = 3.0d+2 * lscale
77 REAL,
PARAMETER :: zmin = -1.5d+2 * lscale
78 REAL,
PARAMETER :: zmax = 1.5d+2 * lscale
79 REAL,
PARAMETER :: ztan = 1.0e+2 * lscale
85 REAL,
PARAMETER :: z0 = 3.0e+02 * lscale
86 REAL,
PARAMETER :: a0 = 1.0e-01
87 REAL,
PARAMETER :: s0 = 1.5e+02 * lscale
88 REAL,
PARAMETER :: b0 = -1.0
89 REAL,
PARAMETER :: rho_inf = 1.0e-15
90 REAL,
PARAMETER :: p_inf = 1.4e-10
91 INTEGER,
PARAMETER :: phys = euler
93 INTEGER,
PARAMETER :: mgeo = cylindrical
94 INTEGER,
PARAMETER :: rres = 100
95 INTEGER,
PARAMETER :: phires = 1
96 INTEGER,
PARAMETER :: zres = 101
98 INTEGER,
PARAMETER :: onum = 1000
99 CHARACTER(LEN=256),
PARAMETER :: odir &
101 CHARACTER(LEN=256),
PARAMETER &
102 :: ofname =
'agndisk'
104 CLASS(
fosite),
ALLOCATABLE :: sim
109 CALL sim%InitFosite()
116 CALL initdata(sim%Mesh, sim%Physics, sim%Fluxes, sim%Timedisc,sim%Sources)
127 TYPE(dict_typ),
POINTER :: config
130 TYPE(dict_typ),
POINTER :: mesh, physics, boundary, datafile, logfile, sources, &
131 timedisc,fluxes,stemp, grav
137 "meshtype" / midpoint, &
153 "western" / custom, &
158 "eastern" / farfield, &
159 "bottomer" / farfield, &
160 "topper" / farfield, &
162 "southern" / periodic, &
164 "northern" / periodic)
176 "meshtype" / midpoint, &
192 "western" / custom, &
197 "eastern" / farfield, &
198 "bottomer" / periodic, &
199 "topper" / periodic, &
228 "variables" / primitive, &
229 "limiter" / minmod, &
233 grav => dict(
"stype" / gravity, &
234 "output/accel" / 1, &
235 "pmass/gtype" / pointmass, &
236 "pmass/potential" / newton, &
237 "pmass/outbound" / 1, &
268 "dtlimit" / 1.0e-3, &
269 "maxiter" / 1000000000, &
270 "output/density" / 1, &
271 "output/xvelocity" / 1, &
272 "output/yvelocity" / 1, &
273 "output/zvelocity" / 1, &
274 "output/pressure" / 1, &
275 "output/fluxes" / 1, &
277 "output/sgspressure" / 1, &
278 "output/geometrical_sources" / 1, &
279 "output/external_sources" / 1)
293 "fileformat" / vtk, &
296 "filename" / (trim(odir) // trim(ofname)), &
299 "filecycles" / (onum+1))
301 config => dict(
"mesh" / mesh, &
302 "physics" / physics, &
303 "boundary" / boundary, &
305 "sources" / sources, &
306 "timedisc" / timedisc, &
308 "datafile" / datafile)
312 SUBROUTINE initdata(Mesh,Physics,Fluxes,Timedisc,Sources)
317 CLASS(mesh_base),
INTENT(IN) :: Mesh
318 CLASS(physics_base),
INTENT(INOUT) :: Physics
319 CLASS(fluxes_base),
INTENT(INOUT) :: Fluxes
320 CLASS(timedisc_base),
INTENT(INOUT) :: Timedisc
321 CLASS(sources_list),
ALLOCATABLE,
INTENT(INOUT) :: Sources
326 REAL,
DIMENSION(Mesh%IGMIN:Mesh%IGMAX,Mesh%JGMIN:Mesh%JGMAX,2) :: dv
327 REAL,
DIMENSION(Mesh%IGMIN:Mesh%IGMAX,Mesh%JGMIN:Mesh%JGMAX) :: velo_periodic
328 REAL :: s,z,tau,lambda,tau1,tau2,dPs,dPz,s1,s2,z1,z2,P0,rho0,cs_inf
332 REAL :: mdisk_local,am_local
335 CHARACTER(LEN=20) :: mdisk_str,am_str
338 IF (
ALLOCATED(sources)) &
339 sp => sources%GetSourcesPointer(gravity)
340 IF (.NOT.
ASSOCIATED(sp)) &
341 CALL physics%Error(
"agndisk::InitData",
"no gravity term initialized")
347 timedisc%pvar%data2d(:,physics%XVELOCITY) = 0.0
348 timedisc%pvar%data2d(:,physics%ZVELOCITY) = 0.0
351 cs_inf = sqrt(gamma*p_inf/rho_inf)
353 DO k=mesh%KGMIN,mesh%KGMAX
354 DO j=mesh%JGMIN,mesh%JGMAX
355 DO i=mesh%IGMIN,mesh%IGMAX
356 s = abs(mesh%bccart(i,j,k,1))
357 z = mesh%bccart(i,j,k,3)
362 timedisc%pvar%data4d(i,j,k,physics%PRESSURE) =
pc(tau,rho0,cs_inf,gamma)
363 z1 = mesh%cart%faces(i,j,k,5,3)
364 z2 = mesh%cart%faces(i,j,k,6,3)
367 dpz = (
pc(tau2,rho0,cs_inf,gamma)-
pc(tau1,rho0,cs_inf,gamma))/mesh%dlz%data3d(i,j,k)
368 IF (abs(z).GT.0.0)
THEN
370 timedisc%pvar%data4d(i,j,k,physics%DENSITY) = dpz /
gz(s,z)
372 timedisc%pvar%data4d(i,j,k,physics%DENSITY) = rho0 * (s/s0)**b0
374 s1 = abs(mesh%cart%faces(i,j,k,1,1))
375 s2 = mesh%cart%faces(i,j,k,2,1)
378if (tau1 .eq. 0.0 .or. tau2 .eq. 0.0) print *,tau1,tau2,i,j
379 dps = (
pc(tau2,rho0,cs_inf,gamma)-
pc(tau1,rho0,cs_inf,gamma))/mesh%dlx%data3d(i,j,k)
382 timedisc%pvar%data4d(i,j,k,physics%YVELOCITY) = sqrt(max(0.0,-s*
gs(s,z) / lambda))
391 mdisk = sum(timedisc%pvar%data4d(mesh%IMIN:mesh%IMAX,mesh%JMIN:mesh%JMAX,mesh%KMIN:mesh%KMAX,physics%DENSITY) &
392 * mesh%volume%data3d(mesh%IMIN:mesh%IMAX,mesh%JMIN:mesh%JMAX,mesh%KMIN:mesh%KMAX))
394 am = sum(mesh%hy%bcenter(mesh%IMIN:mesh%IMAX,mesh%JMIN:mesh%JMAX,mesh%KMIN:mesh%KMAX) &
395 * timedisc%pvar%data4d(mesh%IMIN:mesh%IMAX,mesh%JMIN:mesh%JMAX,mesh%KMIN:mesh%KMAX,physics%YVELOCITY) &
396 * timedisc%pvar%data4d(mesh%IMIN:mesh%IMAX,mesh%JMIN:mesh%JMAX,mesh%KMIN:mesh%KMAX,physics%DENSITY) &
397 * mesh%volume%data3d(mesh%IMIN:mesh%IMAX,mesh%JMIN:mesh%JMAX,mesh%KMIN:mesh%KMAX))
401 CALL mpi_allreduce(mdisk_local,mdisk,1,default_mpi_real,mpi_sum, &
402 mesh%comm_cart,ierror)
404 CALL mpi_allreduce(am_local,am,1,default_mpi_real,mpi_sum, &
405 mesh%comm_cart,ierror)
410 rho0 = m_disk / mdisk * rho0
423 SELECT TYPE(beast => timedisc%Boundary%Boundary(east)%p)
424 CLASS IS (boundary_farfield)
425 DO k=mesh%KMIN,mesh%KMAX
426 DO j=mesh%JMIN,mesh%JMAX
428 beast%data(i,j,k,:) = timedisc%pvar%data4d(mesh%IMAX+i,j,k,:)
435 SELECT TYPE(bwest => timedisc%Boundary%Boundary(west)%p)
436 CLASS IS(boundary_custom)
437 CALL bwest%SetCustomBoundaries(mesh,physics, &
438 (/custom_nograd,custom_outflow,custom_kepler,custom_nograd,custom_nograd/))
442 SELECT TYPE (bbottom => timedisc%Boundary%Boundary(bottom)%p)
443 CLASS IS(boundary_farfield)
445 DO j=mesh%JMIN,mesh%JMAX
446 DO i=mesh%IMIN,mesh%IMAX
447 bbottom%data(i,j,k,:) = timedisc%pvar%data4d(i,j,mesh%KMIN-k,:)
454 SELECT TYPE (btop => timedisc%Boundary%Boundary(top)%p)
455 CLASS IS(boundary_farfield)
457 DO j=mesh%JMIN,mesh%JMAX
458 DO i=mesh%IMIN,mesh%IMAX
459 btop%data(i,j,k,:) = timedisc%pvar%data4d(i,j,mesh%KMAX+k,:)
465 CALL physics%Convert2Conservative(timedisc%pvar,timedisc%cvar)
467 CALL timedisc%Boundary%CenterBoundary(mesh,physics,0.0,timedisc%pvar, &
471 WRITE (mdisk_str,
'(ES10.2)') mdisk/msun
472 WRITE (am_str,
'(ES10.2)') am
473 CALL mesh%Info(
" DATA-----> initial condition: non Keplerian flow")
474 CALL mesh%Info(
" disk mass: " // trim(mdisk_str) //
" M_sun")
475 CALL mesh%Info(
" angular momentum: " // trim(am_str) //
" kg/m^2/s")
481 REAL,
INTENT(IN) :: sc,zc
483 res = sqrt(gn*m_bh/sc) * (1.0d0+(zc/sc)**2)**(-0.75)
486 ELEMENTAL FUNCTION gs(sc,zc)
RESULT(res)
488 REAL,
INTENT(IN) :: sc,zc
490 res = -(gn*m_bh/sc**2) / sqrt(1.0d0+(zc/sc)**2)**3
493 ELEMENTAL FUNCTION gz(sc,zc)
RESULT(res)
495 REAL,
INTENT(IN) :: sc,zc
497 res = -(gn*m_bh/sc**2) * (zc/sc) / sqrt(1.0d0+(zc/sc)**2)**3
502 REAL,
INTENT(IN) :: zc
504 res = 1.0d0 + a0/(1.0d0+(zc/z0)**2)
509 REAL,
INTENT(IN) :: sc,zc
511 res = sc*sqrt(1.0d0 + (zc/sc)**2*(1.0d0 + (1.0d0+0.5d0*(zc/z0)**2)/a0))
514 ELEMENTAL FUNCTION pc(tau,rho_s0,cs_inf,gamma)
RESULT(res)
516 REAL,
INTENT(IN) :: tau,rho_s0,cs_inf,gamma
519 res = p_inf * (1.0 + a0/(1.0+a0) * 1.0/(1.0-b0) * 0.5 * gamma * &
520 (cc/cs_inf)**2 * rho_s0/rho_inf * (rs/s0) * (tau/s0)**(b0-1.0))
elemental real function func_tau(sc, zc)
elemental real function pc(tau, rho_s0, cs_inf, gamma)
elemental real function gz(sc, zc)
elemental real function func_lambda(zc)
elemental real function vkepler(sc, zc)
elemental real function gs(sc, zc)
subroutine initdata(Mesh, Physics, Fluxes, Timedisc)
subroutine makeconfig(Sim, config)
generic source terms module providing functionaly common to all source terms
generic gravity terms module providing functionaly common to all gravity terms