cellDecomposerTemplates.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 "fvMesh.H"
29 #include "emptyFvPatchFields.H"
31 #include "mapPolyMesh.H"
32 //#include "polyPatch.H"
33 //#include "lduSchedule.H"
34 //#include "meshToMesh.H"
35 
36 //template<class Type>
37 //void Foam::functionObjects::cellDecomposer::evaluateConstraintTypes
38 //(
39 // GeometricField<Type, fvPatchField, volMesh>& fld
40 //) const
41 //{
42 // auto& fldBf = fld.boundaryFieldRef();
43 //
44 // const UPstream::commsTypes commsType = UPstream::defaultCommsType;
45 // const label startOfRequests = UPstream::nRequests();
46 //
47 // if
48 // (
49 // commsType == UPstream::commsTypes::blocking
50 // || commsType == UPstream::commsTypes::nonBlocking
51 // )
52 // {
53 // forAll(fldBf, patchi)
54 // {
55 // fvPatchField<Type>& tgtField = fldBf[patchi];
56 //
57 // if
58 // (
59 // tgtField.type() == tgtField.patch().patch().type()
60 // && polyPatch::constraintType(tgtField.patch().patch().type())
61 // )
62 // {
63 // tgtField.initEvaluate(commsType);
64 // }
65 // }
66 //
67 // // Wait for outstanding requests
68 // if (commsType == UPstream::commsTypes::nonBlocking)
69 // {
70 // UPstream::waitRequests(startOfRequests);
71 // }
72 //
73 // forAll(fldBf, patchi)
74 // {
75 // fvPatchField<Type>& tgtField = fldBf[patchi];
76 //
77 // if
78 // (
79 // tgtField.type() == tgtField.patch().patch().type()
80 // && polyPatch::constraintType(tgtField.patch().patch().type())
81 // )
82 // {
83 // tgtField.evaluate(commsType);
84 // }
85 // }
86 // }
87 // else if (commsType == UPstream::commsTypes::scheduled)
88 // {
89 // const lduSchedule& patchSchedule =
90 // fld.mesh().globalData().patchSchedule();
91 //
92 // for (const auto& schedEval : patchSchedule)
93 // {
94 // const label patchi = schedEval.patch;
95 //
96 // fvPatchField<Type>& tgtField = fldBf[patchi];
97 //
98 // if
99 // (
100 // tgtField.type() == tgtField.patch().patch().type()
101 // && polyPatch::constraintType(tgtField.patch().patch().type())
102 // )
103 // {
104 // if (schedEval.init)
105 // {
106 // tgtField.initEvaluate(commsType);
107 // }
108 // else
109 // {
110 // tgtField.evaluate(commsType);
111 // }
112 // }
113 // }
114 // }
115 //}
116 
117 template<class Type>
118 Foam::tmp
119 <
121 >
122 Foam::functionObjects::cellDecomposer::interpolate
123 (
124  const GeometricField<Type, fvPatchField, volMesh>& vf,
125  const fvMesh& sMesh,
126  const labelUList& patchMap,
127  const labelUList& cellMap,
128  const labelUList& faceMap,
129  const bool allowUnmapped
130 ) const
131 {
132  // 1. Create the complete field with dummy patch fields
133  PtrList<fvPatchField<Type>> patchFields(patchMap.size());
134 
135  forAll(patchFields, patchi)
136  {
137  // Set the first one by hand as it corresponds to the
138  // exposed internal faces. Additional interpolation can be put here
139  // as necessary.
140  if (patchMap[patchi] == -1)
141  {
142  patchFields.set
143  (
144  patchi,
145  new emptyFvPatchField<Type>
146  (
147  sMesh.boundary()[patchi],
149  )
150  );
151  }
152  else
153  {
154  patchFields.set
155  (
156  patchi,
158  (
160  sMesh.boundary()[patchi],
162  )
163  );
164  }
165  }
166 
167  auto tresult = tmp<GeometricField<Type, fvPatchField, volMesh>>::New
168  (
169  IOobject
170  (
171  "subset"+vf.name(),
172  sMesh.time().timeName(),
173  sMesh,
176  ),
177  sMesh,
178  vf.dimensions(),
179  Field<Type>(vf.primitiveField(), cellMap),
180  patchFields
181  );
182  auto& result = tresult.ref();
183  result.oriented() = vf.oriented();
184 
185 
186  // 2. Change the fvPatchFields to the correct type using a mapper
187  // constructor (with reference to the now correct internal field)
188 
189  auto& bf = result.boundaryFieldRef();
190 
191  forAll(bf, patchi)
192  {
193  const label basePatchId = patchMap[patchi];
194 
195  if (basePatchId != -1)
196  {
197  // Construct addressing
198  const fvPatch& subPatch = sMesh.boundary()[patchi];
199  const fvPatch& basePatch = vf.mesh().boundary()[basePatchId];
200  const label baseStart = basePatch.start();
201  const label baseSize = basePatch.size();
202 
203  labelList directAddressing(subPatch.size());
204 
205  forAll(directAddressing, i)
206  {
207  const label baseFacei = faceMap[subPatch.start()+i];
208 
209  if (baseFacei >= baseStart && baseFacei < baseStart+baseSize)
210  {
211  directAddressing[i] = baseFacei-baseStart;
212  }
213  else
214  {
215  // Mapped from internal face. Do what? Leave up to
216  // fvPatchField
217  directAddressing[i] = -1;
218  }
219  }
220 
221 
222  directFvPatchFieldMapper mapper(directAddressing);
223 
224  // allowUnmapped : special mode for if we do not want to be
225  // warned for unmapped faces (e.g. from fvMeshDistribute).
226 
227  const bool hasUnmapped = mapper.hasUnmapped();
228  if (allowUnmapped)
229  {
230  mapper.hasUnmapped() = false;
231  }
232 
233  bf.set
234  (
235  patchi,
237  (
238  vf.boundaryField()[basePatchId],
239  subPatch,
240  result.internalField(),
241  mapper
242  )
243  );
244 
245  if (allowUnmapped && hasUnmapped)
246  {
247  // Set unmapped values to zeroGradient. This is the default
248  // action for unmapped fvPatchFields. Note that this bypasses
249  // any special logic for handling unmapped fvPatchFields but
250  // since this is only used inside fvMeshDistribute ...
251 
252  tmp<Field<Type>> tfld(bf[patchi].patchInternalField());
253  const Field<Type>& fld = tfld();
254 
255  Field<Type> value(bf[patchi]);
256  forAll(directAddressing, i)
257  {
258  if (directAddressing[i] == -1)
259  {
260  value[i] = fld[i];
261  }
262  }
263  bf[patchi].fvPatchField<Type>::operator=(value);
264  }
265  }
266  }
267 
268  return tresult;
269 }
270 
271 
272 template<class Type>
273 bool Foam::functionObjects::cellDecomposer::mapFieldType() const
274 {
275  typedef GeometricField<Type, fvPatchField, volMesh> VolFieldType;
276 
277  const fvMesh& mapRegion =
278  this->mesh_.time().lookupObject<fvMesh>(mapRegion_);
279 
280  const labelList patchMap(identity(mapRegion.boundaryMesh().size()));
281 
282  const wordList fieldNames
283  (
284  this->mesh_.sortedNames<VolFieldType>(fieldNames_)
285  );
286 
287  const bool processed = !fieldNames.empty();
288 
289  for (const word& fieldName : fieldNames)
290  {
291  const VolFieldType& field = lookupObject<VolFieldType>(fieldName);
292 
293  auto* mapFieldPtr = mapRegion.getObjectPtr<VolFieldType>(fieldName);
294 
295  if (!mapFieldPtr)
296  {
297  mapFieldPtr = new VolFieldType
298  (
299  IOobject
300  (
301  fieldName,
302  time_.timeName(),
303  mapRegion,
307  ),
308  mapRegion,
309  dimensioned<Type>(field.dimensions(), Zero)
310  );
311 
312  mapFieldPtr->store();
313  }
314 
315  auto& mappedField = *mapFieldPtr;
316 
317  mappedField = interpolate
318  (
319  field,
320  mapRegion,
321  patchMap,
322  mapPtr_().cellMap(),
323  mapPtr_().faceMap(),
324  false //allowUnmapped
325  );
326  Log << " " << fieldName << ": interpolated";
327 
328  //evaluateConstraintTypes(mappedField);
329  }
330 
331  return processed;
332 }
333 
334 
335 template<class Type>
336 bool Foam::functionObjects::cellDecomposer::writeFieldType() const
337 {
338  typedef GeometricField<Type, fvPatchField, volMesh> VolFieldType;
339 
340  const fvMesh& mapRegion =
341  this->mesh_.time().lookupObject<fvMesh>(mapRegion_);
342 
343  const wordList fieldNames
344  (
345  this->mesh_.sortedNames<VolFieldType>(fieldNames_)
346  );
347 
348  const bool processed = !fieldNames.empty();
349 
350  for (const word& fieldName : fieldNames)
351  {
352  const VolFieldType& mappedField =
353  mapRegion.template lookupObject<VolFieldType>(fieldName);
354 
355  mappedField.write();
356 
357  Log << " " << fieldName << ": written";
358  }
359 
360  return processed;
361 }
362 
363 
364 // ************************************************************************* //
rDeltaTY field()
static const DimensionedField< Type, GeoMesh > & null() noexcept
Return a null DimensionedField (reference to a nullObject).
bool interpolate(const vector &p1, const vector &p2, const vector &o, vector &n, scalar l)
Definition: curveTools.C:75
Ignore writing from objectRegistry::writeObject()
DirectFieldMapper< fvPatchFieldMapper > directFvPatchFieldMapper
A fvPatchFieldMapper with direct mapping.
UList< label > labelUList
A UList of labels.
Definition: UList.H:78
Pair< int > faceMap(const label facePi, const face &faceP, const label faceNi, const face &faceN)
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
labelList identity(const label len, label start=0)
Return an identity map of the given length with (map[i] == i), works like std::iota() but returning a...
Definition: labelLists.C:44
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))
List< word > wordList
List of word.
Definition: fileName.H:59
Nothing to be read.
static const word & calculatedType() noexcept
The type name for calculated patch fields.
Definition: fvPatchField.H:204
#define Log
Definition: PDRblock.C:28
static autoPtr< functionObject > New(const word &name, const Time &runTime, const dictionary &dict)
Select from dictionary, based on its "type" entry.
List< label > labelList
A List of labels.
Definition: List.H:62
A class for managing temporary objects.
Definition: HashPtrTable.H:50
Request registration (bool: true)
static tmp< fvPatchField< Type > > New(const word &patchFieldType, const fvPatch &, const DimensionedField< Type, volMesh > &)
Return a pointer to a new patchField created on freestore given.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127