diff options
Diffstat (limited to 'Examples/Advect/src/advect_flux.cc')
-rw-r--r-- | Examples/Advect/src/advect_flux.cc | 145 |
1 files changed, 78 insertions, 67 deletions
diff --git a/Examples/Advect/src/advect_flux.cc b/Examples/Advect/src/advect_flux.cc index dd27e27..1375eb2 100644 --- a/Examples/Advect/src/advect_flux.cc +++ b/Examples/Advect/src/advect_flux.cc @@ -12,50 +12,36 @@ #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 INV(x) ((CCTK_REAL)1.0 / (x)) #define SQR(x) ((x) * (x)) -#define CUB(x) ((x) * (x) * (x)) +#define CUB(x) ((x) * SQR(x)) +#define QAD(x) (SQR(SQR(x))) -static void advect_flux_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 min[3], int const max[3], int const n_subblock_gfs, CCTK_REAL * restrict const subblock_gfs[]) +static void advect_flux_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; - /* Declare the variables used for looping over grid points */ - CCTK_INT i, j, k; - // CCTK_INT index = INITVALUE; - - /* Declare finite differencing variables */ - - if (verbose > 1) - { - CCTK_VInfo(CCTK_THORNSTRING,"Entering advect_flux_Body"); - } - - if (cctk_iteration % advect_flux_calc_every != advect_flux_calc_offset) - { - return; - } - - const char *groups[] = {"Advect::Frho_group","Advect::rho_group","Advect::v_group"}; - GenericFD_AssertGroupStorage(cctkGH, "advect_flux", 3, groups); - - /* 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); @@ -69,24 +55,20 @@ static void advect_flux_Body(cGH const * restrict const cctkGH, int const dir, i CCTK_REAL const hdzi = 0.5 * dzi; /* Initialize predefined quantities */ - CCTK_REAL const p1o1 = 1; CCTK_REAL const p1o12dx = 0.0833333333333333333333333333333*INV(dx); CCTK_REAL const p1o12dy = 0.0833333333333333333333333333333*INV(dy); CCTK_REAL const p1o12dz = 0.0833333333333333333333333333333*INV(dz); - CCTK_REAL const p1o144dxdy = 0.00694444444444444444444444444444*INV(dx)*INV(dy); - CCTK_REAL const p1o144dxdz = 0.00694444444444444444444444444444*INV(dx)*INV(dz); - CCTK_REAL const p1o144dydz = 0.00694444444444444444444444444444*INV(dy)*INV(dz); + CCTK_REAL const p1o144dxdy = 0.00694444444444444444444444444444*INV(dx*dy); + CCTK_REAL const p1o144dxdz = 0.00694444444444444444444444444444*INV(dx*dz); + CCTK_REAL const p1o144dydz = 0.00694444444444444444444444444444*INV(dy*dz); 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 p1o4dxdy = 0.25*INV(dx)*INV(dy); - CCTK_REAL const p1o4dxdz = 0.25*INV(dx)*INV(dz); - CCTK_REAL const p1o4dydz = 0.25*INV(dy)*INV(dz); - CCTK_REAL const p1odx = INV(dx); + CCTK_REAL const p1o4dxdy = 0.25*INV(dx*dy); + CCTK_REAL const p1o4dxdz = 0.25*INV(dx*dz); + CCTK_REAL const p1o4dydz = 0.25*INV(dy*dz); CCTK_REAL const p1odx2 = INV(SQR(dx)); - CCTK_REAL const p1ody = INV(dy); CCTK_REAL const p1ody2 = INV(SQR(dy)); - CCTK_REAL const p1odz = INV(dz); CCTK_REAL const p1odz2 = INV(SQR(dz)); CCTK_REAL const pm1o12dx2 = -0.0833333333333333333333333333333*INV(SQR(dx)); CCTK_REAL const pm1o12dy2 = -0.0833333333333333333333333333333*INV(SQR(dy)); @@ -95,41 +77,47 @@ static void advect_flux_Body(cGH const * restrict const cctkGH, int const dir, i CCTK_REAL const pm1o2dy = -0.5*INV(dy); CCTK_REAL const pm1o2dz = -0.5*INV(dz); + /* Assign local copies of arrays functions */ + + + + /* Calculate temporaries and arrays functions */ + + /* Copy local copies back to grid functions */ + /* Loop over the grid points */ - for (k = min[2]; k < max[2]; k++) + #pragma omp parallel + CCTK_LOOP3(advect_flux, + i,j,k, imin[0],imin[1],imin[2], imax[0],imax[1],imax[2], + cctk_ash[0],cctk_ash[1],cctk_ash[2]) { - for (j = min[1]; j < max[1]; j++) - { - for (i = min[0]; i < max[0]; i++) - { - int const index = CCTK_GFINDEX3D(cctkGH,i,j,k) ; - - /* Assign local copies of grid functions */ - - CCTK_REAL rhoL = rho[index]; - CCTK_REAL v1L = v1[index]; - CCTK_REAL v2L = v2[index]; - CCTK_REAL v3L = v3[index]; - - - /* Include user supplied include files */ - - /* Precompute derivatives */ - - /* Calculate temporaries and grid functions */ - CCTK_REAL Frho1L = rhoL*v1L; - - CCTK_REAL Frho2L = rhoL*v2L; - - CCTK_REAL Frho3L = rhoL*v3L; - - /* Copy local copies back to grid functions */ - Frho1[index] = Frho1L; - Frho2[index] = Frho2L; - Frho3[index] = Frho3L; - } - } + ptrdiff_t const index = di*i + dj*j + dk*k; + + /* Assign local copies of grid functions */ + + CCTK_REAL rhoL = rho[index]; + CCTK_REAL v1L = v1[index]; + CCTK_REAL v2L = v2[index]; + CCTK_REAL v3L = v3[index]; + + + /* Include user supplied include files */ + + /* Precompute derivatives */ + + /* Calculate temporaries and grid functions */ + CCTK_REAL Frho1L = rhoL*v1L; + + CCTK_REAL Frho2L = rhoL*v2L; + + CCTK_REAL Frho3L = rhoL*v3L; + + /* Copy local copies back to grid functions */ + Frho1[index] = Frho1L; + Frho2[index] = Frho2L; + Frho3[index] = Frho3L; } + CCTK_ENDLOOP3(advect_flux); } extern "C" void advect_flux(CCTK_ARGUMENTS) @@ -137,5 +125,28 @@ extern "C" void advect_flux(CCTK_ARGUMENTS) DECLARE_CCTK_ARGUMENTS; DECLARE_CCTK_PARAMETERS; - GenericFD_LoopOverEverything(cctkGH, &advect_flux_Body); + + if (verbose > 1) + { + CCTK_VInfo(CCTK_THORNSTRING,"Entering advect_flux_Body"); + } + + if (cctk_iteration % advect_flux_calc_every != advect_flux_calc_offset) + { + return; + } + + const char *const groups[] = { + "Advect::Frho_group", + "Advect::rho_group", + "Advect::v_group"}; + GenericFD_AssertGroupStorage(cctkGH, "advect_flux", 3, groups); + + + GenericFD_LoopOverEverything(cctkGH, advect_flux_Body); + + if (verbose > 1) + { + CCTK_VInfo(CCTK_THORNSTRING,"Leaving advect_flux_Body"); + } } |