56 const scalar minRatio,
69 for (label facei = 0; facei <
mesh_.
nFaces(); facei++)
74 const vector& fcCorr = faceCorrection[facei];
82 const vector& fn = faceNormals[facei];
83 const point& fc = faceCentres[facei];
140 ownHeight.setSize(mesh_.nFaces());
141 neiHeight.setSize(mesh_.nInternalFaces());
143 const labelList& own = mesh_.faceOwner();
144 const labelList& nei = mesh_.faceNeighbour();
146 for (label facei = 0; facei < mesh_.nInternalFaces(); facei++)
153 ownHeight[facei] = ((fc-ownCc)&
n);
154 neiHeight[facei] = ((neiCc-fc)&
n);
157 for (label facei = mesh_.nInternalFaces(); facei < mesh_.nFaces(); facei++)
163 ownHeight[facei] = ((fc-ownCc)&
n);
183 const labelList& own = mesh_.faceOwner();
184 const labelList& nei = mesh_.faceNeighbour();
187 for (label facei = 0; facei < mesh_.nInternalFaces(); facei++)
193 const vector&
n = faceNormals[facei];
194 const point& fc = faceCentres[facei];
197 const label ownCelli = own[facei];
201 const scalar ownHeight = ((fc-ownCc)&
n);
202 if (ownHeight < minOwnHeight[facei])
214 const label neiCelli = nei[facei];
218 const scalar neiHeight = ((neiCc-fc)&
n);
219 if (neiHeight < minNeiHeight[facei])
232 for (label facei = mesh_.nInternalFaces(); facei < mesh_.nFaces(); facei++)
238 const vector&
n = faceNormals[facei];
239 const point& fc = faceCentres[facei];
242 const label ownCelli = own[facei];
246 const scalar ownHeight = ((fc-ownCc)&
n);
247 if (ownHeight < minOwnHeight[facei])
263Foam::tmp<Foam::pointField>
271 const labelList& own = mesh_.faceOwner();
272 const labelList& nei = mesh_.faceNeighbour();
276 auto& cc = tcc.ref();
281 for (label facei = 0; facei < mesh_.nInternalFaces(); facei++)
286 const vector&
n = faceNormals[facei];
288 const point& ownCc = cellCentres[own[facei]];
289 const point& neiCc = cellCentres[nei[facei]];
303 const scalar w = 0.5*faceWeights[facei];
304 cc[own[facei]] +=
point(w*d);
305 cellWeights[own[facei]] += w;
307 cc[nei[facei]] -=
point(w*d);
308 cellWeights[nei[facei]] += w;
325 const label meshFacei =
pp.start()+i;
326 const label bFacei = meshFacei-mesh_.nInternalFaces();
331 const vector&
n = faceNormals[meshFacei];
334 const point& ownCc = cellCentres[fc[i]];
335 const point& neiCc = neiCellCentres[bFacei];
348 const scalar w = 0.5*faceWeights[meshFacei];
349 cc[fc[i]] +=
point(w*d);
350 cellWeights[fc[i]] += w;
358 if (cellWeights[celli] > VSMALL)
360 cc[celli] = cellCentres[celli] + cc[celli]/cellWeights[celli];
364 cc[celli] = cellCentres[celli];
381 const labelList& own = mesh_.faceOwner();
382 const labelList& nei = mesh_.faceNeighbour();
386 auto& newFc = tnewFc.ref();
389 for (label facei = 0; facei < mesh_.nInternalFaces(); facei++)
395 const vector&
n = faceNormals[facei];
396 const point& oldFc = faceCentres[facei];
416 const solveScalar
f = (deltaFc&
n)/(deltaCc&
n);
428 newFc[facei] = oldFc + d;
446 const label facei =
pp.start()+i;
447 const label bFacei = facei-mesh_.nInternalFaces();
453 const vector&
n = faceNormals[facei];
454 const point& oldFc = faceCentres[facei];
464 const solveScalar
f = (deltaFc&
n)/(deltaCc&
n);
472 newFc[facei] = oldFc + d;
480 const label facei =
pp.start()+i;
486 const vector&
n = faceNormals[facei];
487 const point& oldFc = faceCentres[facei];
497 newFc[facei] = oldFc+d;
523 faceWeights.setSize(cosAngles.size());
530 const scalar cosAngle = cosAngles[facei];
531 if (cosAngle < minCos)
533 faceWeights[facei] = 1.0;
535 else if (cosAngle > maxCos)
537 faceWeights[facei] = 0.0;
542 1.0-(cosAngle-minCos)/(maxCos-minCos);
551Foam::averageNeighbourFvGeometryScheme::averageNeighbourFvGeometryScheme
560 dict.getCheckOrDefault<label>
564 [](label val) {
return val >= 0; }
572 [](scalar val) { return val > 0 && val <= 1; }
581 [](scalar val) { return val >= 0 && val <= 1; }
587 Pout<<
"averageNeighbourFvGeometryScheme :"
588 <<
" nIters:" << nIters_
589 <<
" relax:" << relax_
590 <<
" minRatio:" << minRatio_ <<
endl;
604 Pout<<
"averageNeighbourFvGeometryScheme::movePoints() : "
605 <<
"recalculating primitiveMesh centres" <<
endl;
623 const vectorField faceNormals(faceAreas/magFaceAreas);
630 calcAspectRatioWeights(cellWeight, faceWeight);
633 cellWeight *= relax_;
650 minOwnHeight *= minRatio_;
651 minNeiHeight *= minRatio_;
655 autoPtr<OBJstream> osPtr;
656 autoPtr<surfaceWriter> writerPtr;
663 mesh_.time().timePath()
664 /
"cellCentre_correction.obj"
667 Pout<<
"averageNeighbourFvGeometryScheme::movePoints() :"
668 <<
" writing cell centre path to " << osPtr().
name() <<
endl;
674 mesh_.time().globalPath()
676 / mesh_.pointsInstance()
686 writerPtr->useTimeDir(
true);
688 writerPtr->beginTime(mesh_.time());
694 (outputDir /
"cosAngle"),
698 writerPtr->endTime();
707 for (label iter = 0; iter < nIters_; iter++)
740 writerPtr->beginTime(
instant(scalar(iter)));
741 writerPtr->write(
"cosAngles", cosAngles);
742 writerPtr->endTime();
751 Pout<<
" face:" << facei
752 <<
" at:" << mesh_.faceCentres()[facei]
753 <<
" cosAngle:" << cosAngles[facei]
754 <<
" faceWeight:" << faceWeights[facei]
760 tcc = averageNeighbourCentres
774 const label nClipped = clipPyramids
790 forAll(cellCentres, celli)
792 const point& cc = cellCentres[celli];
818 Pout<<
" iter:" << iter
819 <<
" nClipped:" << nClipped
820 <<
" average displacement:" <<
gAverage(magCorrection)
821 <<
" non-ortho angle : average:" <<
gAverage(nonOrthoAngle)
822 <<
" max:" <<
gMax(nonOrthoAngle) <<
endl;
835 vectorField faceCorrection(faceWeight*(tfc-mesh_.faceCentres()));
845 vectorField faceCentres(mesh_.faceCentres() + faceCorrection);
852 Pout<<
"averageNeighbourFvGeometryScheme::movePoints() :"
853 <<
" averageNeighbour weight"
855 <<
" average:" << avg <<
endl;
858 const fileName tp(mesh_.time().timePath());
860 OBJstream str(tp/
"averageNeighbourCellCentres.obj");
861 Pout<<
"Writing lines from old to new cell centre to " << str.name()
863 forAll(mesh_.cellCentres(), celli)
865 const point& oldCc = mesh_.cellCentres()[celli];
866 const point& newCc = cellCentres[celli];
867 str.writeLine(oldCc, newCc);
873 const fileName tp(mesh_.time().timePath());
874 OBJstream str(tp/
"averageFaceCentres.obj");
875 Pout<<
"Writing lines from old to new face centre to " << str.name()
877 forAll(mesh_.faceCentres(), facei)
879 const point& oldFc = mesh_.faceCentres()[facei];
880 const point& newFc = faceCentres[facei];
881 str.writeLine(oldFc, newFc);
891 std::move(faceCentres),
892 std::move(faceAreas),
893 std::move(cellCentres),
894 std::move(cellVolumes)
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
uindirectPrimitivePatch pp(UIndirectList< face >(mesh.faces(), faceLabels), mesh.points())
const polyBoundaryMesh & pbm
Generic templated field type that is much like a Foam::List except that it is expected to hold numeri...
void setSize(label n)
Alias for resize().
An OFstream that keeps track of vertices and provides convenience output methods for OBJ files.
Ostream & writeLine(const point &p0, const point &p1)
Write line joining two points.
virtual const fileName & name() const override
Read/write access to the name of the stream.
virtual const fileName & name() const override
Get the name of the output serial stream. (eg, the name of the Fstream file name).
void size(const label n)
Older name for setAddressableSize.
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
void reset(T *p=nullptr) noexcept
Delete managed object and set to new given pointer.
Geometry calculation scheme to minimise non-orthogonality/.
const label nIters_
Number of averaging iterations.
label clipFaceTet(const scalar minRatio, const vectorField &faceCentres, const vectorField &faceNormals, vectorField &faceCorrection) const
Clip face-centre correction vector if new triangle area would get below min. Return number of clipped...
const scalar minRatio_
Clipping for pyramid heights - allowable shrinkage as fraction.
tmp< pointField > averageCentres(const pointField &cellCentres, const pointField &faceCentres, const vectorField &faceNormals) const
virtual void movePoints()
Do what is necessary if the mesh has moved.
void makePyrHeights(const pointField &cellCentres, const vectorField &faceCentres, const vectorField &faceNormals, scalarField &ownHeight, scalarField &neiHeight) const
Calculate pyramid heights.
tmp< pointField > averageNeighbourCentres(const pointField &cellCentres, const vectorField &faceNormals, const scalarField &faceWeights) const
Average neighbouring cell centres to minimise non-ortho.
const scalar relax_
Blending between old-iteration cell centres and current average.
label clipPyramids(const pointField &cellCentres, const vectorField &faceCentres, const vectorField &faceNormals, const scalarField &minOwnHeight, const scalarField &minNeiHeight, vectorField &correction) const
Clip correction vector if new pyramid height would get below min. Return number of clipped cells.
void makeNonOrthoWeights(const pointField &cellCentres, const vectorField &faceNormals, scalarField &cosAngles, scalarField &faceWeights) const
Make weights based on non-orthogonality.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
T getCheck(const word &keyword, const Predicate &pred, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a T with additional checking FatalIOError if not found, or if the number of tokens is...
T getCheckOrDefault(const word &keyword, const T &deflt, const Predicate &pred, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a T, or return the given default value. FatalIOError if it is found and the number of...
A face is a list of labels corresponding to mesh vertices.
A class for handling file names.
static bool clean(std::string &str)
Cleanup filename string, possibly applies other transformations such as changing the path separator e...
static word outputPrefix
Directory prefix.
Abstract base class for geometry calculation schemes.
const fvMesh & mesh_
Hold reference to mesh.
const fvMesh & mesh() const
Return mesh reference.
Mesh data needed to do the Finite Volume discretisation.
Geometry calculation scheme with automatic stabilisation for high-aspect ratio cells.
virtual void movePoints()
Do what is necessary if the mesh has moved.
void calcAspectRatioWeights(scalarField &cellWeight, scalarField &faceWeight) const
Calculate cell and face weight. Is 0 for cell < minAspect, 1 for.
An instant of time. Contains the time value and name. Uses Foam::Time when formatting the name.
A polyBoundaryMesh is a polyPatch list with registered IO, a reference to the associated polyMesh,...
virtual const faceList & faces() const
Return raw faces.
virtual const pointField & points() const
Return raw points.
A patch is a list of labels that address the faces in the global face list.
label nFaces() const noexcept
Number of mesh faces.
void resetGeometry(pointField &&faceCentres, pointField &&faceAreas, pointField &&cellCentres, scalarField &&cellVolumes)
Reset the local geometry.
static autoPtr< surfaceWriter > New(const word &writeType)
Select construct a surfaceWriter.
A class for managing temporary objects.
static tmp< T > New(Args &&... args)
Construct tmp with forwarding arguments.
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
Represents 0/1 range or concept. Used for tagged dispatch or clamping.
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Namespace for handling debugging switches.
Type gAverage(const FieldField< Field, Type > &f, const label comm)
The global arithmetic average of a FieldField.
List< label > labelList
A List of labels.
bool mkDir(const fileName &pathName, mode_t mode=0777)
Make a directory and return an error if it could not be created.
dimensionSet clamp(const dimensionSet &a, const dimensionSet &range)
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
T returnReduce(const T &value, BinaryOp bop, const int tag=UPstream::msgType(), const int communicator=UPstream::worldComm)
Perform reduction on a copy, using specified binary operation.
Ostream & endl(Ostream &os)
Add newline and flush stream.
constexpr scalar degToRad() noexcept
Multiplication factor for degrees to radians conversion.
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Field< vector > vectorField
Specialisation of Field<T> for vector.
vector point
Point is a vector.
MinMax< Type > gMinMax(const FieldField< Field, Type > &f)
static constexpr const zero Zero
Global zero (0).
tmp< fvMatrix< Type > > correction(const fvMatrix< Type > &)
Return the correction form of the given matrix by subtracting the matrix multiplied by the current fi...
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
vectorField pointField
pointField is a vectorField.
UList< label > labelUList
A UList of labels.
Vector< solveScalar > solveVector
Type gMax(const FieldField< Field, Type > &f)
dimensionedScalar cos(const dimensionedScalar &ds)
dimensionedScalar acos(const dimensionedScalar &ds)
constexpr scalar radToDeg(const scalar rad) noexcept
Conversion from radians to degrees.
#define forAll(list, i)
Loop across all elements in list.
Unit conversion functions.