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 OpenCFD Ltd.
9 -------------------------------------------------------------------------------
10 License
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"
29 #include "pointConstraints.H"
30 #include "polyMeshTools.H"
31 
32 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
33 
34 namespace Foam
35 {
38 }
39 
40 
41 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
42 
43 bool 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  autoPtr<PackedBoolList> markerPtr
75  (
76  new PackedBoolList(mesh().nPoints(), false)
77  );
78 
79  PackedBoolList& marker(markerPtr());
80 
81  forAll(facesToMove, faceToMoveI)
82  {
83  const label faceI(facesToMove[faceToMoveI]);
84 
85  if (moveInternalFaces || !isInternalOrProcessorFace(faceI))
86  {
87  const face& fPoints(mesh().faces()[faceI]);
88 
89  forAll(fPoints, fPointI)
90  {
91  const label pointI(fPoints[fPointI]);
92 
93  marker[pointI] = true;
94  }
95  }
96  }
97 
99  (
100  mesh(),
101  marker,
102  orEqOp<unsigned int>(),
103  0U
104  );
105 
106  return markerPtr;
107 }
108 
109 
110 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
111 
112 Foam::pointSmoother::pointSmoother
113 (
114  const polyMesh& mesh,
115  const dictionary& dict
116 )
117 :
118  mesh_(mesh)
119 {
120  for (const auto& pp : mesh.boundaryMesh())
121  {
122  if (isA<processorPolyPatch>(pp))
123  {
124  processorPatchIDs_.insert(pp.index());
125  }
126  }
127 }
128 
129 
130 // * * * * * * * * * * * * * * * * Selectors * * * * * * * * * * * * * * * * //
131 
134 (
135  const word& pointSmootherType,
136  const polyMesh& mesh,
137  const dictionary& dict
138 )
139 {
140  Info<< "Selecting pointSmoother type " << pointSmootherType << endl;
141 
142  auto cstrIter = dictionaryConstructorTablePtr_->find(pointSmootherType);
143 
144  if (cstrIter == dictionaryConstructorTablePtr_->end())
145  {
146  FatalErrorIn("pointSmoother::New")
147  << "Unknown " << typeName << " type "
148  << pointSmootherType << endl << endl
149  << "Valid " << typeName << " types are : " << endl
150  << dictionaryConstructorTablePtr_->sortedToc()
151  << exit(FatalError);
152  }
154  return cstrIter()(mesh, dict);
155 }
156 
157 
160 (
161  const polyMesh& mesh,
162  const dictionary& dict
163 )
164 {
165  word pointSmootherType(dict.lookup(typeName));
167  return New(pointSmootherType, mesh, dict);
168 }
169 
170 
171 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
172 
174 {}
175 
176 
177 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
178 
180 (
181  const labelList& facesToMove,
182  const pointField& oldPoints,
183  const pointField& currentPoints,
184  const pointField& faceCentres,
185  const vectorField& faceAreas,
186  const pointField& cellCentres,
187  const scalarField& cellVolumes,
188  pointVectorField& pointDisplacement,
189  const bool correctBCs
190 ) const
191 {
192  // Apply point smoothing (without boundary conditions!)
193  calculate
194  (
195  facesToMove,
196  oldPoints,
197  currentPoints,
198  faceCentres,
199  faceAreas,
200  cellCentres,
201  cellVolumes,
202  pointDisplacement.ref()
203  );
204 
205  // Note: do not want to apply boundary conditions whilst smoothing so
206  // disable for now
208  //for (auto& ppf : pointDisplacement.boundaryFieldRef())
209  //{
210  // ppf.operator==(ppf.patchInternalField());
211  //}
212 
213  if (correctBCs)
214  {
215  // Apply (multi-)patch constraints
217  (
218  pointDisplacement().mesh()
219  ).constrainDisplacement(pointDisplacement);
220  }
221 }
222 
223 
225 (
226  const pointField& points,
227  const pointField& faceCentres,
228  const vectorField& faceAreas,
229  const pointField& cellCentres,
230  const scalarField& cellVolumes
231 ) const
232 {
233  // See e.g.
234  // - src/functionObjects/field/stabilityBlendingFactor/
235  // stabilityBlendingFactor.H
236  // - maybe add function to motionSmootherCheck to return value instead
237  // of yes/no for invalid?
238  // - for now just look at non-ortho
239 
240  tmp<scalarField> tortho
241  (
243  (
244  mesh(),
245  faceAreas,
246  cellCentres
247  )
248  );
250  //return max(tortho, scalar(0.0));
251  return tortho;
252 }
253 
254 
256 (
257  const pointField& points,
258  const pointField& faceCentres,
259  const vectorField& faceAreas,
260  const pointField& cellCentres,
261  const scalarField& cellVolumes
262 ) const
263 {
264  // Get face-based non-ortho
265  tmp<scalarField> tfaceOrtho
266  (
267  faceQuality
268  (
269  points,
270  faceCentres,
271  faceAreas,
272  cellCentres,
273  cellVolumes
274  )
275  );
276  const auto& faceOrtho = tfaceOrtho();
277 
278  // Min over cells
279  tmp<scalarField> tortho(new scalarField(mesh().nCells(), GREAT));
280  auto& ortho = tortho.ref();
281 
282  const auto& own = mesh().faceOwner();
283  const auto& nei = mesh().faceNeighbour();
284 
285  forAll(own, facei)
286  {
287  auto& o = ortho[own[facei]];
288  o = min(o, faceOrtho[facei]);
289  }
290  for (label facei = 0; facei < mesh().nInternalFaces(); facei++)
291  {
292  auto& o = ortho[nei[facei]];
293  o = min(o, faceOrtho[facei]);
294  }
295 
296  return tortho;
297 }
298 
299 
300 // ************************************************************************* //
const polyMesh & mesh() const
ITstream & lookup(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find and return an entry data stream. FatalIOError if not found, or not a stream. ...
Definition: dictionary.C:367
dictionary dict
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:125
error FatalError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL ERROR&#39; header text and sta...
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
Abstract base class for point smoothing methods. Handles parallel communication via reset and average...
Definition: pointSmoother.H:50
virtual const labelList & faceNeighbour() const
Return face neighbour.
Definition: polyMesh.C:1122
static const pointConstraints & New(const pointMesh &mesh, Args &&... args)
Get existing or create MeshObject registered with typeName.
Definition: MeshObject.C:53
T & ref() const
Return non-const reference to the contents of a non-null managed pointer.
Definition: tmpI.H:235
void constrainDisplacement(pointVectorField &displacement, const bool overrideValue=false) const
Apply boundary conditions (single-patch constraints),.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
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.
bool insert(const Key &key)
Insert a new entry, not overwriting existing entries.
Definition: HashSet.H:232
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:38
dynamicFvMesh & mesh
const pointField & points
const polyBoundaryMesh & boundaryMesh() const noexcept
Return boundary mesh.
Definition: polyMesh.H:609
A class for handling words, derived from Foam::string.
Definition: word.H:63
label nPoints
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
static autoPtr< pointSmoother > New(const word &pointSmootherType, const polyMesh &mesh, const dictionary &dict)
Construct given type.
virtual const labelList & faceOwner() const
Return face owner.
Definition: polyMesh.C:1116
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.
label nInternalFaces() const noexcept
Number of internal faces.
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:26
static void syncPointList(const polyMesh &mesh, List< T > &pointValues, const CombineOp &cop, const T &nullValue, const TransformOp &top)
Synchronize values on all mesh points.
virtual ~pointSmoother()
Destructor.
autoPtr< PackedBoolList > pointsToMove(const labelList &facesToMove, const bool moveInternalFaces) const
Get a boolean list of the points to be moved.
Definition: pointSmoother.C:62
defineRunTimeSelectionTable(reactionRateFlameArea, dictionary)
defineTypeNameAndDebug(combustionModel, 0)
static tmp< scalarField > faceOrthogonality(const polyMesh &mesh, const vectorField &fAreas, const vectorField &cellCtrs)
Generate orthogonality field. (1 for fully orthogonal, < 1 for non-orthogonal)
Definition: polyMeshTools.C:30
U
Definition: pEqn.H:72
bool isInternalOrProcessorFace(const label faceI) const
Test if the given face is internal or on a processor boundary.
Definition: pointSmoother.C:36
messageStream Info
Information stream (stdout output on master, null elsewhere)
Internal & ref(const bool updateAccessTime=true)
Same as internalFieldRef()
Field< vector > vectorField
Specialisation of Field<T> for vector.
Pointer management similar to std::unique_ptr, with some additional methods and type checking...
Definition: HashPtrTable.H:48
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:75
List< label > labelList
A List of labels.
Definition: List.H:62
A class for managing temporary objects.
Definition: HashPtrTable.H:50
#define FatalErrorIn(functionName)
Report an error message using Foam::FatalError.
Definition: error.H:600
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...
uindirectPrimitivePatch pp(UIndirectList< face >(mesh.faces(), faceLabels), mesh.points())
Namespace for OpenFOAM.
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, point locations supplied. Move to motionSolver level?
bitSet PackedBoolList