aboutsummaryrefslogtreecommitdiff
diff options
context:
space:
mode:
authorrhaas <rhaas@ac85fae7-cede-4708-beff-ae01c7fa1c26>2013-01-11 15:04:13 +0000
committerrhaas <rhaas@ac85fae7-cede-4708-beff-ae01c7fa1c26>2013-01-11 15:04:13 +0000
commit7c39d4964e8b048681a402347b2dc7c4d57f8fee (patch)
tree5b2958e2f36fdc0366e27dcd56c0e05e3ff1340d
parentc5d0dacc66b08443898399c0eb14ac973bbce544 (diff)
GRHydro_InitData: Add basic vector potential support
Initial Avec constrained to poloidal at the moment. From: Tanja Bode <tanja.bode@physics.gatech.edu> git-svn-id: http://svn.einsteintoolkit.org/cactus/EinsteinInitialData/GRHydro_InitData/trunk@195 ac85fae7-cede-4708-beff-ae01c7fa1c26
-rw-r--r--param.ccl5
-rw-r--r--schedule.ccl2
-rw-r--r--src/GRHydro_PoloidalMagFieldM.F9018
3 files changed, 23 insertions, 2 deletions
diff --git a/param.ccl b/param.ccl
index 276cbb5..51c4c9d 100644
--- a/param.ccl
+++ b/param.ccl
@@ -35,6 +35,11 @@ EXTENDS KEYWORD initial_Bvec
"magnetized Bondi" :: "radial magnetic field appropriate for Bondi test"
}
+EXTENDS KEYWORD initial_Avec
+{
+ "poloidalmagfield" :: "Poloidal Magnetic Field"
+}
+
EXTENDS KEYWORD initial_entropy
{
"magnetized Bondi" :: "Initial entropy for a radial magnetic field appropriate for Bondi test"
diff --git a/schedule.ccl b/schedule.ccl
index c12093a..76563c9 100644
--- a/schedule.ccl
+++ b/schedule.ccl
@@ -250,7 +250,7 @@ if (CCTK_EQUALS(initial_hydro, "magnetized_bondi_solution"))
} "setup GRHydro vars for the magnetized Bondi solution"
}
-if (CCTK_EQUALS(initial_Bvec, "poloidalmagfield"))
+if (CCTK_EQUALS(initial_Bvec, "poloidalmagfield") || CCTK_EQUALS(initial_Avec, "poloidalmagfield"))
{
# SCHEDULE GRHydro_PoloidalMagFieldM AT CCTK_INITIAL AFTER IN HydroBase_Initial AFTER rnsid_init AFTER TOV_Initial_Data after CCCC_StarMapper_InitialData
SCHEDULE GRHydro_PoloidalMagFieldM AT CCTK_POSTINITIAL
diff --git a/src/GRHydro_PoloidalMagFieldM.F90 b/src/GRHydro_PoloidalMagFieldM.F90
index 52e9ca8..37b4377 100644
--- a/src/GRHydro_PoloidalMagFieldM.F90
+++ b/src/GRHydro_PoloidalMagFieldM.F90
@@ -29,6 +29,9 @@ Keisuke Taniguchi - Phys. Rev. D 78, 024012 (2008)
#define Bconsx(i,j,k) Bcons(i,j,k,1)
#define Bconsy(i,j,k) Bcons(i,j,k,2)
#define Bconsz(i,j,k) Bcons(i,j,k,3)
+#define Avecx(i,j,k) Avec(i,j,k,1)
+#define Avecy(i,j,k) Avec(i,j,k,2)
+#define Avecz(i,j,k) Avec(i,j,k,3)
/*@@
@routine GRHydro_PoloidalMagFieldM
@@ -56,7 +59,7 @@ subroutine GRHydro_PoloidalMagFieldM(CCTK_ARGUMENTS)
DECLARE_CCTK_PARAMETERS
DECLARE_CCTK_FUNCTIONS
- CCTK_INT :: i, j, k, nx, ny, nz
+ CCTK_INT :: i, j, k, nx, ny, nz, set_Avec
CCTK_REAL :: det
CCTK_REAL :: sdet
CCTK_REAL :: dx,dy,dz
@@ -75,8 +78,15 @@ subroutine GRHydro_PoloidalMagFieldM(CCTK_ARGUMENTS)
dy = CCTK_DELTA_SPACE(2)
dz = CCTK_DELTA_SPACE(3)
+ set_Avec = 0
+ if ( CCTK_EQUALS(Bvec_evolution_method,"GRHydro_Avec") ) then
+ set_Avec = 1
+ end if
+
write(*,*)'GRHydro_InitData: Setting up initial poloidal magnetic field'
+
+
do i=2,nx-1
do j=2,ny-1
do k=2,nz-1
@@ -123,6 +133,12 @@ subroutine GRHydro_PoloidalMagFieldM(CCTK_ARGUMENTS)
Bvecy(i,j,k) = - Ax_dz/sdet
Bvecz(i,j,k) = (Ax_dy-Ay_dx)/sdet
+ if ( set_Avec.gt.0 ) then
+ Avecx(i,j,k) = Ax
+ Avecy(i,j,k) = Ay
+ Avecz(i,j,k) = Az
+ end if
+
!Bvecx(i,j,k) = 0.0d0
!Bvecy(i,j,k) = 0.0d0
!Bvecz(i,j,k) = 0.00000001/sdet