motionSmootherAlgoTemplates.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) 2011-2017 OpenFOAM Foundation
9 Copyright (C) 2016-2018 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
12 This file is part of OpenFOAM.
13
14 OpenFOAM is free software: you can redistribute it and/or modify it
15 under the terms of the GNU General Public License as published by
16 the Free Software Foundation, either version 3 of the License, or
17 (at your option) any later version.
18
19 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22 for more details.
23
24 You should have received a copy of the GNU General Public License
25 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26
27\*---------------------------------------------------------------------------*/
28
29#include "motionSmootherAlgo.H"
30#include "meshTools.H"
32#include "pointConstraint.H"
33#include "pointConstraints.H"
34#include "syncTools.H"
35
36// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
37
38template<class Type>
39void Foam::motionSmootherAlgo::checkConstraints
40(
41 GeometricField<Type, pointPatchField, pointMesh>& pf
42)
43{
44 typedef GeometricField<Type, pointPatchField, pointMesh> FldType;
45
46 const polyMesh& mesh = pf.mesh();
47
48 const polyBoundaryMesh& bm = mesh.boundaryMesh();
49
50 // First count the total number of patch-patch points
51
52 label nPatchPatchPoints = 0;
53
54 for (const polyPatch& pp : bm)
55 {
56 if (!isA<emptyPolyPatch>(pp))
57 {
58 nPatchPatchPoints += pp.boundaryPoints().size();
59 }
60 }
61
62
63 typename FldType::Boundary& bFld = pf.boundaryFieldRef();
64
65
66 // Evaluate in reverse order
67
68 forAllReverse(bFld, patchi)
69 {
70 bFld[patchi].initEvaluate(Pstream::commsTypes::blocking); // buffered
71 }
72
73 forAllReverse(bFld, patchi)
74 {
75 bFld[patchi].evaluate(Pstream::commsTypes::blocking);
76 }
77
78
79 // Save the values
80
81 Field<Type> boundaryPointValues(nPatchPatchPoints);
82 nPatchPatchPoints = 0;
83
84 forAll(bm, patchi)
85 {
86 if (!isA<emptyPolyPatch>(bm[patchi]))
87 {
88 const labelList& bp = bm[patchi].boundaryPoints();
89 const labelList& meshPoints = bm[patchi].meshPoints();
90
91 forAll(bp, pointi)
92 {
93 label ppp = meshPoints[bp[pointi]];
94 boundaryPointValues[nPatchPatchPoints++] = pf[ppp];
95 }
96 }
97 }
98
99
100 // Forward evaluation
101
102 bFld.evaluate();
103
104
105 // Check
106
107 nPatchPatchPoints = 0;
108
109 forAll(bm, patchi)
110 {
111 if (!isA<emptyPolyPatch>(bm[patchi]))
112 {
113 const labelList& bp = bm[patchi].boundaryPoints();
114 const labelList& meshPoints = bm[patchi].meshPoints();
115
116 for (const label pointi : bp)
117 {
118 const label ppp = meshPoints[pointi];
119
120 const Type& savedVal = boundaryPointValues[nPatchPatchPoints++];
121
122 if (savedVal != pf[ppp])
123 {
125 << "Patch fields are not consistent on mesh point "
126 << ppp << " coordinate " << mesh.points()[ppp]
127 << " at patch " << bm[patchi].name() << '.'
128 << endl
129 << "Reverse evaluation gives value " << savedVal
130 << " , forward evaluation gives value " << pf[ppp]
131 << abort(FatalError);
132 }
133 }
134 }
135 }
136}
137
138
139template<class Type>
141Foam::motionSmootherAlgo::avg
142(
144 const scalarField& edgeWeight
145) const
146{
148 (
150 (
151 "avg("+fld.name()+')',
152 fld.time().timeName(),
153 fld.db(),
156 false
157 ),
158 fld.mesh(),
159 dimensioned<Type>(fld.dimensions(), Zero)
160 );
161 auto& res = tres.ref();
162
163 const polyMesh& mesh = fld.mesh()();
164
165
166 // Sum local weighted values and weights
167 // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
168
169 // Note: on coupled edges use only one edge (through isMasterEdge)
170 // This is done so coupled edges do not get counted double.
171
172 scalarField sumWeight(mesh.nPoints(), Zero);
173
174 const edgeList& edges = mesh.edges();
175
176 for (const label edgei : isMasterEdge_)
177 {
178 const edge& e = edges[edgei];
179 const scalar w = edgeWeight[edgei];
180
181 res[e[0]] += w*fld[e[1]];
182 sumWeight[e[0]] += w;
183
184 res[e[1]] += w*fld[e[0]];
185 sumWeight[e[1]] += w;
186 }
187
188
189 // Add coupled contributions
190 // ~~~~~~~~~~~~~~~~~~~~~~~~~
191
193 (
194 mesh,
195 res,
197 Type(Zero) // null value
198 );
200 (
201 mesh,
202 sumWeight,
204 scalar(0) // null value
205 );
206
207
208 // Average
209 // ~~~~~~~
210
211 forAll(res, pointi)
212 {
213 if (mag(sumWeight[pointi]) < VSMALL)
214 {
215 // Unconnected point. Take over original value
216 res[pointi] = fld[pointi];
217 }
218 else
219 {
220 res[pointi] /= sumWeight[pointi];
221 }
222 }
223
224 // Single and multi-patch constraints
225 pointConstraints::New(fld.mesh()).constrain(res, false);
226
227 return tres;
228}
229
230
231template<class Type>
233(
235 const scalarField& edgeWeight,
237) const
238{
239 tmp<pointVectorField> tavgFld = avg(fld, edgeWeight);
240 const pointVectorField& avgFld = tavgFld();
241
242 forAll(fld, pointi)
243 {
244 if (isInternalPoint_.test(pointi))
245 {
246 newFld[pointi] = 0.5*fld[pointi] + 0.5*avgFld[pointi];
247 }
248 }
249
250 // Single and multi-patch constraints
251 pointConstraints::New(fld.mesh()).constrain(newFld, false);
252}
253
254
255template<class Type, class CombineOp>
256void Foam::motionSmootherAlgo::testSyncField
257(
258 const Field<Type>& fld,
259 const CombineOp& cop,
260 const Type& zero,
261 const scalar maxMag
262) const
263{
264 if (debug)
265 {
266 Pout<< "testSyncField : testing synchronisation of Field<Type>."
267 << endl;
268 }
269
270 Field<Type> syncedFld(fld);
271
273 (
274 mesh_,
275 syncedFld,
276 cop,
277 zero
278 );
279
280 forAll(syncedFld, i)
281 {
282 if (mag(syncedFld[i] - fld[i]) > maxMag)
283 {
285 << "On element " << i << " value:" << fld[i]
286 << " synchronised value:" << syncedFld[i]
287 << abort(FatalError);
288 }
289 }
290}
291
292
293template<class Type>
295(
296 const dictionary& dict,
297 const word& keyword,
298 const bool noExit,
299 enum keyType::option matchOpt,
300 const Type& defaultValue
301)
302{
303 Type val(defaultValue);
304
305 if
306 (
308 (
309 keyword,
310 val,
311 matchOpt,
312 !noExit
313 )
314 )
315 {
317 << "Entry '" << keyword << "' not found in dictionary "
318 << dict.name() << endl;
319 }
320 return val;
321}
322
323
324// ************************************************************************* //
Info<< nl<< "Wrote faMesh in vtk format: "<< writer.output().name()<< nl;}{ vtk::lineWriter writer(aMesh.points(), aMesh.edges(), fileName(aMesh.mesh().time().globalPath()/"finiteArea-edges"));writer.writeGeometry();writer.beginCellData(4);writer.writeProcIDs();{ Field< scalar > fld(faMeshTools::flattenEdgeField(aMesh.magLe(), true))
Generic GeometricField class.
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:170
unsigned int get() const
Get value as unsigned, no range-checking.
Definition: PackedListI.H:302
static autoPtr< Time > New()
Construct (dummy) Time - no functionObjects or libraries.
Definition: Time.C:717
label size() const noexcept
The number of elements in the list.
Definition: UPtrListI.H:106
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
const fileName & name() const noexcept
The dictionary name.
Definition: dictionaryI.H:48
bool readEntry(const word &keyword, T &val, enum keyType::option matchOpt=keyType::REGEX, bool mandatory=true) const
Generic dimensioned Type class.
An edge is a list of two point labels. The functionality it provides supports the discretisation on a...
Definition: edge.H:66
option
Enumeration for the data type and search/match modes (bitmask)
Definition: keyType.H:79
const polyMesh & mesh() const
Reference to mesh.
void smooth(const GeometricField< Type, pointPatchField, pointMesh > &fld, const scalarField &edgeWeight, GeometricField< Type, pointPatchField, pointMesh > &newFld) const
Fully explicit smoothing of fields (not positions)
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:81
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:456
virtual const pointField & points() const
Return raw points.
Definition: polyMesh.C:1083
const edgeList & edges() const
Return mesh edges. Uses calcEdges.
label nPoints() const noexcept
Number of mesh points.
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:65
A class for handling words, derived from Foam::string.
Definition: word.H:68
A class representing the concept of 0 (zero) that can be used to avoid manipulating objects known to ...
Definition: zero.H:63
dynamicFvMesh & mesh
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
List< label > labelList
A List of labels.
Definition: List.H:66
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
errorManip< error > abort(error &err)
Definition: errorManip.H:144
IOerror FatalIOError
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh > > &tdf1, const word &name, const dimensionSet &dimensions)
Global function forwards to reuseTmpDimensionedField::New.
error FatalError
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
dictionary dict
volScalarField & e
Definition: createFields.H:11
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333
#define forAllReverse(list, i)
Reverse loop across all elements in list.
Definition: stdFoam.H:346