meshRefinementTemplates.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-2015 OpenFOAM Foundation
9 Copyright (C) 2019-2020 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 "meshRefinement.H"
30#include "fvMesh.H"
31#include "globalIndex.H"
32#include "syncTools.H"
33
34// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
35
36template<class T>
38(
39 const labelList& newToOld,
40 const T& nullValue,
41 List<T>& elems
42)
43{
44 List<T> newElems(newToOld.size(), nullValue);
45
46 forAll(newElems, i)
47 {
48 const label oldI = newToOld[i];
49
50 if (oldI >= 0)
51 {
52 newElems[i] = elems[oldI];
53 }
54 }
55
56 elems.transfer(newElems);
57}
58
59
60template<class T>
62(
63 const bitSet& isMasterElem,
64 const UList<T>& values
65)
66{
67 if (values.size() != isMasterElem.size())
68 {
70 << "Number of elements in list " << values.size()
71 << " does not correspond to number of elements in isMasterElem "
72 << isMasterElem.size()
73 << exit(FatalError);
74 }
75
76 T sum = T(Zero);
77 label n = 0;
78
79 forAll(values, i)
80 {
81 if (isMasterElem.test(i))
82 {
83 sum += values[i];
84 n++;
85 }
86 }
87
90
91 if (n > 0)
92 {
93 return sum/n;
94 }
95 else
96 {
97 return pTraits<T>::max;
98 }
99}
100
101
102// Compare two lists over all boundary faces
103template<class T>
105(
106 const scalar tol,
107 const string& msg,
108 const UList<T>& faceData,
109 const UList<T>& syncedFaceData
110) const
111{
112 const label nBFaces = mesh_.nBoundaryFaces();
113
114 if (faceData.size() != nBFaces || syncedFaceData.size() != nBFaces)
115 {
117 << "Boundary faces:" << nBFaces
118 << " faceData:" << faceData.size()
119 << " syncedFaceData:" << syncedFaceData.size()
120 << abort(FatalError);
121 }
122
123 const polyBoundaryMesh& patches = mesh_.boundaryMesh();
124
125 forAll(patches, patchi)
126 {
127 const polyPatch& pp = patches[patchi];
128
129 label bFacei = pp.start() - mesh_.nInternalFaces();
130
131 forAll(pp, i)
132 {
133 const T& data = faceData[bFacei];
134 const T& syncData = syncedFaceData[bFacei];
135
136 if (mag(data - syncData) > tol)
137 {
138 label facei = pp.start()+i;
139
141 << msg
142 << "patchFace:" << i
143 << " face:" << facei
144 << " fc:" << mesh_.faceCentres()[facei]
145 << " patch:" << pp.name()
146 << " faceData:" << data
147 << " syncedFaceData:" << syncData
148 << " diff:" << mag(data - syncData)
149 << abort(FatalError);
150 }
151
152 bFacei++;
153 }
154 }
155}
156
157
158// Print list sorted by coordinates. Used for comparing non-parallel v.s.
159// parallel operation
160template<class T>
162(
163 const UList<point>& points,
164 const UList<T>& data
165)
166{
168
169 pointField allPoints;
170 globalPoints.gather
171 (
172 points,
173 allPoints,
176 );
177
178 List<T> allData;
179 globalPoints.gather
180 (
181 data,
182 allData,
185 );
186
187
188 scalarField magAllPoints(mag(allPoints-point(-0.317, 0.117, 0.501)));
189
190 labelList visitOrder(sortedOrder(magAllPoints));
191 for (const label allPointi : visitOrder)
192 {
193 Info<< allPoints[allPointi] << " : " << allData[allPointi]
194 << endl;
195 }
196}
197
198
199template<class GeoField>
200void Foam::meshRefinement::addPatchFields
201(
202 fvMesh& mesh,
203 const word& patchFieldType
204)
205{
207 (
208 mesh.objectRegistry::lookupClass<GeoField>()
209 );
210
211 forAllIters(flds, iter)
212 {
213 GeoField& fld = *iter();
214 auto& fldBf = fld.boundaryFieldRef();
215
216 label sz = fldBf.size();
217 fldBf.setSize(sz+1);
218 fldBf.set
219 (
220 sz,
221 GeoField::Patch::New
222 (
223 patchFieldType,
224 mesh.boundary()[sz],
225 fld()
226 )
227 );
228 }
229}
230
231
232template<class GeoField>
233void Foam::meshRefinement::reorderPatchFields
234(
235 fvMesh& mesh,
236 const labelList& oldToNew
237)
238{
239 HashTable<GeoField*> flds
240 (
241 mesh.objectRegistry::lookupClass<GeoField>()
242 );
243
244 forAllIters(flds, iter)
245 {
246 iter()->boundaryFieldRef().reorder(oldToNew);
247 }
248}
249
250
251template<class EnumContainer>
253(
254 const EnumContainer& namedEnum,
255 const wordList& words
256)
257{
258 int flags = 0;
259
260 for (const word& w : words)
261 {
262 flags |= namedEnum[w];
263 }
264
265 return flags;
266}
267
268
269template<class Type>
271(
272 const polyMesh& mesh,
273 const bitSet& isMasterEdge,
274 const labelList& meshPoints,
275 const edgeList& edges,
276 const scalarField& edgeWeights,
277 const Field<Type>& pointData,
279)
280{
281 if
282 (
283 edges.size() != isMasterEdge.size()
284 || edges.size() != edgeWeights.size()
285 || meshPoints.size() != pointData.size()
286 )
287 {
289 << "Inconsistent sizes for edge or point data:"
290 << " isMasterEdge:" << isMasterEdge.size()
291 << " edgeWeights:" << edgeWeights.size()
292 << " edges:" << edges.size()
293 << " pointData:" << pointData.size()
294 << " meshPoints:" << meshPoints.size()
295 << abort(FatalError);
296 }
297
298 sum.setSize(meshPoints.size());
299 sum = Type(Zero);
300
301 forAll(edges, edgeI)
302 {
303 if (isMasterEdge.test(edgeI))
304 {
305 const edge& e = edges[edgeI];
306
307 scalar eWeight = edgeWeights[edgeI];
308
309 label v0 = e[0];
310 label v1 = e[1];
311
312 sum[v0] += eWeight*pointData[v1];
313 sum[v1] += eWeight*pointData[v0];
314 }
315 }
316
318 (
319 mesh,
320 meshPoints,
321 sum,
323 Type(Zero) // null value
324 );
325}
326
327
328template<class Type>
330(
331 const dictionary& dict,
332 const word& keyword,
333 const bool noExit,
334 enum keyType::option matchOpt,
335 const Type& deflt
336)
337{
338 Type val(deflt);
339
340 if (!dict.readEntry(keyword, val, matchOpt, !noExit))
341 {
343 << "Entry '" << keyword << "' not found in dictionary "
344 << dict.name() << nl;
345 }
346
347 return val;
348}
349
350
351// ************************************************************************* //
label n
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))
A HashTable similar to std::unordered_map.
Definition: HashTable.H:123
void transfer(List< T > &list)
Definition: List.C:447
unsigned int get() const
Get value as unsigned, no range-checking.
Definition: PackedListI.H:302
label size() const noexcept
Number of entries.
Definition: PackedListI.H:377
A 1D vector of objects of type <T>, where the size of the vector is known and can be used for subscri...
Definition: UList.H:94
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:114
static int & msgType() noexcept
Message tag of standard messages.
Definition: UPstream.H:556
A bitSet stores bits (elements with only two states) in packed internal format and supports a variety...
Definition: bitSet.H:66
bool test(const label pos) const
Test value at specified position, never auto-vivify entries.
Definition: bitSetI.H:521
Database for solution data, solver performance and other reduced data.
Definition: data.H:58
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
An edge is a list of two point labels. The functionality it provides supports the discretisation on a...
Definition: edge.H:66
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:91
const fvBoundaryMesh & boundary() const
Return reference to boundary mesh.
Definition: fvMesh.C:712
Calculates a unique integer (label so might not have enough room - 2G max) for processor + local inde...
Definition: globalIndex.H:68
Calculates points shared by more than two processor patches or cyclic patches.
Definition: globalPoints.H:103
option
Enumeration for the data type and search/match modes (bitmask)
Definition: keyType.H:79
static T gAverage(const bitSet &isMasterElem, const UList< T > &values)
Helper: calculate average.
static void collectAndPrint(const UList< point > &points, const UList< T > &data)
Print list according to (collected and) sorted coordinate.
static int readFlags(const EnumContainer &namedEnum, const wordList &words)
Helper: convert wordList into bit pattern using provided Enum.
static void weightedSum(const polyMesh &mesh, const bitSet &isMasterEdge, const labelList &meshPoints, const edgeList &edges, const scalarField &edgeWeights, const Field< Type > &data, Field< Type > &sum)
Helper: weighted sum (over all subset of mesh points) by.
void testSyncBoundaryFaceList(const scalar mergeDistance, const string &, const UList< T > &, const UList< T > &) const
Compare two lists over all boundary faces.
static void updateList(const labelList &newToOld, const T &nullValue, List< T > &elems)
Helper: reorder list according to map.
A traits class, which is primarily used for primitives.
Definition: pTraits.H:59
const word & name() const noexcept
The patch name.
A polyBoundaryMesh is a polyPatch list with additional search methods and registered IO.
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:81
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:75
label start() const
Return start label of this patch in the polyMesh face list.
Definition: polyPatch.H:364
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 handling words, derived from Foam::string.
Definition: word.H:68
const volScalarField & T
const polyBoundaryMesh & patches
dynamicFvMesh & mesh
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:473
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
const pointField & points
List< label > labelList
A List of labels.
Definition: List.H:66
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &df)
messageStream Info
Information stream (stdout output on master, null elsewhere)
vector point
Point is a vector.
Definition: point.H:43
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
void reduce(const List< UPstream::commsStruct > &comms, T &value, const BinaryOp &bop, const int tag, const label comm)
labelList sortedOrder(const UList< T > &input)
Return the (stable) sort order for the list.
errorManip< error > abort(error &err)
Definition: errorManip.H:144
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
error FatalError
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
constexpr char nl
The newline '\n' character (0x0a)
Definition: Ostream.H:53
dictionary dict
volScalarField & e
Definition: createFields.H:11
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333
#define forAllIters(container, iter)
Iterate across all elements in the container object.
Definition: stdFoam.H:260