Loading...
Searching...
No Matches
pointSmoother.C
Go to the documentation of this file.
1/*---------------------------------------------------------------------------*\
2 ========= |
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4 \\ / O peration |
5 \\ / A nd | www.openfoam.com
6 \\/ M anipulation |
7-------------------------------------------------------------------------------
8 Copyright (C) 2024-2025 OpenCFD Ltd.
9-------------------------------------------------------------------------------
10License
11 This file is part of OpenFOAM.
12
13 OpenFOAM is free software: you can redistribute it and/or modify it
14 under the terms of the GNU General Public License as published by
15 the Free Software Foundation, either version 3 of the License, or
16 (at your option) any later version.
17
18 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21 for more details.
22
23 You should have received a copy of the GNU General Public License
24 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25
26\*---------------------------------------------------------------------------*/
27
28#include "pointSmoother.H"
31
32// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
33
34namespace Foam
35{
36 defineTypeNameAndDebug(pointSmoother, 0);
38}
39
40
41// * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
42
43bool Foam::pointSmoother::isInternalOrProcessorFace(const label faceI) const
44{
45 if (mesh().isInternalFace(faceI))
46 {
47 return true;
48 }
49
50 if
51 (
52 processorPatchIDs_
53 [
54 mesh().boundaryMesh().patchID()
55 [
56 faceI - mesh().nInternalFaces()
57 ]
58 ]
59 )
60 {
61 return true;
62 }
63
64 return false;
65}
66
67
69(
70 const labelList& facesToMove,
71 const bool moveInternalFaces
72) const
73{
74 bitSet marker(mesh().nPoints());
75
76 for (const label facei : facesToMove)
77 {
78 if (moveInternalFaces || !isInternalOrProcessorFace(facei))
79 {
80 marker.set(mesh().faces()[facei]);
81 }
82 }
83
85 (
86 mesh(),
87 marker,
89 0U
90 );
92 return marker;
93}
94
95
96// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
97
98Foam::pointSmoother::pointSmoother
99(
100 const polyMesh& mesh,
101 const dictionary& dict
102)
103:
104 mesh_(mesh)
105{
106 for (const auto& pp : mesh.boundaryMesh())
107 {
108 if (isA<processorPolyPatch>(pp))
109 {
110 processorPatchIDs_.insert(pp.index());
111 }
113}
114
115
116// * * * * * * * * * * * * * * * * Selectors * * * * * * * * * * * * * * * * //
117
120(
121 const word& pointSmootherType,
122 const polyMesh& mesh,
123 const dictionary& dict
124)
125{
126 Info<< "Selecting pointSmoother type " << pointSmootherType << endl;
127
128 auto* ctorPtr = dictionaryConstructorTable(pointSmootherType);
129
130 if (!ctorPtr)
131 {
133 (
134 dict,
135 typeName,
136 pointSmootherType,
137 *dictionaryConstructorTablePtr_
138 ) << exit(FatalIOError);
140
141 return autoPtr<pointSmoother>(ctorPtr(mesh, dict));
142}
143
144
147(
148 const polyMesh& mesh,
149 const dictionary& dict
150)
151{
152 word pointSmootherType(dict.get<word>(typeName));
154 return New(pointSmootherType, mesh, dict);
155}
156
157
158// * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
159
161(
162 const labelList& facesToMove,
163 const pointField& oldPoints,
164 const pointField& currentPoints,
165 const pointField& faceCentres,
166 const vectorField& faceAreas,
167 const pointField& cellCentres,
168 const scalarField& cellVolumes,
169 pointVectorField& pointDisplacement,
170 const bool correctBCs
171) const
172{
173 // Apply point smoothing (without boundary conditions!)
174 calculate
175 (
176 facesToMove,
177 oldPoints,
178 currentPoints,
179 faceCentres,
180 faceAreas,
181 cellCentres,
182 cellVolumes,
183 pointDisplacement.ref()
184 );
185
186 // Note: do not want to apply boundary conditions whilst smoothing so
187 // disable for now
189 //for (auto& ppf : pointDisplacement.boundaryFieldRef())
190 //{
191 // ppf.operator==(ppf.patchInternalField());
192 //}
193
194 if (correctBCs)
195 {
196 // Apply (multi-)patch constraints
199 pointDisplacement().mesh()
200 ).constrainDisplacement(pointDisplacement);
201 }
202}
203
204
206(
207 const pointField& points,
208 const pointField& faceCentres,
209 const vectorField& faceAreas,
210 const pointField& cellCentres,
211 const scalarField& cellVolumes
212) const
213{
214 // See e.g.
215 // - src/functionObjects/field/stabilityBlendingFactor/
216 // stabilityBlendingFactor.H
217 // - maybe add function to motionSmootherCheck to return value instead
218 // of yes/no for invalid?
219 // - for now just look at non-ortho
220
221 tmp<scalarField> tortho
222 (
224 (
225 mesh(),
226 faceAreas,
227 cellCentres
228 )
229 );
230
231 // tortho.ref().clamp_min(0);
232 return tortho;
233}
234
235
237(
238 const pointField& points,
239 const pointField& faceCentres,
240 const vectorField& faceAreas,
241 const pointField& cellCentres,
242 const scalarField& cellVolumes
243) const
244{
245 // Get face-based non-ortho
246 tmp<scalarField> tfaceOrtho
247 (
248 faceQuality
249 (
250 points,
251 faceCentres,
252 faceAreas,
253 cellCentres,
254 cellVolumes
255 )
256 );
257 const auto& faceOrtho = tfaceOrtho();
258
259 // Min over cells
260 auto tortho = tmp<scalarField>::New(mesh().nCells(), GREAT);
261 auto& ortho = tortho.ref();
262
263 const auto& own = mesh().faceOwner();
264 const auto& nei = mesh().faceNeighbour();
265
266 forAll(own, facei)
267 {
268 auto& o = ortho[own[facei]];
269 o = min(o, faceOrtho[facei]);
270 }
271 for (label facei = 0; facei < mesh().nInternalFaces(); facei++)
272 {
273 auto& o = ortho[nei[facei]];
274 o = min(o, faceOrtho[facei]);
275 }
276
277 return tortho;
278}
279
280
281// ************************************************************************* //
uindirectPrimitivePatch pp(UIndirectList< face >(mesh.faces(), faceLabels), mesh.points())
Internal & ref(const bool updateAccessTime=true)
Same as internalFieldRef().
static FOAM_NO_DANGLING_REFERENCE const pointConstraints & New(const pointMesh &mesh, Args &&... args)
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition autoPtr.H:65
A bitSet stores bits (elements with only two states) in packed internal format and supports a variety...
Definition bitSet.H:61
void set(const bitSet &bitset)
Set specified bits from another bitset.
Definition bitSetI.H:502
Addressing for all faces on surface of mesh. Can either be read from polyMesh or from triSurface....
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
void constrainDisplacement(pointVectorField &displacement, const bool overrideValue=false) const
Apply boundary conditions (single-patch constraints),.
Abstract base class for point smoothing methods. Handles parallel communication via reset and average...
virtual void calculate(const labelList &facesToMove, const pointField &oldPoints, const pointField &currentPoints, const pointField &faceCentres, const vectorField &faceAreas, const pointField &cellCentres, const scalarField &cellVolumes, vectorField &pointDisplacement) const =0
Calculate the point displacement.
bool isInternalOrProcessorFace(const label faceI) const
Test if the given face is internal or on a processor boundary.
static autoPtr< pointSmoother > New(const word &pointSmootherType, const polyMesh &mesh, const dictionary &dict)
Construct given type.
bitSet pointsToMove(const labelList &facesToMove, const bool moveInternalFaces) const
Get list of the points to be moved.
virtual tmp< scalarField > cellQuality(const pointField &points, const pointField &faceCentres, const vectorField &faceAreas, const pointField &cellCentres, const scalarField &cellVolumes) const
Check element quality: 1 = best, 0 = invalid. Topology from mesh, point locations supplied....
const polyMesh & mesh() const noexcept
Access the mesh.
void update(const labelList &facesToMove, const pointField &oldPoints, const pointField &currentPoints, const pointField &faceCentres, const vectorField &faceAreas, const pointField &cellCentres, const scalarField &cellVolumes, pointVectorField &pointDisplacement, const bool correctBCs=true) const
Update the point displacements and apply constraints.
virtual tmp< scalarField > faceQuality(const pointField &points, const pointField &faceCentres, const vectorField &faceAreas, const pointField &cellCentres, const scalarField &cellVolumes) const
Check element quality: 1 = best, 0 = invalid. (also negative?) Topology from mesh,...
static tmp< scalarField > faceOrthogonality(const polyMesh &mesh, const vectorField &fAreas, const vectorField &cellCtrs)
Generate orthogonality field. (1 for fully orthogonal, < 1 for non-orthogonal).
Mesh consisting of general polyhedral cells.
Definition polyMesh.H:79
virtual const labelList & faceOwner() const
Return face owner.
Definition polyMesh.C:1101
virtual const labelList & faceNeighbour() const
Return face neighbour.
Definition polyMesh.C:1107
label nInternalFaces() const noexcept
Number of internal faces.
static void syncPointList(const polyMesh &mesh, List< T > &pointValues, const CombineOp &cop, const T &nullValue, const TransformOp &top)
Synchronize values on all mesh points.
A class for managing temporary objects.
Definition tmp.H:75
static tmp< T > New(Args &&... args)
Construct tmp with forwarding arguments.
Definition tmp.H:215
A class for handling words, derived from Foam::string.
Definition word.H:66
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
U
Definition pEqn.H:72
dynamicFvMesh & mesh
#define FatalIOErrorInLookup(ios, lookupTag, lookupName, lookupTable)
Report an error message using Foam::FatalIOError.
Definition error.H:637
const pointField & points
label nPoints
Namespace for OpenFOAM.
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh > > &tf1, const word &name, const dimensionSet &dimensions, const bool initCopy=false)
Global function forwards to reuseTmpDimensionedField::New.
List< label > labelList
A List of labels.
Definition List.H:62
messageStream Info
Information stream (stdout output on master, null elsewhere).
GeometricField< vector, pointPatchField, pointMesh > pointVectorField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
const word GlobalIOList< Tuple2< scalar, vector > >::typeName("scalarVectorTable")
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:26
Field< vector > vectorField
Specialisation of Field<T> for vector.
IOerror FatalIOError
Error stream (stdout output on all processes), with additional 'FOAM FATAL IO ERROR' header text and ...
vectorField pointField
pointField is a vectorField.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition errorManip.H:125
#define defineRunTimeSelectionTable(baseType, argNames)
Define run-time selection table.
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299