aboutsummaryrefslogtreecommitdiff
path: root/src
diff options
context:
space:
mode:
authorjthorn <jthorn@f88db872-0e4f-0410-b76b-b9085cfa78c5>2002-04-07 14:52:21 +0000
committerjthorn <jthorn@f88db872-0e4f-0410-b76b-b9085cfa78c5>2002-04-07 14:52:21 +0000
commitf067f5833d27c1449557183e4bac8c0d939670a5 (patch)
treeae32877e4af1323cf0d576f82d08ab1a251d80ea /src
parente8237486091a9fe61d970422d5ecbe15dc07119b (diff)
fix a bunch more compiler-found bugs (mostly syntax errors)
git-svn-id: http://svn.einsteintoolkit.org/cactus/EinsteinAnalysis/AHFinderDirect/trunk@460 f88db872-0e4f-0410-b76b-b9085cfa78c5
Diffstat (limited to 'src')
-rw-r--r--src/patch/coords.cc17
-rw-r--r--src/patch/patch_system.cc183
-rw-r--r--src/patch/patch_system.hh32
-rw-r--r--src/patch/patch_system_info.hh44
4 files changed, 146 insertions, 130 deletions
diff --git a/src/patch/coords.cc b/src/patch/coords.cc
index aa20418..d152b57 100644
--- a/src/patch/coords.cc
+++ b/src/patch/coords.cc
@@ -20,7 +20,6 @@
#include "jt/stdc.h"
#include "jt/util.hh"
-using jtutil::fuzzy;
#include "fp.hh"
#include "coords.hh"
@@ -45,12 +44,12 @@ namespace local_coords
bool fuzzy_EQ_ang(fp ang1, fp ang2)
{
-return fuzzy<fp>::is_integer( (ang2 - ang1)/(2.0*PI) );
+return jtutil::fuzzy<fp>::is_integer( (ang2 - ang1)/(2.0*PI) );
}
bool fuzzy_EQ_dang(fp dang1, fp dang2)
{
-return fuzzy<fp>::is_integer( (dang2 - dang1)/360.0 );
+return jtutil::fuzzy<fp>::is_integer( (dang2 - dang1)/360.0 );
}
}
@@ -95,12 +94,12 @@ void xyz_of_r_mu_nu(fp r, fp mu, fp nu, fp& x, fp& y, fp& z)
{
const fp sign_y = signum(sin(mu));
const fp sign_z_via_mu = signum(cos(mu));
-assert( fuzzy<fp>::NE(cos(mu), 0.0) );
+assert( jtutil::fuzzy<fp>::NE(cos(mu), 0.0) );
const fp y_over_z = tan(mu);
const fp sign_x = signum(sin(nu));
const fp sign_z_via_nu = signum(cos(nu));
-assert( fuzzy<fp>::NE(cos(nu), 0.0) );
+assert( jtutil::fuzzy<fp>::NE(cos(nu), 0.0) );
const fp x_over_z = tan(nu);
// failure of next assert() ==> inconsistent input (mu,nu)
@@ -131,12 +130,12 @@ const fp phi_bar = 0.5*PI - phi;
const fp sign_z = signum(sin(mu_bar));
const fp sign_y_via_mu_bar = signum(cos(mu_bar));
-assert( fuzzy<fp>::NE(cos(mu_bar), 0.0) );
+assert( jtutil::fuzzy<fp>::NE(cos(mu_bar), 0.0) );
const fp z_over_y = tan(mu_bar);
const fp sign_x = signum(sin(phi_bar));
const fp sign_y_via_phi_bar = signum(cos(phi_bar));
-assert( fuzzy<fp>::NE(cos(phi_bar), 0.0) );
+assert( jtutil::fuzzy<fp>::NE(cos(phi_bar), 0.0) );
const fp x_over_y = tan(phi_bar);
// failure of next assert() ==> inconsistent input (mu,phi)
@@ -166,12 +165,12 @@ const fp nu_bar = 0.5*PI - nu ;
const fp sign_z = signum(sin(nu_bar));
const fp sign_x_via_nu_bar = signum(cos(nu_bar));
-assert( fuzzy<fp>::NE(cos(nu_bar), 0.0) );
+assert( jtutil::fuzzy<fp>::NE(cos(nu_bar), 0.0) );
const fp z_over_x = tan(nu_bar);
const fp sign_y = signum(sin(phi));
const fp sign_x_via_phi = signum(cos(phi));
-assert( fuzzy<fp>::NE(cos(phi), 0.0) );
+assert( jtutil::fuzzy<fp>::NE(cos(phi), 0.0) );
const fp y_over_x = tan(phi );
// failure of next assert() ==> inconsistent input (nu,phi)
diff --git a/src/patch/patch_system.cc b/src/patch/patch_system.cc
index 2412d75..32a17fe 100644
--- a/src/patch/patch_system.cc
+++ b/src/patch/patch_system.cc
@@ -22,10 +22,13 @@
// patch_system::print_gridfn
//
-#include <stdio.h>
+#include <cstdio>
#include <assert.h>
-#include <math.h>
-#include <string.h>
+#include <cmath>
+#include <cstring>
+#include <vector>
+using std::printf;
+using std::strcmp;
#include "cctk.h"
@@ -76,8 +79,8 @@ patch_system::patch_system(fp origin_x_in, fp origin_y_in, fp origin_z_in,
enum patch_system_type type_in,
int N_ghost_points, int N_overlap_points,
fp delta_drho_dsigma,
- int min_gfn, int max_gfn,
- int ghosted_min_gfn, int ghosted_max_gfn,
+ int min_gfn_in, int max_gfn_in,
+ int ghosted_min_gfn_in, int ghosted_max_gfn_in,
int interp_handle, int interp_par_table_handle)
: global_coords_(origin_x_in, origin_y_in, origin_z_in),
@@ -85,7 +88,7 @@ patch_system::patch_system(fp origin_x_in, fp origin_y_in, fp origin_z_in,
N_patches_(N_patches_of_type(type_in)),
all_patches_(N_patches_),
gridfn_storage_(NULL), // set in setup_gridfn_storage()
- ghosted_gridfn_storage(NULL) // set in setup_gridfn_storage()
+ ghosted_gridfn_storage_(NULL) // set in setup_gridfn_storage()
{
if (! jtutil::is_odd(N_overlap_points))
then error_exit(ERROR_EXIT,
@@ -102,57 +105,52 @@ case full_sphere_patch_system:
construct_patches(patch_system_info::full_sphere::patch_info_array,
N_ghost_points, N_extend_points,
delta_drho_dsigma);
- setup_gridfn_storage(min_gfn, max_gfn,
- ghosted_min_gfn, ghosted_max_gfn);
- interlink_full_sphere_patch_system
- (N_ghost_points, N_extend_points,
- N_overlap_points, delta_drho_dsigma,
- interp_handle, interp_par_table_handle);
+ setup_gridfn_storage(min_gfn_in, max_gfn_in,
+ ghosted_min_gfn_in, ghosted_max_gfn_in);
+ interlink_full_sphere_patch_system(N_overlap_points,
+ interp_handle,
+ interp_par_table_handle);
break;
case plus_z_hemisphere_patch_system:
construct_patches(patch_system_info::plus_z_hemisphere
::patch_info_array,
N_ghost_points, N_extend_points,
delta_drho_dsigma);
- setup_gridfn_storage(min_gfn, max_gfn,
- ghosted_min_gfn, ghosted_max_gfn);
- interlink_plus_z_hemisphere_patch_system
- (N_ghost_points, N_extend_points,
- N_overlap_points, delta_drho_dsigma,
- interp_handle, interp_par_table_handle);
+ setup_gridfn_storage(min_gfn_in, max_gfn_in,
+ ghosted_min_gfn_in, ghosted_max_gfn_in);
+ interlink_plus_z_hemisphere_patch_system(N_overlap_points,
+ interp_handle,
+ interp_par_table_handle);
break;
case plus_xy_quadrant_patch_system:
construct_patches(patch_system_info::plus_xy_quadrant::patch_info_array,
N_ghost_points, N_extend_points,
delta_drho_dsigma);
- setup_gridfn_storage(min_gfn, max_gfn,
- ghosted_min_gfn, ghosted_max_gfn);
- interlink_plus_xy_quadrant_patch_system
- (N_ghost_points, N_extend_points,
- N_overlap_points, delta_drho_dsigma,
- interp_handle, interp_par_table_handle);
+ setup_gridfn_storage(min_gfn_in, max_gfn_in,
+ ghosted_min_gfn_in, ghosted_max_gfn_in);
+ interlink_plus_xy_quadrant_patch_system(N_overlap_points,
+ interp_handle,
+ interp_par_table_handle);
break;
case plus_xz_quadrant_patch_system:
construct_patches(patch_system_info::plus_xz_quadrant::patch_info_array,
N_ghost_points, N_extend_points,
delta_drho_dsigma);
- setup_gridfn_storage(min_gfn, max_gfn,
- ghosted_min_gfn, ghosted_max_gfn);
- interlink_plus_xz_quadrant_patch_system
- (N_ghost_points, N_extend_points,
- N_overlap_points, delta_drho_dsigma,
- interp_handle, interp_par_table_handle);
+ setup_gridfn_storage(min_gfn_in, max_gfn_in,
+ ghosted_min_gfn_in, ghosted_max_gfn_in);
+ interlink_plus_xz_quadrant_patch_system(N_overlap_points,
+ interp_handle,
+ interp_par_table_handle);
break;
case plus_xyz_octant_patch_system:
construct_patches(patch_system_info::plus_xyz_octant::patch_info_array,
N_ghost_points, N_extend_points,
delta_drho_dsigma);
- setup_gridfn_storage(min_gfn, max_gfn,
- ghosted_min_gfn, ghosted_max_gfn);
- interlink_plus_xyz_octant_patch_system
- (N_ghost_points, N_extend_points,
- N_overlap_points, delta_drho_dsigma,
- interp_handle, interp_par_table_handle);
+ setup_gridfn_storage(min_gfn_in, max_gfn_in,
+ ghosted_min_gfn_in, ghosted_max_gfn_in);
+ interlink_plus_xyz_octant_patch_system(N_overlap_points,
+ interp_handle,
+ interp_par_table_handle);
break;
default:
error_exit(ERROR_EXIT,
@@ -176,8 +174,6 @@ delete gridfn_storage_;
{
delete &ith_patch(pn);
}
-
-delete all_patches_;
}
//******************************************************************************
@@ -192,7 +188,7 @@ delete all_patches_;
// This function does *NOT* create any of the ghost zones or frontiers,
// and does *NOT* set up any gridfns.
//
-void patch_system::construct_patches(const struct patch_info patch_info[],
+void patch_system::construct_patches(const struct patch_info patch_info_in[],
int N_ghost_points, int N_extend_points,
fp delta_drho_dsigma)
{
@@ -200,7 +196,7 @@ N_grid_points_ = 0;
ghosted_N_grid_points_ = 0;
for (int pn = 0 ; pn < N_patches() ; ++pn)
{
- const struct patch_info& pi = patch_info[pn];
+ const struct patch_info& pi = patch_info_in[pn];
struct patch *p;
switch (pi.ctype)
@@ -233,7 +229,7 @@ ghosted_N_grid_points_ = 0;
default:
error_exit(ERROR_EXIT,
"***** patch_system::construct_patches():\n"
-" unknown patch_info[pn=%d].ctype=0x%02d='%c'!\n"
+" unknown patch_info_in[pn=%d].ctype=0x%02d='%c'!\n"
,
pn, pi.ctype, pi.ctype); /*NOTREACHED*/
}
@@ -278,37 +274,39 @@ ghosted_N_grid_points_ = 0;
// FIXME: We should pad the gridfn storage as necessary to avoid cache
// conflicts, but we don't do this at present.
//
-void setup_gridfn_storage(int min_gfn, int max_gfn,
- int ghosted_min_gfn, int ghosted_max_gfn)
+void patch_system::setup_gridfn_storage
+ (int min_gfn_in, int max_gfn_in,
+ int ghosted_min_gfn_in, int ghosted_max_gfn_in)
{
-const int N_gridfns = jtutil::how_many_in_range(min_gfn, max_gfn);
-const int ghosted_N_gridfns
- = jtutil::how_many_in_range(ghosted_min_gfn, ghosted_max_gfn);
+const int N_gridfns_in = jtutil::how_many_in_range(min_gfn_in, max_gfn_in);
+const int ghosted_N_gridfns_in
+ = jtutil::how_many_in_range(ghosted_min_gfn_in, ghosted_max_gfn_in);
const int gfn_stride = N_grid_points();
const int ghosted_gfn_stride = ghosted_N_grid_points();
-const int N_storage = gfn_stride * N_gridfns;
-const int ghosted_N_storage = ghosted_gfn_stride * ghosted_N_gridfns;
+const int N_storage = gfn_stride * N_gridfns_in;
+const int ghosted_N_storage = ghosted_gfn_stride * ghosted_N_gridfns_in;
// storage arrays for all gridfns
gridfn_storage_ = new fp[N_storage];
ghosted_gridfn_storage_ = new fp[ghosted_N_storage];
// divide up the storage array among the patches
+// and set up the storage in the individual patches themselves
int posn = 0;
int ghosted_posn = 0;
- for (pn = 0 ; pn < N_patches() ; ++pn)
+ for (int pn = 0 ; pn < N_patches() ; ++pn)
{
const struct grid_arrays::gridfn_pars gridfn_pars
= {
- min_gfn, max_gfn,
+ min_gfn_in, max_gfn_in,
& gridfn_storage_[posn],
gfn_stride, 0, 0
};
const struct grid_arrays::gridfn_pars ghosted_gridfn_pars
= {
- ghosted_min_gfn, ghosted_max_gfn,
+ ghosted_min_gfn_in, ghosted_max_gfn_in,
& ghosted_gridfn_storage_[ghosted_posn],
ghosted_gfn_stride, 0, 0
};
@@ -325,19 +323,19 @@ assert(ghosted_posn == ghosted_N_grid_points());
// check to make sure storage for distinct gridfns
// forms a partition of the overall storage array
const patch& pfirst = ith_patch(0);
-const patcn& plast = ith_patch(N_patches()-1);
+const patch& plast = ith_patch(N_patches()-1);
for (int gfn = min_gfn() ; gfn+1 < max_gfn() ; ++gfn)
{
// range of storage occupied by gridfns:
// gfn --> [gfn_first, gfn_last]
// gfn+1 --> [gfn1_first, gfn1_last]
- const fp* const gfn_last
+ const fp* const gfn_last_ptr
= & plast.gridfn(gfn, plast.max_irho(),
plast.max_isigma());
- const fp* const gfn1_first
+ const fp* const gfn1_first_ptr
= & pfirst.gridfn(gfn+1, pfirst.min_irho(),
pfirst.min_isigma());
- assert(gfn1_first == gfn_last + 1);
+ assert(gfn1_first_ptr == gfn_last_ptr + 1);
}
for (int gfn = ghosted_min_gfn() ; gfn+1 < ghosted_max_gfn() ; ++gfn)
@@ -345,15 +343,15 @@ const patcn& plast = ith_patch(N_patches()-1);
// range of storage occupied by ghosted gridfns:
// gfn --> [gfn_first, gfn_last]
// gfn+1 --> [gfn1_first, gfn1_last]
- const fp* const gfn1_last
+ const fp* const gfn_last_ptr
= & plast.ghosted_gridfn(gfn,
plast.ghosted_max_irho(),
plast.ghosted_max_isigma());
- const fp* const gfn1_first
+ const fp* const gfn1_first_ptr
= & pfirst.ghosted_gridfn(gfn+1,
pfirst.ghosted_min_irho(),
pfirst.ghosted_min_isigma());
- assert(gnosted_gfn1_first == gfn_last + 1);
+ assert(gfn1_first_ptr == gfn_last_ptr + 1);
}
// check to make sure storage for distinct patches
@@ -368,11 +366,11 @@ const patcn& plast = ith_patch(N_patches()-1);
// range of storage occupied by gridfn:
// p --> [p_first, p_last]
// p1 --> [p1_first, p1_last]
- const fp* const p_last = & p.gridfn(gfn, p.max_irho(),
- p.max_isigma());
- const fp* const p1_first = & p1.gridfn(gfn, p1.min_irho(),
- p1.min_isigma());
- assert(p1_first == p_last + 1);
+ const fp* const p_last_ptr
+ = & p.gridfn(gfn, p.max_irho(), p.max_isigma());
+ const fp* const p1_first_ptr
+ = & p1.gridfn(gfn, p1.min_irho(), p1.min_isigma());
+ assert(p1_first_ptr == p_last_ptr + 1);
}
}
@@ -386,13 +384,13 @@ const patcn& plast = ith_patch(N_patches()-1);
// range of storage occupied by ghosted gridfn:
// p --> [p_first, p_last]
// p1 --> [p1_first, p1_last]
- const fp* const p_last
+ const fp* const p_last_ptr
= & p.ghosted_gridfn(gfn, p.ghosted_max_irho(),
p.ghosted_max_isigma());
- const fp* const p1_first
+ const fp* const p1_first_ptr
= & p1.ghosted_gridfn(gfn, p1.ghosted_min_irho(),
p1.ghosted_min_isigma());
- assert(p1_first == p_last + 1);
+ assert(p1_first_ptr == p_last_ptr + 1);
}
}
}
@@ -404,10 +402,16 @@ const patcn& plast = ith_patch(N_patches()-1);
// frontiers for a full-sphere patch system.
//
void patch_system::interlink_full_sphere_patch_system
- (int N_ghost_points, int N_extend_points,
- int N_overlap_points, fp delta_drho_dsigma,
+ (int N_overlap_points,
int interp_handle, int interp_par_table_handle)
{
+patch& pz = ith_patch(patch_number_of_name("+z"));
+patch& px = ith_patch(patch_number_of_name("+x"));
+patch& py = ith_patch(patch_number_of_name("+y"));
+patch& mx = ith_patch(patch_number_of_name("-x"));
+patch& my = ith_patch(patch_number_of_name("-y"));
+patch& mz = ith_patch(patch_number_of_name("-z"));
+
// create the ghost zones
setup_adjacent_ghost_zones(pz, px, N_overlap_points);
setup_adjacent_ghost_zones(pz, py, N_overlap_points);
@@ -470,10 +474,14 @@ assert_all_ghost_zones_fully_setup();
// frontiers for a +z hemisphere patch system.
//
void patch_system::interlink_plus_z_hemisphere_patch_system
- (int N_ghost_points, int N_extend_points,
- int N_overlap_points, fp delta_drho_dsigma,
+ (int N_overlap_points,
int interp_handle, int interp_par_table_handle)
{
+patch& pz = ith_patch(patch_number_of_name("+z"));
+patch& px = ith_patch(patch_number_of_name("+x"));
+patch& py = ith_patch(patch_number_of_name("+y"));
+patch& mx = ith_patch(patch_number_of_name("-x"));
+patch& my = ith_patch(patch_number_of_name("-y"));
// create the ghost zones
setup_adjacent_ghost_zones(pz, px, N_overlap_points);
@@ -525,10 +533,13 @@ assert_all_ghost_zones_fully_setup();
// frontiers for a +xy quadrant patch system.
//
void patch_system::interlink_plus_xy_quadrant_patch_system
- (int N_ghost_points, int N_extend_points,
- int N_overlap_points, fp delta_drho_dsigma,
+ (int N_overlap_points,
int interp_handle, int interp_par_table_handle)
{
+patch& pz = ith_patch(patch_number_of_name("+z"));
+patch& px = ith_patch(patch_number_of_name("+x"));
+patch& py = ith_patch(patch_number_of_name("+y"));
+patch& mz = ith_patch(patch_number_of_name("-z"));
// create the ghost zones
setup_adjacent_ghost_zones(pz, px, N_overlap_points);
@@ -573,10 +584,14 @@ assert_all_ghost_zones_fully_setup();
// frontiers for a +xz quadrant patch system.
//
void patch_system::interlink_plus_xz_quadrant_patch_system
- (int N_ghost_points, int N_extend_points,
- int N_overlap_points, fp delta_drho_dsigma,
+ (int N_overlap_points,
int interp_handle, int interp_par_table_handle)
{
+patch& pz = ith_patch(patch_number_of_name("+z"));
+patch& px = ith_patch(patch_number_of_name("+x"));
+patch& py = ith_patch(patch_number_of_name("+y"));
+patch& my = ith_patch(patch_number_of_name("-y"));
+
// create the ghost zones
setup_adjacent_ghost_zones(pz, px, N_overlap_points);
setup_adjacent_ghost_zones(pz, py, N_overlap_points);
@@ -620,10 +635,13 @@ assert_all_ghost_zones_fully_setup();
// frontiers for a +xyz octant patch system.
//
void patch_system::interlink_plus_xyz_octant_patch_system
- (int N_ghost_points, int N_extend_points,
- int N_overlap_points, fp delta_drho_dsigma,
+ (int N_overlap_points,
int interp_handle, int interp_par_table_handle)
{
+patch& pz = ith_patch(patch_number_of_name("+z"));
+patch& px = ith_patch(patch_number_of_name("+x"));
+patch& py = ith_patch(patch_number_of_name("+y"));
+
// create the ghost zones
setup_adjacent_ghost_zones(pz, px, N_overlap_points);
setup_adjacent_ghost_zones(pz, py, N_overlap_points);
@@ -846,8 +864,8 @@ error_exit(ERROR_EXIT,
// This process is described in detail in the header comments in
// "ghost_zone.hh".
//
-void patch_system::synchronize_ghost_zones(int ghosted_min_gfn,
- int ghosted_max_gfn)
+void patch_system::synchronize_ghost_zones(int ghosted_min_gfn_to_sync,
+ int ghosted_max_gfn_to_sync)
{
//
// Phase 1:
@@ -866,7 +884,8 @@ void patch_system::synchronize_ghost_zones(int ghosted_min_gfn,
{
ghost_zone& gz = p.minmax_ang_ghost_zone(want_min, want_rho);
if (gz.is_symmetry())
- then gz.synchronize(ghosted_min_gfn, ghosted_max_gfn,
+ then gz.synchronize(ghosted_min_gfn_to_sync,
+ ghosted_max_gfn_to_sync,
false, // want min-par corner?
true, // want non-corner?
false); // want max-par corner?
@@ -890,7 +909,8 @@ void patch_system::synchronize_ghost_zones(int ghosted_min_gfn,
{
ghost_zone& gz = p.minmax_ang_ghost_zone(want_min, want_rho);
if (gz.is_interpatch())
- then gz.synchronize(ghosted_min_gfn, ghosted_max_gfn,
+ then gz.synchronize(ghosted_min_gfn_to_sync,
+ ghosted_max_gfn_to_sync,
true, // want min-par corner?
true, // want non-corner?
true); // want max-par corner?
@@ -915,7 +935,8 @@ void patch_system::synchronize_ghost_zones(int ghosted_min_gfn,
{
ghost_zone& gz = p.minmax_ang_ghost_zone(want_min, want_rho);
if (gz.is_symmetry())
- then gz.synchronize(ghosted_min_gfn, ghosted_max_gfn,
+ then gz.synchronize(ghosted_min_gfn_to_sync,
+ ghosted_max_gfn_to_sync,
true, // want min-par corner?
false, // want non-corner?
true); // want max-par corner?
diff --git a/src/patch/patch_system.hh b/src/patch/patch_system.hh
index 1db72db..140675e 100644
--- a/src/patch/patch_system.hh
+++ b/src/patch/patch_system.hh
@@ -94,7 +94,7 @@ public:
// get patches by patch number
patch &ith_patch(int pn) const
- { return * (*all_patches_)(pn); }
+ { return * all_patches_[pn]; }
// find a patch by name, return patch number; error_exit() if not found
int patch_number_of_name(const char* name) const;
@@ -125,8 +125,8 @@ public:
// i.e. update the ghost-zone values of the specified gridfns
// via the appropriate sequence of symmetry operations
// and interpatch interpolations
- void synchronize_ghost_zones(int ghosted_min_gfn,
- int ghosted_max_gfn);
+ void synchronize_ghost_zones(int ghosted_min_gfn_to_sync,
+ int ghosted_max_gfn_to_sync);
// print a gridfn (via C stdio) in ASCII format
// output is assumed to already be open
@@ -150,8 +150,8 @@ public:
enum patch_system_type type_in,
int N_ghost_points, int N_overlap_points,
fp delta_drho_dsigma,
- int min_gfn, max_gfn,
- int ghosted_min_gfn, int ghosted_max_gfn);
+ int min_gfn_in, int max_gfn_in,
+ int ghosted_min_gfn_in, int ghosted_max_gfn_in,
int interp_handle_in, int interp_par_table_handle_in);
~patch_system();
@@ -164,34 +164,30 @@ private:
// and link them into the patch system
// does *NOT* create ghost zones or frontiers
// does *NOT* set up gridfns
- void construct_patches(const struct patch_info patch_info[],
+ void construct_patches(const struct patch_info patch_info_in[],
int N_ghost_points, int N_extend_points,
fp delta_drho_dsigma);
// setup all gridfns with contiguous-across-patches storage
- void setup_gridfn_storage(int min_gfn, int max_gfn,
- int ghosted_min_gfn, int ghosted_max_gfn);
+ void setup_gridfn_storage
+ (int min_gfn_in, int max_gfn_in,
+ int ghosted_min_gfn_in, int ghosted_max_gfn_in);
// create/interlink ghost zones and frontiers
void interlink_full_sphere_patch_system
- (int N_ghost_points, int N_extend_points,
- int N_overlap_points, fp delta_drho_dsigma,
+ (int N_overlap_points,
int interp_handle, int interp_par_table_handle);
void interlink_plus_z_hemisphere_patch_system
- (int N_ghost_points, int N_extend_points,
- int N_overlap_points, fp delta_drho_dsigma,
+ (int N_overlap_points,
int interp_handle, int interp_par_table_handle);
void interlink_plus_xy_quadrant_patch_system
- (int N_ghost_points, int N_extend_points,
- int N_overlap_points, fp delta_drho_dsigma,
+ (int N_overlap_points,
int interp_handle, int interp_par_table_handle);
void interlink_plus_xz_quadrant_patch_system
- (int N_ghost_points, int N_extend_points,
- int N_overlap_points, fp delta_drho_dsigma,
+ (int N_overlap_points,
int interp_handle, int interp_par_table_handle);
void interlink_plus_xyz_octant_patch_system
- (int N_ghost_points, int N_extend_points,
- int N_overlap_points, fp delta_drho_dsigma,
+ (int N_overlap_points,
int interp_handle, int interp_par_table_handle);
// create/interlink a pair of periodic-symmetry ghost zones
diff --git a/src/patch/patch_system_info.hh b/src/patch/patch_system_info.hh
index a4924b3..72b0191 100644
--- a/src/patch/patch_system_info.hh
+++ b/src/patch/patch_system_info.hh
@@ -14,22 +14,22 @@ namespace full_sphere
static const struct patch_info patch_info_array[]
= {
// +z patch (90 x 90 degrees): dmu [ -45, 45], dnu [ -45, 45]
- {"+z", patch_is_plus, 'z', -45.0, 45.0, -45.0, 45.0},
+ {"+z", patch::patch_is_plus, 'z', -45.0, 45.0, -45.0, 45.0},
// +x patch (90 x 90 degrees): dnu [ 45, 135], dphi [ -45, 45]
- {"+x", patch_is_plus, 'x', 45.0, 135.0, -45.0, 45.0},
+ {"+x", patch::patch_is_plus, 'x', 45.0, 135.0, -45.0, 45.0},
// +y patch (90 x 90 degrees): dmu [ 45, 135], dphi [ 45, 135]
- {"+y", patch_is_plus, 'y', 45.0, 135.0, 45.0, 135.0},
+ {"+y", patch::patch_is_plus, 'y', 45.0, 135.0, 45.0, 135.0},
// -x patch (90 x 90 degrees): dnu [-135, -45], dphi [ 135, 225]
- {"-x", patch_is_minus, 'x', -135.0, -45.0, 135.0, 225.0},
+ {"-x", patch::patch_is_minus, 'x', -135.0, -45.0, 135.0, 225.0},
// -y patch (90 x 90 degrees): dmu [-135, -45], dphi [-135, -45]
- {"-y", patch_is_minus, 'y', -135.0, -45.0, -135.0, -45.0},
+ {"-y", patch::patch_is_minus, 'y', -135.0, -45.0, -135.0, -45.0},
// -z patch (90 x 90 degrees): dmu [ 135, 225], dnu [ 135, 225]
- {"-z", patch_is_minus, 'z', 135.0, 225.0, 135.0, 225.0},
+ {"-z", patch::patch_is_minus, 'z', 135.0, 225.0, 135.0, 225.0},
};
static const int N_patches = sizeof(patch_info_array) / sizeof(patch_info);
}; // namespace patch_system_info::full_sphere
@@ -43,19 +43,19 @@ namespace plus_z_hemisphere
static const struct patch_info patch_info_array[]
= {
// +z patch (90 x 90 degrees): dmu [ -45, 45], dnu [ -45, 45]
- {"+z", patch_is_plus, 'z', -45.0, 45.0, -45.0, 45.0},
+ {"+z", patch::patch_is_plus, 'z', -45.0, 45.0, -45.0, 45.0},
// +x patch (45 x 90 degrees): dnu [ 45, 90], dphi [ -45, 45]
- {"+x", patch_is_plus, 'x', 45.0, 90.0, -45.0, 45.0},
+ {"+x", patch::patch_is_plus, 'x', 45.0, 90.0, -45.0, 45.0},
// +y patch (45 x 90 degrees): dmu [ 45, 90], dphi [ 45, 135]
- {"+y", patch_is_plus, 'y', 45.0, 90.0, 45.0, 135.0},
+ {"+y", patch::patch_is_plus, 'y', 45.0, 90.0, 45.0, 135.0},
// -x patch (45 x 90 degrees): dnu [ -90, -45], dphi [ 135, 225]
- {"-x", patch_is_minus, 'x', -90.0, -45.0, 135.0, 225.0},
+ {"-x", patch::patch_is_minus, 'x', -90.0, -45.0, 135.0, 225.0},
// -y patch (45 x 90 degrees): dmu [ -90, -45], dphi [-135, -45]
- {"-y", patch_is_minus, 'y', -90.0, -45.0, -135.0, -45.0},
+ {"-y", patch::patch_is_minus, 'y', -90.0, -45.0, -135.0, -45.0},
};
static const int N_patches = sizeof(patch_info_array) / sizeof(patch_info);
}; // namespace patch_system_info::plus_z_hemisphere
@@ -69,16 +69,16 @@ namespace plus_xy_quadrant
static const struct patch_info patch_info_array[]
= {
// +z patch (45 x 45 degrees): dmu [ 0, 45], dnu [ 0, 45]
- {"+z", patch_is_plus, 'z', 0.0, 45.0, 0.0, 45.0},
+ {"+z", patch::patch_is_plus, 'z', 0.0, 45.0, 0.0, 45.0},
// +x patch (90 x 45 degrees): dnu [ 45, 135], dphi [ 0, 45]
- {"+x", patch_is_plus, 'x', 45.0, 135.0, 0.0, 45.0},
+ {"+x", patch::patch_is_plus, 'x', 45.0, 135.0, 0.0, 45.0},
// +y patch (90 x 45 degrees): dmu [ 45, 135], dphi [ 45, 90]
- {"+y", patch_is_plus, 'y', 45.0, 135.0, 45.0, 90.0},
+ {"+y", patch::patch_is_plus, 'y', 45.0, 135.0, 45.0, 90.0},
// -z patch (45 x 45 degrees): dmu [ 135, 180], dnu [ 135, 180]
- {"-z", patch_is_minus, 'z', 135.0, 180.0, 135.0, 180.0},
+ {"-z", patch::patch_is_minus, 'z', 135.0, 180.0, 135.0, 180.0},
};
static const int N_patches = sizeof(patch_info_array) / sizeof(patch_info);
}; // namespace patch_system_info::plus_xy_quadrant
@@ -93,16 +93,16 @@ namespace plus_xz_quadrant
static const struct patch_info patch_info_array[]
= {
// +z patch (90 x 45 degrees): dmu [ -45, 45], dnu [ 0, 45]
- {"+z", patch_is_plus, 'z', -45.0, 45.0, 0.0, 45.0},
+ {"+z", patch::patch_is_plus, 'z', -45.0, 45.0, 0.0, 45.0},
// +x patch (45 x 90 degrees): dnu [ 45, 90], dphi [ -45, 45]
- {"+x", patch_is_plus, 'x', 45.0, 90.0, -45.0, 45.0},
+ {"+x", patch::patch_is_plus, 'x', 45.0, 90.0, -45.0, 45.0},
// +y patch (45 x 45 degrees): dmu [ 45, 90], dphi [ 45, 90]
- {"+y", patch_is_plus, 'y', 45.0, 90.0, 45.0, 90.0},
+ {"+y", patch::patch_is_plus, 'y', 45.0, 90.0, 45.0, 90.0},
// -y patch (45 x 45 degrees): dmu [ -90, -45], dphi [ -90, -45]
- {"-y", patch_is_minus, 'y', -90.0, -45.0, -90.0, -45.0},
+ {"-y", patch::patch_is_minus, 'y', -90.0, -45.0, -90.0, -45.0},
};
static const int N_patches = sizeof(patch_info_array) / sizeof(patch_info);
}; // namespace patch_system_info::plus_xz_quadrant
@@ -117,11 +117,11 @@ namespace plus_xyz_octant
static const struct patch_info patch_info_array[]
= {
// +z patch (45 x 45 degrees): dmu [ 0, 45], dnu [ 0, 45]
- {"+z", patch_is_plus, 'z', 0.0, 45.0, 0.0, 45.0},
+ {"+z", patch::patch_is_plus, 'z', 0.0, 45.0, 0.0, 45.0},
// +x patch (45 x 45 degrees): dnu [ 45, 90], dphi [ 0, 45]
- {"+x", patch_is_plus, 'x', 45.0, 90.0, 0.0, 45.0},
+ {"+x", patch::patch_is_plus, 'x', 45.0, 90.0, 0.0, 45.0},
// +y patch (45 x 45 degrees): dmu [ 45, 90], dphi [ 45, 90]
- {"+y", patch_is_plus, 'y', 45.0, 90.0, 45.0, 90.0},
+ {"+y", patch::patch_is_plus, 'y', 45.0, 90.0, 45.0, 90.0},
};
static const int N_patches = sizeof(patch_info_array) / sizeof(patch_info);
}; // namespace patch_system_info::octant