Go to the documentation of this file.
39 namespace patchDistMethods
47 Foam::patchDistMethods::exact::patchSurface()
const
49 if (!patchSurfPtr_.valid())
61 treeBoundBox(localBb).extend(
rndGen, 1E-3)
77 dict.add(
"mergeDistance", 1
e-6*localBb.mag());
80 Info<<
"Triangulating local patch faces" <<
nl <<
endl;
86 new distributedTriSurfaceMesh
108 Info<<
"Redistributing surface" <<
nl <<
endl;
109 autoPtr<mapDistribute>
faceMap;
110 autoPtr<mapDistribute> pointMap;
111 patchSurfPtr_().distribute
121 return patchSurfPtr_();
127 Foam::patchDistMethods::exact::exact
138 Foam::patchDistMethods::exact::exact
176 if (info[cellI].hit())
178 const point& cc = mesh_.cellCentres()[cellI];
179 y[cellI] =
mag(cc-info[cellI].hitPoint());
187 y.correctBoundaryConditions();
191 OBJstream str(mesh_.time().timePath()/
"wallPoint.obj");
192 Info<<
type() <<
": dumping nearest wall point to " << str.name()
194 forAll(mesh_.cellCentres(), cellI)
196 const point& cc = mesh_.cellCentres()[cellI];
206 n.correctBoundaryConditions();
int debug
Static debugging option.
List< label > labelList
A List of labels.
virtual const pointField & points() const
Return raw points.
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
static const Enum< distributionType > distributionTypeNames_
Pair< int > faceMap(const label facePi, const face &faceP, const label faceNi, const face &faceN)
OFstream that keeps track of vertices.
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Ostream & endl(Ostream &os)
Add newline and flush stream.
List< treeBoundBox > meshBb(1, treeBoundBox(boundBox(coarseMesh.points(), false)).extend(rndGen, 1e-3))
defineTypeNameAndDebug(advectionDiffusion, 0)
#define forAll(list, i)
Loop across all elements in list.
virtual void findNearest(const pointField &sample, const scalarField &nearestDistSqr, List< pointIndexHit > &) const
messageStream Info
Information stream (uses stdout - output is on the master only)
bool notNull(const T *ptr)
True if ptr is not a pointer (of type T) to the nullObject.
virtual bool correct(volScalarField &y)
Correct the given distance-to-patch field.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Specialisation of patchDist for wall distance calculation.
const fvMesh & mesh_
Reference to the mesh.
Macros for easy insertion into run-time selection tables.
Mesh data needed to do the Finite Volume discretisation.
virtual void getNormal(const List< pointIndexHit > &, vectorField &normal) const
From a set of points and indices get the normal.
const labelHashSet patchIDs_
Set of patch IDs.
addToRunTimeSelectionTable(patchDistMethod, advectionDiffusion, dictionary)
line< point, const point & > linePointRef
Line using referred points.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
IOoject and searching on distributed triSurface. All processor hold (possibly overlapping) part of th...
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
const dimensionedScalar e
Elementary charge.
const Time & time() const
Return the top-level database.
static const GeometricField< vector, fvPatchField, volMesh > & null()
Return a null geometric field.
const word & constant() const
Return constant name.