aboutsummaryrefslogtreecommitdiff
path: root/src
diff options
context:
space:
mode:
authorjthorn <jthorn@f88db872-0e4f-0410-b76b-b9085cfa78c5>2002-04-01 14:34:29 +0000
committerjthorn <jthorn@f88db872-0e4f-0410-b76b-b9085cfa78c5>2002-04-01 14:34:29 +0000
commit87120814659e7fc82c81e16559f92bfdd2e0acb9 (patch)
treeee05f79f9737322f945ceb8734056630a8b0ddbf /src
parent33f5afc59239cfcac04ae6c310e5517832ab354f (diff)
a few more mods -- looks ok now (ready to try compiling)
git-svn-id: http://svn.einsteintoolkit.org/cactus/EinsteinAnalysis/AHFinderDirect/trunk@414 f88db872-0e4f-0410-b76b-b9085cfa78c5
Diffstat (limited to 'src')
-rw-r--r--src/patch/patch_interp.cc71
-rw-r--r--src/patch/patch_interp.hh20
2 files changed, 62 insertions, 29 deletions
diff --git a/src/patch/patch_interp.cc b/src/patch/patch_interp.cc
index 4e6ef58..0238aa4 100644
--- a/src/patch/patch_interp.cc
+++ b/src/patch/patch_interp.cc
@@ -4,6 +4,7 @@
//
// patch_frontier::patch_frontier
// patch_frontier::~patch_frontier
+// patch_frontier::[min,max]_ipar_for_gridfn_data
// patch_frontier::interpolate
//
@@ -55,16 +56,8 @@ patch_frontier::patch_frontier(const patch_edge& my_edge_in,
: my_patch_(my_edge_in.my_patch()),
my_edge_(my_edge_in),
min_iperp_(min_iperp_in), max_iperp_(max_iperp_in),
- min_ipar_(my_patch()
- .ghost_zone_on_edge(my_edge_in.min_par_adjacent_edge())
- .is_symmetry()
- ? my_edge_in.min_ipar_with_corners()
- : my_edge_in.min_ipar_without_corners()),
- max_ipar_(my_patch()
- .ghost_zone_on_edge(my_edge_in.max_par_adjacent_edge())
- .is_symmetry()
- ? my_edge_in.max_ipar_with_corners()
- : my_edge_in.max_ipar_without_corners())
+ min_ipar_(min_ipar_for_gridfn_data()),
+ max_ipar_(max_ipar_for_gridfn_data()),
min_parindex_used_(min_parindex_used_in),
max_parindex_used_(max_parindex_used_in),
interp_par_(interp_par_in),
@@ -73,8 +66,10 @@ patch_frontier::patch_frontier(const patch_edge& my_edge_in,
interp_handle_(interp_handle_in),
interp_par_table_handle_(Util_TableClone(interp_par_table_handle_in)),
gridfn_coord_origin_(my_edge().par_of_ipar(min_ipar())),
- gridfn_coord_delta_(my_edge().par_map().delta_fp()),
- gridfn_type_codes(ghosted_min_gfn_in, ghosted_max_gfn_in)
+ gridfn_coord_delta_ (my_edge().par_map().delta_fp()),
+ gridfn_type_codes_ (ghosted_min_gfn_in, ghosted_max_gfn_in)
+ gridfn_data_ptrs_ (ghosted_min_gfn_in, ghosted_max_gfn_in)
+ result_buffer_ptrs_(ghosted_min_gfn_in, ghosted_max_gfn_in)
{
int status;
@@ -89,15 +84,7 @@ if (status < 0)
interp_par_table_handle_in,
status); /*NOTREACHED*/
-// set up CCTK_VARIABLE_* type codes for gridfns
- for (int gfn = ghosted_min_gfn() ; gfn <= ghosted_max_gfn() ; ++gfn)
- {
- gridfn_type_codes_(gfn) = CCTK_VARIABLE_REAL; // fp == CCTK_REAL
- }
-
-//
// set up parameter table for interpolator
-//
const int N_dims = 1;
const CCTK_INT stride = my_patch().ghosted_iang_stride(my_edge().is_rho());
status = Util_TableSetIntArray(interp_par_table_handle_,
@@ -111,6 +98,12 @@ if (status < 0)
,
interp_par_table_handle_,
status); /*NOTREACHED*/
+
+// set up gridfn type codes
+ for (int gfn = ghosted_min_gfn() ; gfn <= ghosted_max_gfn() ; ++gfn)
+ {
+ gridfn_type_codes_(gfn) = CCTK_VARIABLE_REAL; // fp == CCTK_REAL
+ }
}
//*****************************************************************************
@@ -126,13 +119,43 @@ Util_TableDestroy(interp_par_table_handle_);
//*****************************************************************************
//
+// These functions compute the [min,max] range of ipar in this patch
+// from which we may use gridfn data in interpolation. This computation
+// depends on our adjacent ghost zones' types as described in the comments
+// to our constructor, but needs only my_patch_ and my_edge_ members of
+// this object to be valid.
+//
+int patch_frontier::min_ipar_for_gridfn_data() const
+{
+const ghost_zone& min_par_adjacent_ghost_zone
+ = my_patch()
+ .ghost_zone_on_edge(my_edge().min_par_adjacent_edge());
+return min_par_adjacent_ghost_zone.is_symmetry()
+ ? my_edge_in.min_ipar_with_corners()
+ : my_edge_in.min_ipar_without_corners();
+}
+
+int patch_frontier::max_ipar_for_gridfn_data() const
+{
+const ghost_zone& max_par_adjacent_ghost_zone
+ = my_patch()
+ .ghost_zone_on_edge(my_edge().max_par_adjacent_edge());
+return max_par_adjacent_ghost_zone.is_symmetry()
+ ? my_edge_in.max_ipar_with_corners()
+ : my_edge_in.max_ipar_without_corners();
+}
+
+//*****************************************************************************
+
+//
// This function interpolates the specified range of ghosted gridfns
// at all the coordinates specified when we were constructed.
// It stores the interpolated data in result_buffer(gfn, iperp,parindex) .
// It calls error_exit() if the point is out of range or any other
// interpolator error is detected.
//
-void patch_frontier::interpolate(int ghosted_min_gfn, ghosted_max_gfn,
+void patch_frontier::interpolate(int ghosted_min_gfn_to_interp,
+ int ghosted_max_gfn_to_interp,
jtutil::array3d<fp>& result_buffer) const
{
//
@@ -153,9 +176,9 @@ void patch_frontier::interpolate(int ghosted_min_gfn, ghosted_max_gfn,
const int N_dims = 1;
const int N_gridfns = jtutil::how_many_in_range(ghosted_min_gfn_to_interp,
ghosted_max_gfn_to_interp);
+const CCTK_INT N_gridfn_data_points = jtutil::how_many_in_range(min_ipar(),
+ max_ipar());
const int interp_coords_type_code = CCTK_VARIABLE_REAL; // fp == CCTK_REAL
-const CCTK_INT gridfn_data_dims[N_dims]
- = { jtutil::how_many_in_range(min_ipar(), max_ipar()) };
const CCTK_INT *gridfn_type_codes_ptr
= & gridfn_type_codes_(ghosted_min_gfn_to_interp);
@@ -207,7 +230,7 @@ const CCTK_INT *gridfn_type_codes_ptr
interp_coords_type_code,
interp_coords,
N_gridfns,
- gridfn_data_dims,
+ &N_gridfn_data_points,
gridfn_type_codes_ptr,
gridfn_data,
N_gridfns,
diff --git a/src/patch/patch_interp.hh b/src/patch/patch_interp.hh
index 0422f21..4644c5b 100644
--- a/src/patch/patch_interp.hh
+++ b/src/patch/patch_interp.hh
@@ -75,10 +75,22 @@ public:
// interpolate specified range of ghosted gridfns
// at all the coordinates specified when we were constructed,
// store interpolated data in result_buffer(gfn, iperp,parindex)
- void interpolate(int ghosted_min_gfn, ghosted_max_gfn,
+ void interpolate(int ghosted_min_gfn_to_interp,
+ int ghosted_max_gfn_to_interp,
jtutil::array3d<fp>& result_buffer);
+
+ //
+ // ***** internal functions *****
+ //
private:
+ // range of ipar in this patch from which we may use gridfn data
+ // in interpolation -- depends on adjacent ghost zones' types as
+ // described in comments to constructor (above), but needs only
+ // my_patch_ and my_edge_ members of this object to be valid
+ int min_ipar_for_gridfn_data() const;
+ int max_ipar_for_gridfn_data() const;
+
//
// ***** constructor, destructor *****
@@ -143,11 +155,9 @@ private:
const patch& my_patch_;
const patch_edge& my_edge_;
+ // range of (iperp,ipar) in this patch from which
+ // we will use gridfn data in interpolation
int min_iperp_, max_iperp_;
-
- // range of ipar in this patch from which we may use gridfn
- // data in interpolation -- depends on adjacent ghost zones'
- // type as described in comments to constructor (above)
int min_ipar_, max_ipar_;
// [min,max] parindex at each iperp