diff options
author | miguel <miguel@a678b1cf-93e1-4b43-a69d-d43939e66649> | 2001-04-12 10:51:59 +0000 |
---|---|---|
committer | miguel <miguel@a678b1cf-93e1-4b43-a69d-d43939e66649> | 2001-04-12 10:51:59 +0000 |
commit | 6fd09e4d86b3a13278c83874f976c3ba2aa7cd42 (patch) | |
tree | e8575326bc17006dcdcfccb0c4b3851131a34b04 | |
parent | 65c6f4ada5ed70fb476c58de5866b98eb7c911e8 (diff) |
Cleaning up.
git-svn-id: http://svn.einsteintoolkit.org/cactus/EinsteinInitialData/IDBrillData/trunk@36 a678b1cf-93e1-4b43-a69d-d43939e66649
-rw-r--r-- | interface.ccl | 7 | ||||
-rw-r--r-- | par/brilldata.par | 50 | ||||
-rw-r--r-- | schedule.ccl | 5 | ||||
-rw-r--r-- | src/brilldata.F | 25 | ||||
-rw-r--r-- | src/finishbrilldata.F | 15 | ||||
-rw-r--r-- | src/setupbrilldata2D.F | 5 |
6 files changed, 53 insertions, 54 deletions
diff --git a/interface.ccl b/interface.ccl index 887b578..2c45e45 100644 --- a/interface.ccl +++ b/interface.ccl @@ -1,3 +1,6 @@ +# Interface definition for thorn BrillData +# $Header$ + implements: IDBrillData inherits: einstein ellbase @@ -5,13 +8,13 @@ inherits: einstein ellbase private: -real linearelliptic TYPE=GF +CCTK_REAL linearelliptic TYPE=GF { Mlinear, Nsource } "Coefficients for linear elliptic equation" -real brillconf TYPE=GF +CCTK_REAL brillconf TYPE=GF { brillpsi } "Conformal factor for Brill data" diff --git a/par/brilldata.par b/par/brilldata.par index 47f15af..5ce8e7f 100644 --- a/par/brilldata.par +++ b/par/brilldata.par @@ -1,50 +1,42 @@ -################################################################# +################################################ # # brilldata.par # -# Simply Brill wave data of the Hoz et.al. type +# Simple Brill wave data of the Holz et.al. type # -################################################################# +################################################ -ActiveThorns = "iobasic ellbase ellsor IDBrillData pugh pughslab pughreduce einstein cartgrid3d ioutil ioascii" +ActiveThorns = "iobasic ellbase ellsor idbrilldata pugh pughreduce einstein cartgrid3d ioutil ioascii" -# Brill wave initial data +# General -Einstein::initial_data = "brilldata" +driver::global_nx = 32 +driver::global_ny = 32 +driver::global_nz = 32 -# Brill wave parameters [Holz et.al. type data] +grid::type = "byspacing" +grid::dxyz = 0.2 +grid::domain = "octant" -idbrilldata::brill_q = 2 # Form of function q [0,1,2] -idbrilldata::brill_a = 1.0 # Amplitude (default 0) +cactus::cctk_itlast = 0 -idbrilldata::outputpsi = "yes" - -# Elliptic solver. +# Brill wave initial data -#idbrilldata::solver = "sor" #[petsc,sor,bam] +Einstein::initial_data = "brilldata" -# General +# Brill wave parameters [Holz et.al. type data] -driver::global_nx = 32 -driver::global_ny = 32 -driver::global_nz = 32 +idbrilldata::brill_q = 2 # Form of function q [0,1,2] +idbrilldata::brill_a = 1.0 # Amplitude (default 0) -grid::type = "byspacing" -grid::dxyz = 0.2 -grid::domain = "octant" +# Elliptic solver. -cactus::cctk_itlast = 0.0 +idbrilldata::brill_solver = "sor" #[petsc,sor,bam] # Output. -IO::outdir = "BrillData" - IO::out_every = 1 -IOASCII::out1D_vars = "einstein::gxx einstein::gyy einstein::gzz - einstein::gxy einstein::gxz einstein::gyz - IDBrillData::brillpsi" +IOASCII::out1D_vars = "einstein::gxx einstein::gyy einstein::gzz einstein::gxy einstein::gxz einstein::gyz idbrilldata::brillpsi" -IOBasic::outScalar_vars = "einstein::gxx einstein::gyy einstein::gzz - einstein::gxy einstein::gxz einstein::gyz - IDBrillData::brillpsi" +IOBasic::outScalar_vars = "einstein::gxx einstein::gyy einstein::gzz einstein::gxy einstein::gxz einstein::gyz idbrilldata::brillpsi" diff --git a/schedule.ccl b/schedule.ccl index b6f879f..34c3f74 100644 --- a/schedule.ccl +++ b/schedule.ccl @@ -3,6 +3,10 @@ if (CCTK_Equals(initial_data,"brilldata")) { + + STORAGE: brillconf + STORAGE: linearelliptic + schedule BrilData_InitSymBound at CCTK_BASEGRID { LANG: C @@ -10,7 +14,6 @@ if (CCTK_Equals(initial_data,"brilldata")) schedule brilldata at CCTK_INITIAL { - STORAGE: linearelliptic,brillconf LANG: Fortran } "Construct Brill wave initial data" } diff --git a/src/brilldata.F b/src/brilldata.F index 24a745e..253ab15 100644 --- a/src/brilldata.F +++ b/src/brilldata.F @@ -1,12 +1,11 @@ + #include "cctk.h" #include "cctk_Parameters.h" #include "cctk_Arguments.h" - subroutine brilldata(CCTK_FARGUMENTS) +#include "CactusEinstein/Einstein/src/Einstein.h" -c Author: Carsten Gundlach. -c -c Driver routine for calculating Brill wave initial data. + subroutine brilldata(CCTK_FARGUMENTS) implicit none @@ -14,12 +13,11 @@ c Driver routine for calculating Brill wave initial data. DECLARE_CCTK_PARAMETERS DECLARE_CCTK_FUNCTIONS -c Declare arrays for the linear elliptic solver call: - - integer Mlin_index,Nsrc_index,field_index,ierr - integer metpsi_index(7) + integer Mlin_index,Nsrc_index,field_index,ierr + integer metpsi_index(7) + integer ierr - CCTK_REAL AbsTol(3), RelTol(3) + CCTK_REAL AbsTol(3),RelTol(3) c Set up background metric and coefficients for linear solve. @@ -41,8 +39,8 @@ c conformal factor. First some preliminaries. call CCTK_VarIndex(metpsi_index(7),"einstein::psi") call CCTK_VarIndex(field_index,"IDBrillData::brillpsi") - call CCTK_VarIndex(Mlin_index,"IDBrillData::Mlinear") - call CCTK_VarIndex(Nsrc_index,"IDBrillData::Nsource") + call CCTK_VarIndex(Mlin_index, "IDBrillData::Mlinear") + call CCTK_VarIndex(Nsrc_index, "IDBrillData::Nsource") AbsTol(1)= brill_thresh AbsTol(2)= -1 @@ -83,9 +81,10 @@ c Elliptic solver. . field_index,Mlin_index,Nsrc_index,AbsTol,RelTol,"bam") end if -c Synchronize conformal factor. +c Synchronization and symmetry boundaries. - call CCTK_SyncGroup(cctkGH,"einstein::confac") + call CCTK_SyncGroup(cctkGH,"IDBrillData::brillconf") + call CartSymGN(ierr,cctkGH,"IDBrillData::brillconf") c Reconstruct physical metric. diff --git a/src/finishbrilldata.F b/src/finishbrilldata.F index 4cd4f18..ee49b57 100644 --- a/src/finishbrilldata.F +++ b/src/finishbrilldata.F @@ -1,14 +1,11 @@ + #include "cctk.h" #include "cctk_Parameters.h" #include "cctk_Arguments.h" -#include "../../../CactusEinstein/Einstein/src/Einstein.h" - subroutine finishbrilldata(CCTK_FARGUMENTS) +#include "CactusEinstein/Einstein/src/Einstein.h" -C Author: Carsten Gundlach. -C -C Set up 3-metric, extrinsic curvature and BM variables -C once the conformal factor has been found. + subroutine finishbrilldata(CCTK_FARGUMENTS) implicit none @@ -66,7 +63,7 @@ C C e^2q (drho^2 + dz^2) + rho^2 dphi^2 = C e^2q (dx^2 + dy^2 + dz^2) + (1-e^2q) (xdy-ydx)^2/rho^2 C -C The individual coefficients can be read off as +C The individual coefficients can be read off as if (rho1 .gt. rhofudge) then @@ -77,8 +74,8 @@ C The individual coefficients can be read off as else -C This fudge assumes that q = O(rho^2) near the axis. Which -C it should be, or the data will be singular. +C This fudge assumes that q = O(rho^2) near the axis. Which +C it should be, or the data will be singular. gzz(i,j,k) = psi4 gxx(i,j,k) = psi4 diff --git a/src/setupbrilldata2D.F b/src/setupbrilldata2D.F index 9b00de9..5643cf0 100644 --- a/src/setupbrilldata2D.F +++ b/src/setupbrilldata2D.F @@ -58,9 +58,14 @@ C [Delta(g) + Mlinear] psi = Nsource. psi = one + psix = zero + psiy = zero + psiz = zero + gxx = one gyy = one gzz = one + gxy = zero gxz = zero gyz = zero |