aboutsummaryrefslogtreecommitdiff
path: root/interface.ccl
diff options
context:
space:
mode:
authorrhaas <rhaas@ac85fae7-cede-4708-beff-ae01c7fa1c26>2012-09-04 12:28:40 +0000
committerrhaas <rhaas@ac85fae7-cede-4708-beff-ae01c7fa1c26>2012-09-04 12:28:40 +0000
commit91270e519a2413c00e39aee4688765bd0e4d9ab1 (patch)
tree9af5d4df8ca80a6d53c2a810ef360de8696f55c3 /interface.ccl
parent0f066b0ae931fffa90bc79303072891357b4f333 (diff)
Changes to make GRHydro MHD work with nuclear/hot equation of state:
1) Switch to EOSOmni pointwise C2P routine and modify where necessary. 2) Modify Con2PrimM.F90 to allow for the evolution of temperature and adjust the wrapper routine. 3) Create EigenProblemM_hot pointwise routine and call that from HLLEM.F90 when temperature is evolved. Additionally adjust HLLEM where necessary. 4) Adjust InterfacesM.h to incorporate the newly created functions. 5) Fix a loop problem (not taking into account constraint transport) in PPM reconstruction of Y_e 6) Introduce Prim2ConM_hot and call this pointwise routine from Prim2ConM.F90 when temperature is evolved. Additionally also make this routine available to initial data routine in GRHydro_InitData 7) Adjust loops in GRHydro_PoloidalMagFieldM.F90 to not set boundary points it cannot set but instead call boundary group afterwards! Pay attention as this will not work with boundary conditions set to "none" in MHD case anymore but is the correct thing to do. 8) Allow StarMapper to extend HydroBase::initial_hydro = "starmapper". 9) Smaller fixes. From: Philipp Moesta <pmoesta@tapir.caltech.edu> git-svn-id: http://svn.einsteintoolkit.org/cactus/EinsteinInitialData/GRHydro_InitData/trunk@153 ac85fae7-cede-4708-beff-ae01c7fa1c26
Diffstat (limited to 'interface.ccl')
-rw-r--r--interface.ccl17
1 files changed, 17 insertions, 0 deletions
diff --git a/interface.ccl b/interface.ccl
index 98fcd8d..ec711c4 100644
--- a/interface.ccl
+++ b/interface.ccl
@@ -66,6 +66,22 @@ SUBROUTINE Prim2ConGenM(CCTK_INT IN handle, \
CCTK_REAL OUT press, CCTK_REAL IN Bvecx, CCTK_REAL IN Bvecy, \
CCTK_REAL IN Bvecz, CCTK_REAL OUT w_lorentz)
+SUBROUTINE Prim2ConGenM_hot(CCTK_INT IN handle, CCTK_INT IN GRHydro_reflevel, CCTK_INT IN i, CCTK_INT IN j, CCTK_INT IN k, \
+ CCTK_REAL IN x, CCTK_REAL IN y, CCTK_REAL IN z, \
+ CCTK_REAL IN gxx, CCTK_REAL IN gxy, \
+ CCTK_REAL IN gxz, CCTK_REAL IN gyy, \
+ CCTK_REAL IN gyz, CCTK_REAL IN gzz, \
+ CCTK_REAL IN det, CCTK_REAL OUT dens, \
+ CCTK_REAL OUT sx, CCTK_REAL OUT sy, \
+ CCTK_REAL OUT sz, CCTK_REAL OUT tau, \
+ CCTK_REAL OUT Bconsx, CCTK_REAL OUT Bconsy, \
+ CCTK_REAL OUT Bconsz, CCTK_REAL IN rho, CCTK_REAL IN velx, \
+ CCTK_REAL IN vely, \
+ CCTK_REAL IN velz, CCTK_REAL IN epsilon, \
+ CCTK_REAL OUT press, CCTK_REAL IN Bvecx, CCTK_REAL IN Bvecy, \
+ CCTK_REAL IN Bvecz, CCTK_REAL OUT w_lorentz, \
+ CCTK_REAL INOUT temperature, CCTK_REAL IN Y_e)
+
SUBROUTINE Con2PrimPoly(CCTK_INT IN handle, CCTK_REAL OUT dens, \
CCTK_REAL OUT sx, CCTK_REAL OUT sy, \
@@ -154,6 +170,7 @@ USES FUNCTION Prim2ConPoly
USES FUNCTION Prim2ConGen
USES FUNCTION Prim2ConPolyM
USES FUNCTION Prim2ConGenM
+USES FUNCTION Prim2ConGenM_hot
USES FUNCTION Con2PrimPoly
USES FUNCTION Con2PrimGen
USES FUNCTION Con2PrimGenM