diff options
author | jthorn <jthorn@f88db872-0e4f-0410-b76b-b9085cfa78c5> | 2002-10-31 18:30:50 +0000 |
---|---|---|
committer | jthorn <jthorn@f88db872-0e4f-0410-b76b-b9085cfa78c5> | 2002-10-31 18:30:50 +0000 |
commit | 1f22316779d658286289f583363095a31df67290 (patch) | |
tree | c6484f6b516c4a61c6440409e6677e0da58f4f93 /param.ccl | |
parent | fea351e3f95fd7e7f525c018ee989db296e3f7ad (diff) |
*** empty log message ***
git-svn-id: http://svn.einsteintoolkit.org/cactus/EinsteinAnalysis/AHFinderDirect/trunk@874 f88db872-0e4f-0410-b76b-b9085cfa78c5
Diffstat (limited to 'param.ccl')
-rw-r--r-- | param.ccl | 41 |
1 files changed, 27 insertions, 14 deletions
@@ -455,42 +455,55 @@ keyword patch_system_type "what type of patch system should we use?" # ... 1 if FINITE_DIFF_ORDER is set to 2 in "src/include/config.hh" # The code checks for this being too small, and reports a fatal error if so. # -int N_ghost_points "number of ghost zones on each side of a patch" +int ghost_zone_width "number of ghost zones on each side of a patch" { 0:* :: "any integer >= 0" } 2 -int N_overlap_points \ +# +# Our code that computes surface integrals over patches (used for +# computing BH diagnostics like centroids, areas, masses, etc) silently +# assumes that this parameter is == 1, so you should probably leave +# it at that setting. +# +int patch_overlap_width \ "number of grid points that nominally-just-touching patches should overlap" { 1:*:2 :: "any integer >= 0; current implementation requires that it be odd" } 1 # +# This parameter sets the angular resolution of all the patch systems: +# the angular grid spacing in degrees is 90.0/N_zones_per_right_angle. +# # In practice the error in the horizon position is usually dominated # by the errors from interpolating the Cactus gij and Kij to the horizon # position, not by the angular finite differencing or interpatch interpolation -# errors. Thus this parameter can be made quite large (low resolution) +# errors. Thus this parameter can be made fairly small (low resolution) # for better performance, without seriously affecting the accuracy # with which we can locate the horizon. # -# This spacing must integrally divide 90 degrees for full sphere patch -# systems, or 45 degrees for other patch systems. +# For any patch system type other than "full sphere", there are patches +# with 45 degree widths, so this parameter must be even. # -# Also, the computation of surface integrals over the horizon (for horizon -# areas, masses, and centroids) is most accurate if the number of angular -# grid zones (i.e. the patch size divided by this parameter) is even or -# >= 7 (or both) for both coordinates of all patches. +# Normally we use Simpson's rule (in a variant which doesn't require the +# number of zones to be even) for angular integrations over the patch +# system. However, if the number of angular zones in a patch is very +# small and odd (i.e. 3 or 5), then we have to use the trapezoid rule +# instead, so the integrations are less accurate. This occurs for +# N_zones_per_right_angle = 3 or 5 for a full sphere patch system, +# or N_zones_per_right_angle = 6 or 10 for any other patch system type. # -# If you are thinking of setting this to a small value (high resolution), +# If you are thinking of setting this to a large value (high resolution), # note also that with the current dense-matrix storage of the Jacobian, # the memory/running time of the LAPACK linear system solve scales as -# the inverse 4th/6th power of this parameter! +# the 4th/6th power of this parameter! For example, doubling the resolution +# takes 16 times as much memory, and 64 times as long to run! # -real delta_drho_dsigma "angular grid spacing of patches, in degrees" +int N_zones_per_right_angle "sets angular resolution of patch systems" { -(0.0:* :: "any real number > 0.0" -} 7.5 +1:* :: "any integer >= 1; must be even for patch systems other than full-sphere" +} 12 ################################################################################ |