aboutsummaryrefslogtreecommitdiff
path: root/src
diff options
context:
space:
mode:
authorallen <allen@90eea020-d82d-4da5-bf6e-4ee79ff7632f>1999-09-14 11:37:26 +0000
committerallen <allen@90eea020-d82d-4da5-bf6e-4ee79ff7632f>1999-09-14 11:37:26 +0000
commit00f252048d442fc26430f0b6833c63464cefa81c (patch)
tree78c7f6ae3652c24aad8bc7c0742f17487e60d17c /src
parent9b0e42d3d60ca20b45b0394b3c59b8a58541a676 (diff)
Design changes
git-svn-id: http://svn.cactuscode.org/arrangements/CactusWave/WaveToyC/trunk@8 90eea020-d82d-4da5-bf6e-4ee79ff7632f
Diffstat (limited to 'src')
-rw-r--r--src/WaveToy.c127
1 files changed, 63 insertions, 64 deletions
diff --git a/src/WaveToy.c b/src/WaveToy.c
index 7bb2db7..652198e 100644
--- a/src/WaveToy.c
+++ b/src/WaveToy.c
@@ -17,52 +17,6 @@
#include "CactusBase/Boundary/src/Boundary.h"
- /*@@
- @routine WaveToyC_Boundaries
- @date
- @author Tom Goodale
- @desc
- Boundary conditions for the wave equation
- @enddesc
- @calls ApplyFlatBC,ApplyRadiativeBC
- @calledby
- @history
-
- @endhistory
-
-@@*/
-
-void WaveToyC_Boundaries(CCTK_CARGUMENTS)
-{
- DECLARE_CCTK_CARGUMENTS
- DECLARE_CCTK_PARAMETERS
-
- int ierr;
- int sw[3];
-
- /* Set the stencil width */
- sw[0]=1;
- sw[1]=1;
- sw[2]=1;
-
- ApplySymmetry(cctkGH,"wavetoyc::scalartmps");
-
- if (CCTK_EQUALS(bound,"flat"))
- {
- ierr = ApplyFlatBC(cctkGH,sw,"wavetoyc::phi_new");
- }
- else if (CCTK_Equals(bound,"radiation"))
- {
- ierr = ApplyRadiativeBC(cctkGH,0,sw,"wavetoyc::phi_new","wavetoy::phi");
- }
-
- if (ierr < 0)
- {
- CCTK_WARN(0,"Boundary conditions not applied - giving up!");
- }
-
-}
-
/*@@
@routine WaveToyC_Evolution
@@ -71,7 +25,7 @@ void WaveToyC_Boundaries(CCTK_CARGUMENTS)
@desc
Evolution for the wave equation
@enddesc
- @calls CCTK_SyncGroup, wavetoy_boundaries
+ @calls CCTK_SyncGroup, WaveToyC_Boundaries
@calledby
@history
@@ -81,18 +35,19 @@ void WaveToyC_Boundaries(CCTK_CARGUMENTS)
void WaveToyC_Evolution(CCTK_CARGUMENTS)
{
+
DECLARE_CCTK_CARGUMENTS
- DECLARE_CCTK_PARAMETERS
- int i,j,k;
+ int i,j,k;
+ int index;
int istart, jstart, kstart, iend, jend, kend;
CCTK_REAL dx,dy,dz,dt,dx2,dy2,dz2,dt2;
/* Set up shorthands */
- dx = cctk_delta_space[0];
- dy = cctk_delta_space[1];
- dz = cctk_delta_space[2];
- dt = cctk_delta_time;
+ dx = CCTK_DELTA_SPACE(0);
+ dy = CCTK_DELTA_SPACE(1);
+ dz = CCTK_DELTA_SPACE(2);
+ dt = CCTK_DELTA_TIME;
dx2=dx*dx;
dy2=dy*dy;
@@ -115,11 +70,12 @@ void WaveToyC_Evolution(CCTK_CARGUMENTS)
{
for (i=istart; i<iend; i++)
{
+
phi_new[CCTK_GFINDEX3D(cctkGH,i,j,k)] =
2*(1 - (dt2)*(1/dx2 + 1/dy2 + 1/dz2))*
- phi[CCTK_GFINDEX3D(cctkGH,i,j,k)] -
- phi_old[CCTK_GFINDEX3D(cctkGH,i,j,k)]
- + (dt2) *
+ phi[CCTK_GFINDEX3D(cctkGH,i,j,k)] -
+ phi_old[CCTK_GFINDEX3D(cctkGH,i,j,k)]
+ + (dt2) *
( ( phi[CCTK_GFINDEX3D(cctkGH,i+1,j,k)]
+phi[CCTK_GFINDEX3D(cctkGH,i-1,j,k)] )/dx2
+( phi[CCTK_GFINDEX3D(cctkGH,i,j+1,k)]
@@ -131,14 +87,12 @@ void WaveToyC_Evolution(CCTK_CARGUMENTS)
}
- /* Synchronize */
- CCTK_SyncGroup(cctkGH,"wavetoy::scalarevolve");
-
+ /* Synchronize before applying boundary conditions */
+ CCTK_SyncGroup(cctkGH,"wavetoyc::scalartmps");
/* Apply boundary conditions */
WaveToyC_Boundaries(CCTK_PASS_CTOC);
-
/* Update timeslices */
for (k=0; k<cctk_lsh[2]; k++)
{
@@ -146,10 +100,10 @@ void WaveToyC_Evolution(CCTK_CARGUMENTS)
{
for (i=0; i<cctk_lsh[0]; i++)
{
- phi_old[CCTK_GFINDEX3D(cctkGH,i,j,k)]
- = phi[CCTK_GFINDEX3D(cctkGH,i,j,k)];
- phi[CCTK_GFINDEX3D(cctkGH,i,j,k)]
- = phi_new[CCTK_GFINDEX3D(cctkGH,i,j,k)];
+
+ index = CCTK_GFINDEX3D(cctkGH,i,j,k);
+ phi_old[index] = phi[index];
+ phi[index] = phi_new[index];
}
}
}
@@ -157,6 +111,51 @@ void WaveToyC_Evolution(CCTK_CARGUMENTS)
}
+ /*@@
+ @routine WaveToyC_Boundaries
+ @date
+ @author Tom Goodale
+ @desc
+ Boundary conditions for the wave equation
+ @enddesc
+ @calls ApplySymmetry,ApplyFlatBC,ApplyRadiativeBC
+ @calledby
+ @history
+
+ @endhistory
+
+@@*/
+
+void WaveToyC_Boundaries(CCTK_CARGUMENTS)
+{
+ DECLARE_CCTK_CARGUMENTS
+ DECLARE_CCTK_PARAMETERS
+
+ int ierr;
+ int sw[3];
+
+ /* Set the stencil width */
+ sw[0]=1;
+ sw[1]=1;
+ sw[2]=1;
+
+ ApplySymmetry(cctkGH,"wavetoyc::scalartmps");
+
+ if (CCTK_EQUALS(bound,"flat"))
+ {
+ ierr = ApplyFlatBC(cctkGH,sw,"wavetoyc::phi_new");
+ }
+ else if (CCTK_Equals(bound,"radiation"))
+ {
+ ierr = ApplyRadiativeBC(cctkGH,0,sw,"wavetoyc::phi_new","wavetoy::phi");
+ }
+
+ if (ierr < 0)
+ {
+ CCTK_WARN(0,"Boundary conditions not applied - giving up!");
+ }
+
+}