Appearance
Constants
This page documents the named integer constants re-exported by use fossil. They are passed as optional arguments to control the behaviour of three things:
- Signed-distance queries — the
sign_algorithmargument ondistance,compute_distance, andis_point_insidepicks one of theSIGN_*codes. - The AABB acceleration structure — the
aabb_tree_kindargument onload_from_fileandanalyzepicks one ofAABB_TREE_SAH_BVH/AABB_TREE_OCTREE;aabb_refinement_levelsacceptsAABB_AUTO_REFINEMENTto delegate the octree depth to the auto-tune heuristic. - Error reporting — the optional
statusargument on every mutating procedure returns one of theSTATUS_*codes.
Each constant is a small integer(I4P). You never need to know its numeric value — refer to it by name and you stay portable.
Sign algorithms
The sign of a signed-distance query (whether the query point lies inside or outside the closed surface) is decided by one of three algorithms. The default is SIGN_PSEUDO_NORMAL; the other two are kept for benchmarking and as fallbacks for meshes whose orientation cannot be sanitised.
SIGN_PSEUDO_NORMAL (default)
Purpose. Bærentzen–Aanæs angle-weighted pseudo-normal test. At the closest point on the surface, the sign of dot(point - closest, pseudo_normal) classifies the query point as inside or outside. The pseudo-normal is fused into the existing distance traversal so signed distance costs the same as unsigned distance.
When to use. Almost always. This is the default for a reason: it is the fastest of the three, robust on closed meshes, and produces a continuous sign field (no flickering on flat regions, no ray-grazing artefacts).
Precondition. Normals must be consistently outward-oriented. sanitize produces this state; sanitize_normals alone is enough if you already trust the rest of the mesh.
Example.
fortran
program ex_sign_pseudo_normal
use fossil
use penf, only : I4P, R8P
use vecfor_R8P, only : vector_R8P, ex_R8P
implicit none
type(surface_stl_object) :: surface
real(R8P) :: d
call surface%load_from_file(file_name='cube.stl', guess_format=.true.)
call surface%sanitize
d = surface%distance(point=2._R8P * ex_R8P, is_signed=.true., &
sign_algorithm=SIGN_PSEUDO_NORMAL, is_square_root=.true.)
print '(A, ES12.5)', 'signed distance (pseudo-normal) = ', d
end program ex_sign_pseudo_normalSee also. distance, is_point_inside, sanitize_normals.
Reference. Bærentzen & Aanæs, Signed Distance Computation Using the Angle Weighted Pseudo-normal, IEEE TVCG 11(3), 2005.
SIGN_RAY_INTERSECTIONS
Purpose. Axis-aligned ray cast from the query point: count surface intersections along the +x axis; odd means inside, even means outside. The sign is decided after the unsigned distance has been computed (two passes total).
When to use. Meshes with mixed normal orientation that you cannot or do not want to repair (open shells, scan fragments). Also useful as a regression oracle against SIGN_PSEUDO_NORMAL on hand-crafted test surfaces.
Cost. ~3× slower than SIGN_PSEUDO_NORMAL because of the second pass and the per-facet ray-triangle test. Ray-grazing tie-breaks are handled but introduce a small failure rate near edge-on geometry.
Example.
fortran
program ex_sign_ray
use fossil
use penf, only : I4P, R8P
use vecfor_R8P, only : vector_R8P, ex_R8P
implicit none
type(surface_stl_object) :: surface
real(R8P) :: d
call surface%load_from_file(file_name='cube.stl', guess_format=.true.)
d = surface%distance(point=2._R8P * ex_R8P, is_signed=.true., &
sign_algorithm=SIGN_RAY_INTERSECTIONS, is_square_root=.true.)
print '(A, ES12.5)', 'signed distance (ray) = ', d
end program ex_sign_raySee also. SIGN_PSEUDO_NORMAL, SIGN_SOLID_ANGLE.
SIGN_SOLID_ANGLE
Purpose. Sum the projected solid angles subtended by each facet from the query point; the sum approaches ±4π inside and 0 outside.
When to use. Reference / regression. Not accelerated — every query is O(N) in the facet count regardless of the AABB tree state. Use on small surfaces or low-frequency queries only.
Example.
fortran
program ex_sign_solid_angle
use fossil
use penf, only : I4P, R8P
use vecfor_R8P, only : vector_R8P, ex_R8P
implicit none
type(surface_stl_object) :: surface
real(R8P) :: d
call surface%load_from_file(file_name='cube.stl', guess_format=.true.)
d = surface%distance(point=2._R8P * ex_R8P, is_signed=.true., &
sign_algorithm=SIGN_SOLID_ANGLE, is_square_root=.true.)
print '(A, ES12.5)', 'signed distance (solid angle) = ', d
end program ex_sign_solid_angleSee also. SIGN_PSEUDO_NORMAL, facet%solid_angle.
sign_algorithm_from_string
Helper for CLI / config parsing. Maps a lower-case name to the matching SIGN_* constant, or raises error stop on an unknown string.
fortran
algo = sign_algorithm_from_string('pseudo_normal') ! == SIGN_PSEUDO_NORMAL
algo = sign_algorithm_from_string('ray_intersections') ! == SIGN_RAY_INTERSECTIONS
algo = sign_algorithm_from_string('solid_angle') ! == SIGN_SOLID_ANGLEAABB tree kind
Each surface_stl_object owns one acceleration structure for distance and inside-test queries. The kind is chosen at load / analyse time and cannot be changed without rebuilding (which is what set_tree_kind followed by analyze does for you).
AABB_TREE_SAH_BVH (default)
Purpose. Binary bounding-volume hierarchy built top-down by partitioning triangles along the longest axis of their centroid bounding box, using the bucketed surface-area heuristic (16 buckets per axis). The tree adapts to local triangle density: leaves cluster tightly where the surface is dense and spread thinly where it is not.
When to use. Almost always. On dragon-fine.stl (24 k facets, 32³ query grid) it is ~110× faster than the octree. Build cost is O(N log N).
Note. Ignores the aabb_refinement_levels argument — depth is a function of the data, not of a user knob.
See also. AABB_TREE_OCTREE, AABB_AUTO_REFINEMENT, aabb_tree_object.
AABB_TREE_OCTREE
Purpose. Eight-way space-partitioning octree. Each non-leaf node is split into eight axis-aligned children regardless of triangle distribution. Depth is fixed up front via aabb_refinement_levels (or AABB_AUTO_REFINEMENT).
When to use. Benchmarking the SAH BVH; rare meshes where the BVH's top-down split happens to perform poorly; teaching / debugging because the geometry is easier to visualise. Not the default — the BVH wins on every realistic input we have measured.
Example. (switching kind on a freshly-loaded surface)
fortran
program ex_octree
use fossil
use penf, only : I4P, R8P
implicit none
type(surface_stl_object) :: surface
call surface%load_from_file(file_name='cube.stl', guess_format=.true., &
aabb_tree_kind=AABB_TREE_OCTREE, &
aabb_refinement_levels=AABB_AUTO_REFINEMENT)
print '(A, I0)', 'octree refinement levels = ', surface%aabb%get_refinement_levels()
end program ex_octreeSee also. AABB_TREE_SAH_BVH, aabb%set_tree_kind.
AABB_AUTO_REFINEMENT
Purpose. Sentinel passed in place of an integer depth: ask the octree to pick a depth from the facet count, ceil(log8(N/64)), clamped to [1, 6].
When to use. Whenever you choose AABB_TREE_OCTREE and have no specific depth in mind. Has no effect when AABB_TREE_SAH_BVH is active (the BVH ignores it).
Status codes
Every mutating procedure on surface_stl_object and a few on facet_object accepts an optional intent(out) status argument. When you supply it, the procedure never calls error stop — failures are reported through the status code instead. When you omit it, the procedure aborts on unrecoverable errors.
| Constant | Meaning |
|---|---|
STATUS_OK | Success. Numeric value 0. |
STATUS_ALLOC_FAIL | An internal allocate failed (usually out of memory). |
STATUS_AMBIGUOUS_ARGS | Conflicting optionals — e.g. both delta= and x=/y=/z= passed to translate. |
STATUS_FILE_NOT_FOUND | load_from_file could not open the input path for reading. |
STATUS_FILE_OPEN_FAIL | save_into_file could not open the output path for writing. |
STATUS_INVALID_INPUT | The input file contains NaN / Inf vertex coordinates and was refused. The mesh is left untouched. |
Idiom for fail-soft loading.
fortran
program ex_status
use fossil
use penf, only : I4P, R8P
implicit none
type(surface_stl_object) :: surface
integer(I4P) :: stat
call surface%load_from_file(file_name='maybe-missing.stl', guess_format=.true., status=stat)
select case (stat)
case (STATUS_OK)
print '(A,I0,A)', 'loaded ', surface%get_facets_number(), ' facets'
case (STATUS_FILE_NOT_FOUND)
print '(A)', 'file not found'
case (STATUS_INVALID_INPUT)
print '(A)', 'file contains NaN or Inf coords; load refused'
case default
print '(A,I0)', 'load failed with status ', stat
end select
end program ex_statusSee also. load_from_file, save_into_file, translate, resize.