aboutsummaryrefslogtreecommitdiff
path: root/src
diff options
context:
space:
mode:
authorpollney <pollney@5c0f84ea-6048-4d6e-bfa4-55cd5f2e0dd7>2002-05-22 20:04:04 +0000
committerpollney <pollney@5c0f84ea-6048-4d6e-bfa4-55cd5f2e0dd7>2002-05-22 20:04:04 +0000
commitc0e4e7577cbd914e04a8770e9112253b905adf82 (patch)
tree861294dca7c279fc916a8d6be98c20d99658bf5c /src
parentc5430d6f7d52fc5c702b019367b21c0a6dcbe472 (diff)
Adding Eric's standing plane wave solution.
git-svn-id: http://svn.einsteintoolkit.org/cactus/EinsteinInitialData/IDLinearWaves/trunk@76 5c0f84ea-6048-4d6e-bfa4-55cd5f2e0dd7
Diffstat (limited to 'src')
-rw-r--r--src/make.code.defn1
-rw-r--r--src/standing_planewaves.F7797
2 files changed, 98 insertions, 0 deletions
diff --git a/src/make.code.defn b/src/make.code.defn
index 9969001..ff57d67 100644
--- a/src/make.code.defn
+++ b/src/make.code.defn
@@ -4,6 +4,7 @@
# Source files in this directory
SRCS = ParamCheck.c\
planewaves.F77\
+ standing_planewaves.F77\
teukwaves.F77
# Subdirectories containing source files
diff --git a/src/standing_planewaves.F77 b/src/standing_planewaves.F77
new file mode 100644
index 0000000..c4af49b
--- /dev/null
+++ b/src/standing_planewaves.F77
@@ -0,0 +1,97 @@
+c Written 2002-05-22 by Erik Schnetter <schnetter@uni-tuebingen.de>
+c $Header$
+
+#include "cctk.h"
+#include "cctk_Arguments.h"
+#include "cctk_Parameters.h"
+
+ subroutine IDLinearWaves_Standing_PlaneWaves (CCTK_ARGUMENTS)
+
+ implicit none
+
+ DECLARE_CCTK_ARGUMENTS
+ DECLARE_CCTK_FUNCTIONS
+ DECLARE_CCTK_PARAMETERS
+
+ CCTK_REAL pi
+ parameter (pi = 3.14159265358979d0)
+ CCTK_REAL aa, kk
+ CCTK_REAL zz, tt
+ CCTK_REAL bb, bbdot
+ integer i,j,k
+
+ aa = amplitude
+ kk = 2*pi / wavelength
+
+ tt = cctk_time
+
+ do k=1,cctk_lsh(3)
+ do j=1,cctk_lsh(2)
+ do i=1,cctk_lsh(1)
+
+ zz = z(i,j,k)
+
+ bb = aa * sin(kk*zz) * cos(kk*tt)
+ bbdot = -kk*aa * sin(kk*zz) * sin(kk*tt)
+
+ gxx(i,j,k) = 1 + bb
+ gxy(i,j,k) = 0
+ gxz(i,j,k) = 0
+ gyy(i,j,k) = 1 - bb
+ gyz(i,j,k) = 0
+ gzz(i,j,k) = 1
+
+ kxx(i,j,k) = bbdot / (-2)
+ kxy(i,j,k) = 0
+ kxz(i,j,k) = 0
+ kyy(i,j,k) = - bbdot / (-2)
+ kyz(i,j,k) = 0
+ kzz(i,j,k) = 0
+
+ end do
+ end do
+ end do
+
+c initialise the conformal factor
+
+ if (CCTK_EQUALS(metric_type, "static conformal")) then
+
+ conformal_state = 1
+
+ if (CCTK_EQUALS(conformal_storage,"factor+derivs")) then
+
+ conformal_state = 2
+
+ else if (CCTK_EQUALS(conformal_storage,"factor+derivs+2nd derivs")) then
+ conformal_state = 3
+
+ end if
+
+ do k=1,cctk_lsh(3)
+ do j=1,cctk_lsh(2)
+ do i=1,cctk_lsh(1)
+
+ psi(i,j,k) = 1
+
+ if (conformal_state .gt. 1) then
+ psix(i,j,k) = 0
+ psiy(i,j,k) = 0
+ psiz(i,j,k) = 0
+ end if
+
+ if (conformal_state .gt. 2) then
+ psixy(i,j,k) = 0
+ psixz(i,j,k) = 0
+ psiyz(i,j,k) = 0
+ psixx(i,j,k) = 0
+ psiyy(i,j,k) = 0
+ psizz(i,j,k) = 0
+ end if
+
+ end do
+ end do
+ end do
+
+ end if
+
+ end