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 -------------------------------------------------------------------------------
11 License
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 
36 template<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 
60 template<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 
88  reduce(sum, sumOp<T>());
89  reduce(n, sumOp<label>());
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
103 template<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
160 template<class T>
162 (
163  const UList<point>& points,
164  const UList<T>& data
165 )
166 {
167  const globalIndex globalPoints(points.size());
168 
170  globalPoints.gather
171  (
172  points,
173  allPoints,
174  UPstream::msgType(),
175  Pstream::commsTypes::blocking
176  );
177 
178  List<T> allData;
179  globalPoints.gather
180  (
181  data,
182  allData,
183  UPstream::msgType(),
184  Pstream::commsTypes::blocking
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 
199 template<class GeoField>
200 void 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,
222  (
223  patchFieldType,
224  mesh.boundary()[sz],
225  fld()
226  )
227  );
228  }
229 }
230 
231 
232 template<class GeoField>
233 void 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 
251 template<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 
269 template<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,
322  plusEqOp<Type>(),
323  Type(Zero) // null value
324  );
325 }
326 
327 
328 template<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 // ************************************************************************* //
Foam::labelList
List< label > labelList
A List of labels.
Definition: List.H:67
Foam::meshRefinement::get
static Type get(const dictionary &dict, const word &keyword, const bool noExit, enum keyType::option matchOpt=keyType::REGEX, const Type &deflt=Zero)
Wrapper around dictionary::get which does not exit.
Definition: meshRefinementTemplates.C:330
Foam::meshRefinement::weightedSum
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.
Definition: meshRefinementTemplates.C:271
Foam::meshRefinement::gAverage
static T gAverage(const bitSet &isMasterElem, const UList< T > &values)
Helper: calculate average.
Definition: meshRefinementTemplates.C:62
Foam::globalPoints
Calculates points shared by more than two processor patches or cyclic patches.
Definition: globalPoints.H:102
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:65
Foam::bitSet
A bitSet stores bits (elements with only two states) in packed internal format and supports a variety...
Definition: bitSet.H:63
Foam::polyBoundaryMesh
A polyBoundaryMesh is a polyPatch list with additional search methods and registered IO.
Definition: polyBoundaryMesh.H:63
Foam::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
Foam::HashTableOps::values
List< T > values(const HashTable< T, Key, Hash > &tbl, const bool doSort=false)
List of values from HashTable, optionally sorted.
Definition: HashOps.H:149
globalIndex.H
Foam::edge
An edge is a list of two point labels. The functionality it provides supports the discretisation on a...
Definition: edge.H:63
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:369
Foam::bitSet::test
bool test(const label pos) const
Test value at specified position, never auto-vivify entries.
Definition: bitSetI.H:520
Foam::dictionary::name
const fileName & name() const noexcept
The dictionary name.
Definition: dictionaryI.H:48
Foam::meshRefinement::mesh
const fvMesh & mesh() const
Reference to mesh.
Definition: meshRefinement.H:944
syncTools.H
Foam::polyMesh
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:77
Foam::meshRefinement::updateList
static void updateList(const labelList &newToOld, const T &nullValue, List< T > &elems)
Helper: reorder list according to map.
Definition: meshRefinementTemplates.C:38
Foam::syncTools::syncPointList
static void syncPointList(const polyMesh &mesh, List< T > &pointValues, const CombineOp &cop, const T &nullValue, const TransformOp &top)
Synchronize values on all mesh points.
Definition: syncToolsTemplates.C:721
Foam::sumOp
Definition: ops.H:213
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
n
label n
Definition: TABSMDCalcMethod2.H:31
Foam::Field< vector >
Foam::meshRefinement::readFlags
static int readFlags(const EnumContainer &namedEnum, const wordList &words)
Helper: convert wordList into bit pattern using provided Enum.
Definition: meshRefinementTemplates.C:253
Foam::Info
messageStream Info
Information stream (stdout output on master, null elsewhere)
Foam::polyPatch
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:68
Foam::List::transfer
void transfer(List< T > &list)
Definition: List.C:456
Foam::dictionary::readEntry
bool readEntry(const word &keyword, T &val, enum keyType::option matchOpt=keyType::REGEX, bool mandatory=true) const
Definition: dictionaryTemplates.C:302
fld
gmvFile<< "tracers "<< particles.size()<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().x()<< ' ';}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().y()<< ' ';}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().z()<< ' ';}gmvFile<< nl;for(const word &name :lagrangianScalarNames){ IOField< scalar > fld(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
Definition: gmvOutputLagrangian.H:23
Foam::meshRefinement::collectAndPrint
static void collectAndPrint(const UList< point > &points, const UList< T > &data)
Print list according to (collected and) sorted coordinate.
Definition: meshRefinementTemplates.C:162
forAllIters
#define forAllIters(container, iter)
Iterate across all elements in the container object.
Definition: stdFoam.H:223
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::FatalError
error FatalError
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:123
reduce
reduce(hasMovingMesh, orOp< bool >())
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
meshRefinement.H
Foam::meshRefinement::testSyncBoundaryFaceList
void testSyncBoundaryFaceList(const scalar mergeDistance, const string &, const UList< T > &, const UList< T > &) const
Compare two lists over all boundary faces.
Definition: meshRefinementTemplates.C:105
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:85
fvMesh.H
Foam::abort
errorManip< error > abort(error &err)
Definition: errorManip.H:144
T
const volScalarField & T
Definition: createFieldRefs.H:2
Foam::globalIndex
Calculates a unique integer (label so might not have enough room - 2G max) for processor + local inde...
Definition: globalIndex.H:68
Foam::polyPatch::start
label start() const
Return start label of this patch in the polyMesh face list.
Definition: polyPatch.H:361
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
Foam::HashTable
A HashTable similar to std::unordered_map.
Definition: HashTable.H:105
Foam::New
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
Global function forwards to reuseTmpDimensionedField::New.
Definition: DimensionedFieldReuseFunctions.H:105
Foam::fvMesh::boundary
const fvBoundaryMesh & boundary() const
Return reference to boundary mesh.
Definition: fvMesh.C:685
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
Foam::nl
constexpr char nl
Definition: Ostream.H:404
Foam::List< label >
Foam::pTraits
A traits class, which is primarily used for primitives.
Definition: pTraits.H:56
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Foam::PackedList::size
label size() const noexcept
Number of entries.
Definition: PackedListI.H:377
Foam::UList
A 1D vector of objects of type <T>, where the size of the vector is known and can be used for subscri...
Definition: HashTable.H:103
points
const pointField & points
Definition: gmvOutputHeader.H:1
Foam::constant::electromagnetic::e
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
Foam::sum
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &df)
Definition: DimensionedFieldFunctions.C:327
patches
const polyBoundaryMesh & patches
Definition: convertProcessorPatches.H:65
FatalIOErrorInFunction
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:473
Foam::DelaunayMeshTools::allPoints
tmp< pointField > allPoints(const Triangulation &t)
Extract all points in vertex-index order.
Foam::patchIdentifier::name
const word & name() const noexcept
The patch name.
Definition: patchIdentifier.H:135
Foam::plusEqOp
Definition: ops.H:72
Foam::sortedOrder
labelList sortedOrder(const UList< T > &input)
Return the (stable) sort order for the list.
Foam::UList::size
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:114
Foam::point
vector point
Point is a vector.
Definition: point.H:43
Foam::data
Database for solution data, solver performance and other reduced data.
Definition: data.H:55
Foam::keyType::option
option
Enumeration for the data type and search/match modes (bitmask)
Definition: keyType.H:78