Loading...
Searching...
No Matches
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-------------------------------------------------------------------------------
10License
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// fld.boundaryFieldRef().evaluate_if
43// (
44// [](const auto& pfld) -> bool
45// {
46// return
47// (
48// pfld.type() == pfld.patch().patch().type()
49// && polyPatch::constraintType(pfld.patch().patch().type())
50// );
51// }
52// );
53//}
54
55template<class Type>
57<
59>
60Foam::functionObjects::cellDecomposer::interpolate
61(
63 const fvMesh& sMesh,
64 const labelUList& patchMap,
65 const labelUList& cellMap,
66 const labelUList& faceMap,
67 const bool allowUnmapped
68) const
69{
70 // 1. Create the complete field with dummy patch fields
71 PtrList<fvPatchField<Type>> patchFields(patchMap.size());
72
73 forAll(patchFields, patchi)
74 {
75 // Set the first one by hand as it corresponds to the
76 // exposed internal faces. Additional interpolation can be put here
77 // as necessary.
78 if (patchMap[patchi] == -1)
79 {
80 patchFields.set
81 (
82 patchi,
84 (
85 sMesh.boundary()[patchi],
87 )
88 );
89 }
90 else
91 {
92 patchFields.set
93 (
94 patchi,
96 (
98 sMesh.boundary()[patchi],
100 )
101 );
102 }
103 }
104
106 (
108 (
109 "subset"+vf.name(),
110 sMesh.time().timeName(),
111 sMesh,
114 ),
115 sMesh,
116 vf.dimensions(),
117 Field<Type>(vf.primitiveField(), cellMap),
118 patchFields
119 );
120 auto& result = tresult.ref();
121 result.oriented() = vf.oriented();
122
123
124 // 2. Change the fvPatchFields to the correct type using a mapper
125 // constructor (with reference to the now correct internal field)
126
127 auto& bf = result.boundaryFieldRef();
128
129 forAll(bf, patchi)
130 {
131 const label basePatchId = patchMap[patchi];
132
133 if (basePatchId != -1)
134 {
135 // Construct addressing
136 const fvPatch& subPatch = sMesh.boundary()[patchi];
137 const fvPatch& basePatch = vf.mesh().boundary()[basePatchId];
138 const label baseStart = basePatch.start();
139 const label baseSize = basePatch.size();
140
141 labelList directAddressing(subPatch.size());
142
143 forAll(directAddressing, i)
144 {
145 const label baseFacei = faceMap[subPatch.start()+i];
146
147 if (baseFacei >= baseStart && baseFacei < baseStart+baseSize)
148 {
149 directAddressing[i] = baseFacei-baseStart;
150 }
151 else
152 {
153 // Mapped from internal face. Do what? Leave up to
154 // fvPatchField
155 directAddressing[i] = -1;
156 }
157 }
158
159
160 directFvPatchFieldMapper mapper(directAddressing);
161
162 // allowUnmapped : special mode for if we do not want to be
163 // warned for unmapped faces (e.g. from fvMeshDistribute).
164
165 const bool hasUnmapped = mapper.hasUnmapped();
166 if (allowUnmapped)
167 {
168 mapper.hasUnmapped() = false;
169 }
170
171 bf.set
172 (
173 patchi,
175 (
176 vf.boundaryField()[basePatchId],
177 subPatch,
178 result.internalField(),
179 mapper
180 )
181 );
182
183 if (allowUnmapped && hasUnmapped)
184 {
185 // Set unmapped values to zeroGradient. This is the default
186 // action for unmapped fvPatchFields. Note that this bypasses
187 // any special logic for handling unmapped fvPatchFields but
188 // since this is only used inside fvMeshDistribute ...
189
190 tmp<Field<Type>> tfld(bf[patchi].patchInternalField());
191 const Field<Type>& fld = tfld();
192
193 Field<Type> value(bf[patchi]);
194 forAll(directAddressing, i)
195 {
196 if (directAddressing[i] == -1)
197 {
198 value[i] = fld[i];
199 }
200 }
201 bf[patchi].fvPatchField<Type>::operator=(value);
202 }
203 }
204 }
205
206 return tresult;
207}
208
209
210template<class Type>
211bool Foam::functionObjects::cellDecomposer::mapFieldType() const
212{
214
215 const fvMesh& mapRegion =
216 this->mesh_.time().lookupObject<fvMesh>(mapRegion_);
217
218 const labelList patchMap(identity(mapRegion.boundaryMesh().size()));
219
220 const wordList fieldNames
221 (
222 this->mesh_.sortedNames<VolFieldType>(fieldNames_)
223 );
224
225 const bool processed = !fieldNames.empty();
226
227 for (const word& fieldName : fieldNames)
228 {
229 const VolFieldType& field = lookupObject<VolFieldType>(fieldName);
230
231 auto* mapFieldPtr = mapRegion.getObjectPtr<VolFieldType>(fieldName);
232
233 if (!mapFieldPtr)
234 {
235 mapFieldPtr = new VolFieldType
236 (
238 (
239 fieldName,
240 time_.timeName(),
241 mapRegion,
245 ),
246 mapRegion,
247 dimensioned<Type>(field.dimensions(), Zero)
248 );
249
250 mapFieldPtr->store();
251 }
252
253 auto& mappedField = *mapFieldPtr;
254
255 mappedField = interpolate
256 (
257 field,
258 mapRegion,
259 patchMap,
260 mapPtr_().cellMap(),
261 mapPtr_().faceMap(),
262 false //allowUnmapped
263 );
264 Log << " " << fieldName << ": interpolated" << nl;
265
266 //evaluateConstraintTypes(mappedField);
267 }
268
269 return processed;
270}
271
272
273template<class Type>
274bool Foam::functionObjects::cellDecomposer::writeFieldType() const
275{
276 typedef GeometricField<Type, fvPatchField, volMesh> VolFieldType;
277
278 const fvMesh& mapRegion =
279 this->mesh_.time().lookupObject<fvMesh>(mapRegion_);
280
281 const wordList fieldNames
282 (
283 this->mesh_.sortedNames<VolFieldType>(fieldNames_)
284 );
285
286 const bool processed = !fieldNames.empty();
287
288 for (const word& fieldName : fieldNames)
289 {
290 const VolFieldType& mappedField =
291 mapRegion.template lookupObject<VolFieldType>(fieldName);
292
293 mappedField.write();
294
295 Log << " " << fieldName << ": written" << nl;
296 }
297
298 return processed;
299}
300
301
302// ************************************************************************* //
#define Log
Definition PDRblock.C:28
Info<< nl;Info<< "Write faMesh in vtk format:"<< nl;{ vtk::uindirectPatchWriter writer(aMesh.patch(), fileName(aMesh.time().globalPath()/vtkBaseFileName));writer.writeGeometry();globalIndex procAddr(aMesh.nFaces());labelList cellIDs;if(UPstream::master()) { cellIDs.resize(procAddr.totalSize());for(const labelRange &range :procAddr.ranges()) { auto slice=cellIDs.slice(range);slice=identity(range);} } writer.beginCellData(4);writer.writeProcIDs();writer.write("cellID", cellIDs);writer.write("area", aMesh.S().field());writer.write("normal", aMesh.faceAreaNormals());writer.beginPointData(1);writer.write("normal", aMesh.pointAreaNormals());Info<< " "<< writer.output().name()<< nl;}{ vtk::lineWriter writer(aMesh.points(), aMesh.edges(), fileName(aMesh.time().globalPath()/(vtkBaseFileName+"-edges")));writer.writeGeometry();writer.beginCellData(4);writer.writeProcIDs();{ Field< scalar > fld(faMeshTools::flattenEdgeField(aMesh.magLe(), true))
static const this_type & null() noexcept
const Mesh & mesh() const noexcept
Return const reference to mesh.
const dimensionSet & dimensions() const noexcept
Return dimensions.
orientedType oriented() const noexcept
Return oriented type.
virtual bool hasUnmapped() const
Any unmapped values?
Generic templated field type that is much like a Foam::List except that it is expected to hold numeri...
Definition Field.H:172
Generic GeometricField class.
const Internal::FieldType & primitiveField() const noexcept
Return a const-reference to the internal field values.
const Boundary & boundaryField() const noexcept
Return const-reference to the boundary field.
@ REGISTER
Request registration (bool: true).
@ NO_READ
Nothing to be read.
@ NO_WRITE
Ignore writing from objectRegistry::writeObject().
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition IOobject.H:191
const word & name() const noexcept
Return the object name.
Definition IOobjectI.H:205
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition PtrList.H:67
const T * set(const label i) const
Return const pointer to element (can be nullptr), or nullptr for out-of-range access (ie,...
Definition PtrList.H:171
static word timeName(const scalar t, const int precision=precision_)
Return a time name for the given scalar time value formatted with the given precision.
Definition Time.C:714
void size(const label n)
Older name for setAddressableSize.
Definition UList.H:118
label size() const noexcept
The number of entries in the list.
Definition UPtrListI.H:106
Generic dimensioned Type class.
This boundary condition provides an 'empty' condition for reduced dimensions cases,...
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
const Time & time() const
Return the top-level database.
Definition fvMesh.H:360
const fvBoundaryMesh & boundary() const noexcept
Return reference to boundary mesh.
Definition fvMesh.H:395
static const word & calculatedType() noexcept
The type name for calculated patch fields.
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.
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition fvPatch.H:71
virtual label size() const
Patch size is the number of faces, but can be overloaded.
Definition fvPatch.H:242
label start() const noexcept
The patch start within the polyMesh face list.
Definition fvPatch.H:226
const Type & lookupObject(const word &name, const bool recursive=false) const
Lookup and return const reference to the object of the given Type. Fatal if not found or the wrong ty...
Type * getObjectPtr(const word &name, const bool recursive=false) const
Return non-const pointer to the object of the given Type, using a const-cast to have it behave like a...
const polyBoundaryMesh & boundaryMesh() const noexcept
Return boundary mesh.
Definition polyMesh.H:609
A class for managing temporary objects.
Definition tmp.H:75
A class for handling words, derived from Foam::string.
Definition word.H:66
rDeltaTY field()
Pair< int > faceMap(const label facePi, const face &faceP, const label faceNi, const face &faceN)
List< word > wordList
List of word.
Definition fileName.H:60
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.
List< label > labelList
A List of labels.
Definition List.H:62
DirectFieldMapper< fvPatchFieldMapper > directFvPatchFieldMapper
A fvPatchFieldMapper with direct mapping.
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...
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
bool interpolate(const vector &p1, const vector &p2, const vector &o, vector &n, scalar l)
Definition curveTools.C:75
UList< label > labelUList
A UList of labels.
Definition UList.H:75
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299