aboutsummaryrefslogtreecommitdiff
diff options
context:
space:
mode:
authorreisswig <reisswig@40f6ab95-0e4f-0410-8daa-ee8d7420be1d>2007-11-22 11:22:55 +0000
committerreisswig <reisswig@40f6ab95-0e4f-0410-8daa-ee8d7420be1d>2007-11-22 11:22:55 +0000
commit0bb9874d57f57a0c5bbfbf640d995235ceb77932 (patch)
treebf958cf6e92c3da36b4053e00b68d750335107ed
parent7bd4b7c01aa9104b14540bcccf31ccc143dd58f8 (diff)
- introduced warn messages instead of asserts.
- limit resolution by maxntheta/maxnphi if reached git-svn-id: http://svn.cactuscode.org/arrangements/CactusNumerical/SphericalSurface/trunk@32 40f6ab95-0e4f-0410-8daa-ee8d7420be1d
-rw-r--r--src/setup.cc41
1 files changed, 39 insertions, 2 deletions
diff --git a/src/setup.cc b/src/setup.cc
index dc24e82..897efb3 100644
--- a/src/setup.cc
+++ b/src/setup.cc
@@ -101,8 +101,26 @@ extern "C" void SphericalSurface_Setup (CCTK_ARGUMENTS)
sf_nphi[n] += 4 - (sf_nphi[n]-2*nghostsphi[n]) % 4 + 2*nghostsphi[n];
// consistency check
- assert (sf_ntheta[n] >= 3*nghoststheta[n] && sf_ntheta[n] <= maxntheta);
- assert (sf_nphi[n] >= 3*nghostsphi[n] && sf_nphi[n] <= maxnphi);
+ if (sf_ntheta[n] < 3*nghoststheta[n])
+ {
+ CCTK_WARN(1, "AutoRes: Determined sf_ntheta[n] < 3*nghoststheta[n]! Setting sf_ntheta[n] = 3*nghoststheta[n]!");
+ sf_ntheta[n] = 3*nghoststheta[n];
+ }
+ if (sf_ntheta[n] > maxntheta)
+ {
+ CCTK_WARN(1, "AutoRes: Determined sf_ntheta[n] > maxntheta! You should increase maxntheta or auto_res_ratio[n]! Setting sf_ntheta[n] = maxntheta!");
+ sf_ntheta[n] = maxntheta;
+ }
+ if (sf_nphi[n] < 3*nghostsphi[n])
+ {
+ CCTK_WARN(1, "AutoRes: Determined sf_nphi[n] < 3*nghostsphi[n]! Setting sf_nphi[n] = 3*nghostsphi[n]!");
+ sf_nphi[n] = 3*nghostsphi[n];
+ }
+ if (sf_nphi[n] > maxnphi)
+ {
+ CCTK_WARN(1, "AutoRes: Determined sf_nphi[n] > maxnphi! You should increase maxnphi or auto_res_ratio[n]! Setting sf_nphi[n] = maxnphi!");
+ sf_nphi[n] = maxnphi;
+ }
// check, if there are enough maxntheta,maxnphi gridpoints to remove symmetries
// (some thorns require that)
@@ -338,6 +356,25 @@ extern "C" void SphericalSurface_SetupRes (CCTK_ARGUMENTS)
dummy = 0;
MPI_Allreduce(&sf_delta_phi_estimate[n], &dummy, 1, MPI_DOUBLE, MPI_MIN, MPI_COMM_WORLD);
sf_delta_phi_estimate[n] = dummy;
+
+ /*
+ // input arrays and output values
+ const CCTK_INT input_array_variable_indices[2]
+ = { CCTK_VarIndex("SphericalSurface::sf_delta_theta_estimate"),
+ CCTK_VarIndex("SphericalSurface::sf_delta_phi_estimate") };
+ const CCTK_INT output_value_type_codes[2]
+ = { CCTK_VARIABLE_REAL, CCTK_VARIABLE_REAL };
+ void *const output_numbers[2]
+ = { (void *) output_for_real_values,
+ (void *) output_for_real_values };
+
+ const int status
+ = CCTK_ReduceGridArrays(GH,
+ 0,
+ 2, input_array_variable_indices,
+ 2, output_value_type_codes,
+ output_values);
+ */
}
// find maximum value of "sf_minreflevel" over all processors