aboutsummaryrefslogtreecommitdiff
path: root/Examples
diff options
context:
space:
mode:
authorIan Hinder <ian.hinder@aei.mpg.de>2012-02-14 11:33:36 +0100
committerIan Hinder <ian.hinder@aei.mpg.de>2012-02-14 11:33:36 +0100
commitf412f48340dfb00208f3fe87de7b2e8863d9c3cc (patch)
tree20da3540c14fd8828a18680e10ef02d4d3da4b2d /Examples
parent0d99b18cc556b7c4478b2676e312d1c484e1e244 (diff)
Regenerate SimpleWaveCaKernel
Diffstat (limited to 'Examples')
-rw-r--r--Examples/SimpleWaveCaKernel/cakernel.ccl46
-rw-r--r--Examples/SimpleWaveCaKernel/interface.ccl14
-rw-r--r--Examples/SimpleWaveCaKernel/param.ccl24
-rw-r--r--Examples/SimpleWaveCaKernel/schedule.ccl59
-rw-r--r--Examples/SimpleWaveCaKernel/src/CaKernel__calc_bound_rhs.code (renamed from Examples/SimpleWaveCaKernel/src/CaKernel__initial_gaussian.code)46
-rw-r--r--Examples/SimpleWaveCaKernel/src/CaKernel__copy_to_device.code (renamed from Examples/SimpleWaveCaKernel/src/CaKernel__rk1.code)14
-rw-r--r--Examples/SimpleWaveCaKernel/src/RegisterSymmetries.cc5
-rw-r--r--Examples/SimpleWaveCaKernel/src/initial_gaussian.cc131
-rw-r--r--Examples/SimpleWaveCaKernel/src/make.code.defn2
9 files changed, 253 insertions, 88 deletions
diff --git a/Examples/SimpleWaveCaKernel/cakernel.ccl b/Examples/SimpleWaveCaKernel/cakernel.ccl
index ff8d47f..75c394c 100644
--- a/Examples/SimpleWaveCaKernel/cakernel.ccl
+++ b/Examples/SimpleWaveCaKernel/cakernel.ccl
@@ -1,25 +1,4 @@
-CCTK_CUDA_KERNEL initial_gaussian TYPE=gpu_cuda/3dblock STENCIL="0,0,0,0,0,0" TILE="8,8,8" SHARECODE=yes
-{
- CCTK_CUDA_KERNEL_VARIABLE cached=no intent=out
- {
- phi
- }
- "phi"
-
- CCTK_CUDA_KERNEL_VARIABLE cached=no intent=out
- {
- pi
- }
- "pi"
-
- CCTK_CUDA_KERNEL_VARIABLE cached=no intent=in
- {
- x
- }
- "x"
-}
-
-CCTK_CUDA_KERNEL calc_rhs TYPE=gpu_cuda/3dblock STENCIL="1,1,1,1,1,1" TILE="8,8,8" SHARECODE=yes
+CCTK_CUDA_KERNEL calc_rhs TYPE=gpu_cuda/3dblock TILE="8,8,8" SHARECODE=yes STENCIL="1,1,1,1,1,1"
{
CCTK_CUDA_KERNEL_VARIABLE cached=no intent=in
{
@@ -46,7 +25,7 @@ CCTK_CUDA_KERNEL calc_rhs TYPE=gpu_cuda/3dblock STENCIL="1,1,1,1,1,1" TILE="8,8,
"pirhs"
}
-CCTK_CUDA_KERNEL rk1 TYPE=gpu_cuda/3dblock STENCIL="0,0,0,0,0,0" TILE="8,8,8" SHARECODE=yes
+CCTK_CUDA_KERNEL copy_to_device TYPE=gpu_cuda/3dblock TILE="8,8,8" SHARECODE=yes STENCIL="0,0,0,0,0,0"
{
CCTK_CUDA_KERNEL_VARIABLE cached=no intent=inout
{
@@ -54,22 +33,31 @@ CCTK_CUDA_KERNEL rk1 TYPE=gpu_cuda/3dblock STENCIL="0,0,0,0,0,0" TILE="8,8,8" SH
}
"phi"
- CCTK_CUDA_KERNEL_VARIABLE cached=no intent=in
+ CCTK_CUDA_KERNEL_VARIABLE cached=no intent=inout
+ {
+ pi
+ }
+ "pi"
+}
+
+CCTK_CUDA_KERNEL calc_bound_rhs TYPE=gpu_cuda/boundary TILE="8,8,8" SHARECODE=yes
+{
+ CCTK_CUDA_KERNEL_VARIABLE cached=no intent=out
{
phirhs
}
"phirhs"
- CCTK_CUDA_KERNEL_VARIABLE cached=no intent=inout
+ CCTK_CUDA_KERNEL_VARIABLE cached=no intent=out
{
- pi
+ pirhs
}
- "pi"
+ "pirhs"
CCTK_CUDA_KERNEL_VARIABLE cached=no intent=in
{
- pirhs
+ xCopy
}
- "pirhs"
+ "xCopy"
}
diff --git a/Examples/SimpleWaveCaKernel/interface.ccl b/Examples/SimpleWaveCaKernel/interface.ccl
index 0435e8b..710d050 100644
--- a/Examples/SimpleWaveCaKernel/interface.ccl
+++ b/Examples/SimpleWaveCaKernel/interface.ccl
@@ -28,25 +28,31 @@ CCTK_INT FUNCTION Boundary_SelectVarForBC(CCTK_POINTER_TO_CONST IN GH, CCTK_INT
USES FUNCTION Boundary_SelectVarForBC
public:
-CCTK_REAL phi_g type=GF timelevels=1 tags=''
+CCTK_REAL xCopy_g type=GF timelevels=1 tags=''
+{
+ xCopy
+} "xCopy_g"
+
+public:
+CCTK_REAL phi_g type=GF timelevels=2 tags=''
{
phi
} "phi_g"
public:
-CCTK_REAL pi_g type=GF timelevels=1 tags=''
+CCTK_REAL pi_g type=GF timelevels=2 tags=''
{
pi
} "pi_g"
public:
-CCTK_REAL phi_grhs type=GF timelevels=1 tags=''
+CCTK_REAL phi_grhs type=GF timelevels=2 tags=''
{
phirhs
} "phi_grhs"
public:
-CCTK_REAL pi_grhs type=GF timelevels=1 tags=''
+CCTK_REAL pi_grhs type=GF timelevels=2 tags=''
{
pirhs
} "pi_grhs"
diff --git a/Examples/SimpleWaveCaKernel/param.ccl b/Examples/SimpleWaveCaKernel/param.ccl
index ecb8b5e..845e24e 100644
--- a/Examples/SimpleWaveCaKernel/param.ccl
+++ b/Examples/SimpleWaveCaKernel/param.ccl
@@ -31,19 +31,19 @@ CCTK_INT SimpleWaveCaKernel_MaxNumArrayEvolvedVars "Number of Array evolved vari
restricted:
CCTK_INT timelevels "Number of active timelevels" STEERABLE=RECOVER
{
- 0:1 :: ""
-} 1
+ 0:2 :: ""
+} 2
restricted:
CCTK_INT rhs_timelevels "Number of active RHS timelevels" STEERABLE=RECOVER
{
- 0:1 :: ""
+ 0:2 :: ""
} 1
restricted:
CCTK_INT other_timelevels "Number of active timelevels for non-evolved grid functions" STEERABLE=RECOVER
{
- 0:1 :: ""
+ 0:2 :: ""
} 1
restricted:
@@ -59,7 +59,13 @@ CCTK_INT calc_rhs_calc_every "calc_rhs_calc_every" STEERABLE=ALWAYS
} 1
restricted:
-CCTK_INT rk1_calc_every "rk1_calc_every" STEERABLE=ALWAYS
+CCTK_INT copy_to_device_calc_every "copy_to_device_calc_every" STEERABLE=ALWAYS
+{
+ *:* :: ""
+} 1
+
+restricted:
+CCTK_INT calc_bound_rhs_calc_every "calc_bound_rhs_calc_every" STEERABLE=ALWAYS
{
*:* :: ""
} 1
@@ -77,7 +83,13 @@ CCTK_INT calc_rhs_calc_offset "calc_rhs_calc_offset" STEERABLE=ALWAYS
} 0
restricted:
-CCTK_INT rk1_calc_offset "rk1_calc_offset" STEERABLE=ALWAYS
+CCTK_INT copy_to_device_calc_offset "copy_to_device_calc_offset" STEERABLE=ALWAYS
+{
+ *:* :: ""
+} 0
+
+restricted:
+CCTK_INT calc_bound_rhs_calc_offset "calc_bound_rhs_calc_offset" STEERABLE=ALWAYS
{
*:* :: ""
} 0
diff --git a/Examples/SimpleWaveCaKernel/schedule.ccl b/Examples/SimpleWaveCaKernel/schedule.ccl
index 83851b5..5c777d8 100644
--- a/Examples/SimpleWaveCaKernel/schedule.ccl
+++ b/Examples/SimpleWaveCaKernel/schedule.ccl
@@ -1,25 +1,46 @@
# File produced by Kranc
+if (other_timelevels == 1)
+{
+ STORAGE: xCopy_g[1]
+}
+
if (timelevels == 1)
{
STORAGE: phi_g[1]
}
+if (timelevels == 2)
+{
+ STORAGE: phi_g[2]
+}
if (timelevels == 1)
{
STORAGE: pi_g[1]
}
+if (timelevels == 2)
+{
+ STORAGE: pi_g[2]
+}
if (rhs_timelevels == 1)
{
STORAGE: phi_grhs[1]
}
+if (rhs_timelevels == 2)
+{
+ STORAGE: phi_grhs[2]
+}
if (rhs_timelevels == 1)
{
STORAGE: pi_grhs[1]
}
+if (rhs_timelevels == 2)
+{
+ STORAGE: pi_grhs[2]
+}
schedule SimpleWaveCaKernel_Startup at STARTUP
{
@@ -27,49 +48,49 @@ schedule SimpleWaveCaKernel_Startup at STARTUP
OPTIONS: meta
} "create banner"
-schedule SimpleWaveCaKernel_RegisterVars in MoL_Register
-{
- LANG: C
- OPTIONS: meta
-} "Register Variables for MoL"
-
schedule SimpleWaveCaKernel_RegisterSymmetries in SymmetryRegister
{
LANG: C
OPTIONS: meta
} "register symmetries"
-schedule CAKERNEL_Launch_initial_gaussian AT INITIAL
+schedule initial_gaussian AT INITIAL
{
LANG: C
READS: grid::coordinates
WRITES: SimpleWaveCaKernel::phi_g
WRITES: SimpleWaveCaKernel::pi_g
+ WRITES: SimpleWaveCaKernel::xCopy_g
} "initial_gaussian"
-schedule CAKERNEL_Launch_calc_rhs at EVOL
+schedule CAKERNEL_Launch_calc_rhs in MoL_CalcRHS
{
LANG: C
- SYNC: phi_grhs
- SYNC: pi_grhs
+ TAGS: Device=1
READS: SimpleWaveCaKernel::phi_g
READS: SimpleWaveCaKernel::pi_g
WRITES: SimpleWaveCaKernel::phi_grhs
WRITES: SimpleWaveCaKernel::pi_grhs
} "calc_rhs"
-schedule CAKERNEL_Launch_rk1 at EVOL after calc_rhs
+schedule CAKERNEL_Launch_copy_to_device at INITIAL after initial_gaussian
{
LANG: C
- SYNC: phi_g
- SYNC: pi_g
+ TAGS: Device=1
READS: SimpleWaveCaKernel::phi_g
- READS: SimpleWaveCaKernel::phi_grhs
READS: SimpleWaveCaKernel::pi_g
- READS: SimpleWaveCaKernel::pi_grhs
WRITES: SimpleWaveCaKernel::phi_g
WRITES: SimpleWaveCaKernel::pi_g
-} "rk1"
+} "copy_to_device"
+
+schedule CAKERNEL_Launch_calc_bound_rhs in MoL_RHSBoundaries
+{
+ LANG: C
+ TAGS: Device=1
+ READS: SimpleWaveCaKernel::xCopy_g
+ WRITES: SimpleWaveCaKernel::phi_grhs
+ WRITES: SimpleWaveCaKernel::pi_grhs
+} "calc_bound_rhs"
schedule SimpleWaveCaKernel_SelectBoundConds in MoL_PostStep
{
@@ -85,6 +106,12 @@ schedule SimpleWaveCaKernel_CheckBoundaries at BASEGRID
OPTIONS: meta
} "check boundaries treatment"
+schedule SimpleWaveCaKernel_RegisterVars in MoL_Register
+{
+ LANG: C
+ OPTIONS: meta
+} "Register Variables for MoL"
+
schedule group ApplyBCs as SimpleWaveCaKernel_ApplyBCs in MoL_PostStep after SimpleWaveCaKernel_SelectBoundConds
{
# no language specified
diff --git a/Examples/SimpleWaveCaKernel/src/CaKernel__initial_gaussian.code b/Examples/SimpleWaveCaKernel/src/CaKernel__calc_bound_rhs.code
index 40f6037..863ad4f 100644
--- a/Examples/SimpleWaveCaKernel/src/CaKernel__initial_gaussian.code
+++ b/Examples/SimpleWaveCaKernel/src/CaKernel__calc_bound_rhs.code
@@ -14,7 +14,7 @@
#define SQR(x) ((x) * (x))
#define CUB(x) ((x) * (x) * (x))
-CAKERNEL_initial_gaussian_Begin
+CAKERNEL_calc_bound_rhs_Begin
/* Include user-supplied include files */
@@ -51,26 +51,24 @@ CAKERNEL_initial_gaussian_Begin
/* Calculate temporaries and arrays functions */
/* Copy local copies back to grid functions */
- CAKERNEL_initial_gaussian_Computations_Begin
-
- /* Assign local copies of grid functions */
-
- CCTK_REAL xL = I3D(x,0,0,0);
-
-
- /* Include user supplied include files */
-
- /* Precompute derivatives */
-
- /* Calculate temporaries and grid functions */
- CCTK_REAL phiL = exp(-100.*SQR(xL + t));
-
- CCTK_REAL piL = -200.*(xL + t)*exp(-100.*SQR(xL + t));
-
- /* Copy local copies back to grid functions */
- I3D(phi,0,0,0) = phiL;
- I3D(pi,0,0,0) = piL;
-
- CAKERNEL_initial_gaussian_Computations_End
-
-CAKERNEL_initial_gaussian_End
+
+ /* Assign local copies of grid functions */
+
+ CCTK_REAL xCopyL = I3D(xCopy,0,0,0);
+
+
+ /* Include user supplied include files */
+
+ /* Precompute derivatives */
+
+ /* Calculate temporaries and grid functions */
+ CCTK_REAL phirhsL = -200.*(xCopyL + t)*exp(-100.*SQR(xCopyL + t));
+
+ CCTK_REAL pirhsL = exp(-100.*SQR(xCopyL + t))*(-200. +
+ 80000.*xCopyL*t + 40000.*(SQR(xCopyL) + SQR(t)));
+
+ /* Copy local copies back to grid functions */
+ I3D(phirhs,0,0,0) = phirhsL;
+ I3D(pirhs,0,0,0) = pirhsL;
+
+CAKERNEL_calc_bound_rhs_End
diff --git a/Examples/SimpleWaveCaKernel/src/CaKernel__rk1.code b/Examples/SimpleWaveCaKernel/src/CaKernel__copy_to_device.code
index e3feb54..24b4566 100644
--- a/Examples/SimpleWaveCaKernel/src/CaKernel__rk1.code
+++ b/Examples/SimpleWaveCaKernel/src/CaKernel__copy_to_device.code
@@ -14,7 +14,7 @@
#define SQR(x) ((x) * (x))
#define CUB(x) ((x) * (x) * (x))
-CAKERNEL_rk1_Begin
+CAKERNEL_copy_to_device_Begin
/* Include user-supplied include files */
@@ -51,14 +51,12 @@ CAKERNEL_rk1_Begin
/* Calculate temporaries and arrays functions */
/* Copy local copies back to grid functions */
- CAKERNEL_rk1_Computations_Begin
+ CAKERNEL_copy_to_device_Computations_Begin
/* Assign local copies of grid functions */
CCTK_REAL phiL = I3D(phi,0,0,0);
- CCTK_REAL phirhsL = I3D(phirhs,0,0,0);
CCTK_REAL piL = I3D(pi,0,0,0);
- CCTK_REAL pirhsL = I3D(pirhs,0,0,0);
/* Include user supplied include files */
@@ -66,14 +64,14 @@ CAKERNEL_rk1_Begin
/* Precompute derivatives */
/* Calculate temporaries and grid functions */
- phiL = phiL + phirhsL*dt;
+ phiL = phiL;
- piL = piL + pirhsL*dt;
+ piL = piL;
/* Copy local copies back to grid functions */
I3D(phi,0,0,0) = phiL;
I3D(pi,0,0,0) = piL;
- CAKERNEL_rk1_Computations_End
+ CAKERNEL_copy_to_device_Computations_End
-CAKERNEL_rk1_End
+CAKERNEL_copy_to_device_End
diff --git a/Examples/SimpleWaveCaKernel/src/RegisterSymmetries.cc b/Examples/SimpleWaveCaKernel/src/RegisterSymmetries.cc
index 128eb50..c099f5b 100644
--- a/Examples/SimpleWaveCaKernel/src/RegisterSymmetries.cc
+++ b/Examples/SimpleWaveCaKernel/src/RegisterSymmetries.cc
@@ -26,4 +26,9 @@ extern "C" void SimpleWaveCaKernel_RegisterSymmetries(CCTK_ARGUMENTS)
sym[2] = 1;
SetCartSymVN(cctkGH, sym, "SimpleWaveCaKernel::pi");
+ sym[0] = 1;
+ sym[1] = 1;
+ sym[2] = 1;
+ SetCartSymVN(cctkGH, sym, "SimpleWaveCaKernel::xCopy");
+
}
diff --git a/Examples/SimpleWaveCaKernel/src/initial_gaussian.cc b/Examples/SimpleWaveCaKernel/src/initial_gaussian.cc
new file mode 100644
index 0000000..55bd228
--- /dev/null
+++ b/Examples/SimpleWaveCaKernel/src/initial_gaussian.cc
@@ -0,0 +1,131 @@
+/* File produced by Kranc */
+
+#define KRANC_C
+
+#include <assert.h>
+#include <math.h>
+#include <stdio.h>
+#include <stdlib.h>
+#include <string.h>
+#include "cctk.h"
+#include "cctk_Arguments.h"
+#include "cctk_Parameters.h"
+#include "GenericFD.h"
+#include "Differencing.h"
+#include "cctk_Loop.h"
+#include "loopcontrol.h"
+
+/* Define macros used in calculations */
+#define INITVALUE (42)
+#define QAD(x) (SQR(SQR(x)))
+#define INV(x) ((1.0) / (x))
+#define SQR(x) ((x) * (x))
+#define CUB(x) ((x) * (x) * (x))
+
+static void initial_gaussian_Body(cGH const * restrict const cctkGH, int const dir, int const face, CCTK_REAL const normal[3], CCTK_REAL const tangentA[3], CCTK_REAL const tangentB[3], int const imin[3], int const imax[3], int const n_subblock_gfs, CCTK_REAL * restrict const subblock_gfs[])
+{
+ DECLARE_CCTK_ARGUMENTS;
+ DECLARE_CCTK_PARAMETERS;
+
+
+ /* Include user-supplied include files */
+
+ /* Initialise finite differencing variables */
+ ptrdiff_t const di = 1;
+ ptrdiff_t const dj = CCTK_GFINDEX3D(cctkGH,0,1,0) - CCTK_GFINDEX3D(cctkGH,0,0,0);
+ ptrdiff_t const dk = CCTK_GFINDEX3D(cctkGH,0,0,1) - CCTK_GFINDEX3D(cctkGH,0,0,0);
+ ptrdiff_t const cdi = sizeof(CCTK_REAL) * di;
+ ptrdiff_t const cdj = sizeof(CCTK_REAL) * dj;
+ ptrdiff_t const cdk = sizeof(CCTK_REAL) * dk;
+ CCTK_REAL const dx = ToReal(CCTK_DELTA_SPACE(0));
+ CCTK_REAL const dy = ToReal(CCTK_DELTA_SPACE(1));
+ CCTK_REAL const dz = ToReal(CCTK_DELTA_SPACE(2));
+ CCTK_REAL const dt = ToReal(CCTK_DELTA_TIME);
+ CCTK_REAL const t = ToReal(cctk_time);
+ CCTK_REAL const dxi = INV(dx);
+ CCTK_REAL const dyi = INV(dy);
+ CCTK_REAL const dzi = INV(dz);
+ CCTK_REAL const khalf = 0.5;
+ CCTK_REAL const kthird = 1/3.0;
+ CCTK_REAL const ktwothird = 2.0/3.0;
+ CCTK_REAL const kfourthird = 4.0/3.0;
+ CCTK_REAL const keightthird = 8.0/3.0;
+ CCTK_REAL const hdxi = 0.5 * dxi;
+ CCTK_REAL const hdyi = 0.5 * dyi;
+ CCTK_REAL const hdzi = 0.5 * dzi;
+
+ /* Initialize predefined quantities */
+ CCTK_REAL const p1o2dx = 0.5*INV(dx);
+ CCTK_REAL const p1o2dy = 0.5*INV(dy);
+ CCTK_REAL const p1o2dz = 0.5*INV(dz);
+ CCTK_REAL const p1odx2 = INV(SQR(dx));
+ CCTK_REAL const p1ody2 = INV(SQR(dy));
+ CCTK_REAL const p1odz2 = INV(SQR(dz));
+
+ /* Assign local copies of arrays functions */
+
+
+
+ /* Calculate temporaries and arrays functions */
+
+ /* Copy local copies back to grid functions */
+
+ /* Loop over the grid points */
+ #pragma omp parallel
+ CCTK_LOOP3(initial_gaussian,
+ i,j,k, imin[0],imin[1],imin[2], imax[0],imax[1],imax[2],
+ cctk_lsh[0],cctk_lsh[1],cctk_lsh[2])
+ {
+ ptrdiff_t const index = di*i + dj*j + dk*k;
+
+ /* Assign local copies of grid functions */
+
+ CCTK_REAL xL = x[index];
+
+
+ /* Include user supplied include files */
+
+ /* Precompute derivatives */
+
+ /* Calculate temporaries and grid functions */
+ CCTK_REAL phiL = exp(-100.*SQR(xL + t));
+
+ CCTK_REAL piL = -200.*(xL + t)*exp(-100.*SQR(xL + t));
+
+ CCTK_REAL xCopyL = xL;
+
+ /* Copy local copies back to grid functions */
+ phi[index] = phiL;
+ pi[index] = piL;
+ xCopy[index] = xCopyL;
+ }
+ CCTK_ENDLOOP3(initial_gaussian);
+}
+
+extern "C" void initial_gaussian(CCTK_ARGUMENTS)
+{
+ DECLARE_CCTK_ARGUMENTS;
+ DECLARE_CCTK_PARAMETERS;
+
+
+ if (verbose > 1)
+ {
+ CCTK_VInfo(CCTK_THORNSTRING,"Entering initial_gaussian_Body");
+ }
+
+ if (cctk_iteration % initial_gaussian_calc_every != initial_gaussian_calc_offset)
+ {
+ return;
+ }
+
+ const char *groups[] = {"grid::coordinates","SimpleWaveCaKernel::phi_g","SimpleWaveCaKernel::pi_g","SimpleWaveCaKernel::xCopy_g"};
+ GenericFD_AssertGroupStorage(cctkGH, "initial_gaussian", 4, groups);
+
+
+ GenericFD_LoopOverEverything(cctkGH, &initial_gaussian_Body);
+
+ if (verbose > 1)
+ {
+ CCTK_VInfo(CCTK_THORNSTRING,"Leaving initial_gaussian_Body");
+ }
+}
diff --git a/Examples/SimpleWaveCaKernel/src/make.code.defn b/Examples/SimpleWaveCaKernel/src/make.code.defn
index 5e729f5..c43661b 100644
--- a/Examples/SimpleWaveCaKernel/src/make.code.defn
+++ b/Examples/SimpleWaveCaKernel/src/make.code.defn
@@ -1,3 +1,3 @@
# File produced by Kranc
-SRCS = Startup.cc RegisterMoL.cc RegisterSymmetries.cc Boundaries.cc
+SRCS = Startup.cc RegisterSymmetries.cc RegisterMoL.cc initial_gaussian.cc Boundaries.cc