Basilisk CFD
Adaptive Cartesian mesh PDE framework
Loading...
Searching...
No Matches
distance.h File Reference
#include <stdint.h>
#include "PointTriangle.h"
Include dependency graph for distance.h:

Go to the source code of this file.

Data Structures

struct  closest_t
 To compute the minimal distance and its sign, we need to store information on the closest (2 in 2D, 12 in 3D) elements. More...
 

Macros

#define double_to_pointer(d)   (*((void **) &(d)))
 
#define BSIZE   3.
 
#define ND   2
 

Functions

trace coordinput_stl (FILE *fp)
 The 3D triangulated surfaces are defined using the (binary) STL format which can be exported from most CAD modelling programs.
 
trace coordinput_xy (FILE *fp)
 In 2D dimensions, the file format is that used by gnuplot i.e.
 
void bounding_box (coord *p, coord *min, coord *max)
 This function computes the coordinates of the bounding box of a set of segments or triangles.
 
static void update_distance (Point point, coord **i, scalar d)
 
static void refine_distance (Point point, scalar d)
 
static void restriction_distance (Point point, scalar d)
 
static void coarsen_distance (Point point, scalar d)
 
static void delete_distance (scalar d)
 
trace void distance (scalar d, coord *p)
 

Variables

 attribute
 An extra field, holding a pointer to the elements (segments or triangles) intersecting the neighborhood of the cell, is associated with the distance function.
 

Macro Definition Documentation

◆ BSIZE

#define BSIZE   3.

Definition at line 185 of file distance.h.

◆ double_to_pointer

#define double_to_pointer (   d)    (*((void **) &(d)))

Definition at line 184 of file distance.h.

◆ ND

#define ND   2

Definition at line 198 of file distance.h.

Function Documentation

◆ bounding_box()

void bounding_box ( coord p,
coord min,
coord max 
)

This function computes the coordinates of the bounding box of a set of segments or triangles.

Definition at line 159 of file distance.h.

References dimension, HUGE, nodata, p, and x.

◆ coarsen_distance()

static void coarsen_distance ( Point  point,
scalar  d 
)
static

Definition at line 442 of file distance.h.

References d, double_to_pointer, free(), and x.

Referenced by distance().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ delete_distance()

static void delete_distance ( scalar  d)
static

Definition at line 448 of file distance.h.

References d, depth, double_to_pointer, free(), l, and x.

Referenced by distance().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ distance()

trace void distance ( scalar  d,
coord p 
)

Definition at line 459 of file distance.h.

References a, ab, array_append(), array_new(), array_shrink(), boundary_level, coarsen_distance(), d, delete_distance(), depth, dimension, free(), scalar::i, l, n, no_restriction(), nodata, p, point, refine_bilinear(), refine_distance(), restriction_distance(), set_prolongation(), update_distance(), vecdiff, vecdot, vecdotproduct, and x.

Referenced by DistancePointEllipse().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ input_stl()

trace coord * input_stl ( FILE fp)

The 3D triangulated surfaces are defined using the (binary) STL format which can be exported from most CAD modelling programs.

Computation of a distance field from a discretised curve or surface

The goal is to compute the (signed) minimal distance from a set of oriented segments (2d) or triangles (3D) to any point on the grid. This signed distance function is also dynamically updated when the mesh is refined or coarsened.

The scheme is robust and will give results even for inconsistent surface representations (i.e. surfaces with holes, manifold egdes, incompatible orientations etc...). Faces do not need to be properly connected i.e. the scheme works also for "triangle soups".

We need a few utility function such as computation of the minimal distance between a point and a segment or a point and a triangle. The input_stl() function reads such a file a returns an array of triplets of vertex coordinates defining the triangles.

Definition at line 31 of file distance.h.

References a, array_append(), array_new(), array_shrink(), fp, i, j, nodata, p, tag, x, y, and z.

Here is the call graph for this function:

◆ input_xy()

trace coord * input_xy ( FILE fp)

In 2D dimensions, the file format is that used by gnuplot i.e.

pairs of 2D vertex coordinates separated by blank lines. An easy way to create these files is to use a vector graphics editing program (e.g. inkscape) and save the curve as an *.eps* file, then convert it to gnuplot format using

pstoedit -f gnuplot -flat 0.1 file.eps file.gnu

The function below reads such a file and returns an array of pairs of coordinates.

Definition at line 132 of file distance.h.

References a, array_append(), array_new(), array_shrink(), c, fp, nodata, p, and x.

Here is the call graph for this function:

◆ refine_distance()

static void refine_distance ( Point  point,
scalar  d 
)
static

To increase robustness to inconsistent input, we check whether all children are included within the minimum distance sphere. If this is the case then the children and parent must have the same orientation. We enforce this, using the "average" orientation.

Definition at line 406 of file distance.h.

References bilinear(), d, dimension, double_to_pointer, point, s, sign(), update_distance(), and x.

Referenced by distance().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ restriction_distance()

static void restriction_distance ( Point  point,
scalar  d 
)
static

Definition at line 440 of file distance.h.

Referenced by distance().

Here is the caller graph for this function:

◆ update_distance()

static void update_distance ( Point  point,
coord **  i,
scalar  d 
)
static

Definition at line 223 of file distance.h.

References a, ab, array_append(), array_new(), array_shrink(), assert, bilinear(), BSIZE, c, d, dimension, free(), HUGE, i, j, level, n, ND, nn, p, point, PointSegmentDistance(), PointSegmentOrientation(), PointTriangleDistance(), PointTriangleOrientation(), q, s, sign(), sq(), t, type, v, vecdiff, vecdist2, vecdot, vector::x, x, coord::x, vector::y, y, and z.

Referenced by distance(), and refine_distance().

Here is the call graph for this function:
Here is the caller graph for this function:

Variable Documentation

◆ attribute

attribute
Initial value:
{
scalar surface

An extra field, holding a pointer to the elements (segments or triangles) intersecting the neighborhood of the cell, is associated with the distance function.

The neighborhood is a sphere centered on the cell center and with a diameter \(3\Delta\).

Definition at line 180 of file distance.h.