39Foam::reconstructedDistanceFunction::coupledFacesPatch()
const
41 const polyBoundaryMesh&
patches = mesh_.boundaryMesh();
49 nCoupled +=
pp.size();
59 label facei =
pp.start();
63 nCoupledFaces[nCoupled++] = facei++;
79 const label neiRingLevel
84 if (mesh_.topoChanging())
87 if (nextToInterface_.size() != mesh_.nCells())
89 nextToInterface_.resize(mesh_.nCells());
91 coupledBoundaryPoints_ = coupledFacesPatch()().meshPoints();
97 boolList alreadyMarkedPoint(mesh_.nPoints(),
false);
98 nextToInterface_ =
false;
103 for (label level=0;level<=neiRingLevel;level++)
108 forAll(coupledBoundaryPoints_, i)
110 const label
pi = coupledBoundaryPoints_[i];
113 const label celli = cPoints[
pi][j];
114 if (cellDistLevel_[celli] == level-1)
116 syncMap.insert(
pi,
true);
127 const label
pi = iter.key();
129 if (!alreadyMarkedPoint[
pi])
134 const label pCelli = cPoints[
pi][j];
135 if (cellDistLevel_[pCelli] == -1)
137 cellDistLevel_[pCelli] = level;
138 nextToInterface_[pCelli] =
true;
142 alreadyMarkedPoint[
pi] =
true;
147 forAll(cellDistLevel_, celli)
151 if (interfaceCells[celli])
153 cellDistLevel_[celli] = 0;
154 nextToInterface_[celli] =
true;
158 cellDistLevel_[celli] = -1;
163 if (cellDistLevel_[celli] == level-1)
167 const label pI = pCells[celli][i];
169 if (!alreadyMarkedPoint[pI])
173 const label pCelli = cPoints[pI][j];
174 if (cellDistLevel_[pCelli] == -1)
176 cellDistLevel_[pCelli] = level;
177 nextToInterface_[pCelli] =
true;
181 alreadyMarkedPoint[pI] =
true;
211 coupledBoundaryPoints_(coupledFacesPatch()().meshPoints()),
225 nextToInterface_(
mesh.nCells(), false)
242 if (nextToInterface.
size() != centre.
size())
245 <<
"size of nextToInterface: " << nextToInterface.
size()
246 <<
"size of centre:" << centre.
size()
247 <<
"do not match. Did the mesh change?"
249 return reconDistFunc;
255 Map<vector> mapCentres =
257 Map<vector> mapNormal =
263 forAll(nextToInterface, celli)
265 if (nextToInterface[celli])
267 if (
mag(normal[celli]) != 0)
270 scalar dist = (centre[celli] - mesh_.C()[celli]) &
n;
271 reconDistFunc[celli] = dist;
275 scalar averageDist = 0;
276 scalar avgWeight = 0;
277 const point p = mesh_.C()[celli];
279 for (
const label gblIdx : stencil[celli])
305 vector distanceToIntSeg(c -
p);
306 scalar distToSurf = distanceToIntSeg &
n;
309 if (
mag(distanceToIntSeg) != 0)
311 distanceToIntSeg /=
mag(distanceToIntSeg);
312 weight =
sqr(
mag(distanceToIntSeg &
n));
318 averageDist += distToSurf * weight;
325 reconDistFunc[celli] = averageDist / avgWeight;
331 reconDistFunc[celli] = 0;
343 const label pCellI =
pp.faceCells()[i];
345 if (nextToInterface_[pCellI])
347 scalar averageDist = 0;
348 scalar avgWeight = 0;
351 forAll(stencil[pCellI], j)
353 const label gblIdx = stencil[pCellI][j];
383 vector distanceToIntSeg(c -
p);
384 scalar distToSurf = distanceToIntSeg &
n;
387 if (
mag(distanceToIntSeg) != 0)
389 distanceToIntSeg /=
mag(distanceToIntSeg);
390 weight =
sqr(
mag(distanceToIntSeg &
n));
396 averageDist += distToSurf * weight;
403 pRDF[i] = averageDist / avgWeight;
420 return reconDistFunc;
457 1/acap.patch().deltaCoeffs()*
cos(theta)
458 + RDFbf[patchi].patchInternalField();
constexpr scalar pi(M_PI)
uindirectPrimitivePatch pp(UIndirectList< face >(mesh.faces(), faceLabels), mesh.points())
const Mesh & mesh() const noexcept
Return const reference to mesh.
Boundary & boundaryFieldRef(const bool updateAccessTime=true)
Return a reference to the boundary field.
GeometricBoundaryField< vector, fvsPatchField, surfaceMesh > Boundary
void correctBoundaryConditions()
Correct boundary field.
const Boundary & boundaryField() const noexcept
Return const-reference to the boundary field.
bool insert(const Key &key, const T &obj)
Copy insert a new entry, not overwriting existing entries.
@ NO_READ
Nothing to be read.
@ NO_WRITE
Ignore writing from objectRegistry::writeObject().
@ AUTO_WRITE
Automatically write from objectRegistry::writeObject().
Defines the attributes of an object for which implicit objectRegistry management is supported,...
const Time & time() const noexcept
Return Time associated with the objectRegistry.
A List with indirect addressing.
A HashTable to objects of type <T> with a label key.
void size(const label n)
Older name for setAddressableSize.
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
static autoPtr< T > New(Args &&... args)
Construct autoPtr with forwarding arguments.
A fvBoundaryMesh is a fvPatch list with a reference to the associated fvMesh, with additional search ...
Mesh data needed to do the Finite Volume discretisation.
const fvPatch & patch() const noexcept
Return the patch.
const polyPatch & patch() const noexcept
Return the polyPatch.
A patch is a list of labels that address the faces in the global face list.
void markCellsNearSurf(const boolList &interfaceCells, const label neiRingLevel)
const boolList & nextToInterface() const noexcept
reconstructedDistanceFunction(const fvMesh &mesh)
Construct from fvMesh.
const volScalarField & constructRDF(const boolList &nextToInterface, const volVectorField ¢re, const volVectorField &normal, zoneDistribute &distribute, bool updateStencil=true)
void updateContactAngle(const volScalarField &alpha, const volVectorField &U, surfaceVectorField::Boundary &nHatb)
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
Class for parallel communication in a narrow band. It either provides a Map with the neighbouring val...
List< label > getCyclicPatches(const label celli, const label globalIdx, const vector globalIdxCellCentre) const
Finds and returns list of all cyclic patch labels to which celli's.
void setUpCommforZone(const boolList &zone, bool updateStencil=true)
Update stencil with boolList the size has to match mesh nCells.
Type getValue(const VolumeField< Type > &phi, const Map< Type > &valuesFromOtherProc, const label gblIdx) const
Gives patchNumber and patchFaceNumber for a given Geometric volume field.
const labelListList & getStencil() noexcept
Stencil reference.
Map< Type > getDatafromOtherProc(const boolList &zone, const VolumeField< Type > &phi)
Returns stencil and provides a Map with globalNumbering.
vector getPosition(const VolumeField< vector > &positions, const Map< vector > &valuesFromOtherProc, const label gblIdx, const List< label > cyclicPatchID=List< label >()) const
const polyBoundaryMesh & patches
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
const dimensionedScalar c
Speed of light in a vacuum.
Type & refCast(U &obj)
A dynamic_cast (for references) to Type reference.
GeometricField< vector, fvPatchField, volMesh > volVectorField
const dimensionSet dimless
Dimensionless.
List< labelList > labelListList
List of labelList.
List< label > labelList
A List of labels.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
constexpr scalar degToRad() noexcept
Multiplication factor for degrees to radians conversion.
const Type * isA(const U &obj)
Attempt dynamic_cast to Type.
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
vector point
Point is a vector.
List< bool > boolList
A List of bools.
static constexpr const zero Zero
Global zero (0).
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
errorManipArg< error, int > exit(error &err, const int errNo=1)
fvsPatchField< vector > fvsPatchVectorField
dimensionedScalar cos(const dimensionedScalar &ds)
fvPatchField< scalar > fvPatchScalarField
Us boundaryFieldRef().evaluateCoupled< coupledFaPatch >()
#define forAll(list, i)
Loop across all elements in list.
#define forAllConstIters(container, iter)
Iterate across all elements of the container object with const access.
Unit conversion functions.