fvFieldDecomposer.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-2016 OpenFOAM Foundation
9 Copyright (C) 2021-2022 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 "fvFieldDecomposer.H"
30
31// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
32
34
35
36// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
37
39(
40 const labelUList& addressingSlice,
41 const label addressingOffset
42)
43:
44 directAddressing_(addressingSlice)
45{
46 forAll(directAddressing_, i)
47 {
48 // Subtract one to align addressing.
49 directAddressing_[i] -= addressingOffset + 1;
50 }
51}
52
53
54Foam::fvFieldDecomposer::processorVolPatchFieldDecomposer::
55processorVolPatchFieldDecomposer
56(
57 const labelUList& owner, // == mesh.faceOwner()
58 const labelUList& neigh, // == mesh.faceNeighbour()
59 const labelUList& addressingSlice
60)
61:
62 directAddressing_(addressingSlice.size())
63{
64 forAll(directAddressing_, i)
65 {
66 // Subtract one to align addressing.
67 label ai = mag(addressingSlice[i]) - 1;
68
69 if (ai < neigh.size())
70 {
71 // This is a regular face. it has been an internal face
72 // of the original mesh and now it has become a face
73 // on the parallel boundary.
74 // Give face the value of the neighbour.
75
76 if (addressingSlice[i] >= 0)
77 {
78 // I have the owner so use the neighbour value
79 directAddressing_[i] = neigh[ai];
80 }
81 else
82 {
83 directAddressing_[i] = owner[ai];
84 }
85 }
86 else
87 {
88 // This is a face that used to be on a cyclic boundary
89 // but has now become a parallel patch face. I cannot
90 // do the interpolation properly (I would need to look
91 // up the different (face) list of data), so I will
92 // just grab the value from the owner cell
93
94 directAddressing_[i] = owner[ai];
95 }
96 }
97}
98
99
100Foam::fvFieldDecomposer::processorVolPatchFieldDecomposer::
101processorVolPatchFieldDecomposer
102(
103 const polyMesh& mesh,
104 const labelUList& addressingSlice
105)
106:
108 (
109 mesh.faceOwner(),
110 mesh.faceNeighbour(),
111 addressingSlice
112 )
113{}
114
115
116Foam::fvFieldDecomposer::processorSurfacePatchFieldDecomposer::
117processorSurfacePatchFieldDecomposer
118(
119 const labelUList& addressingSlice
120)
121:
122 addressing_(addressingSlice.size()),
123 weights_(addressingSlice.size())
124{
125 forAll(addressing_, i)
126 {
127 addressing_[i].resize(1);
128 weights_[i].resize(1);
129
130 addressing_[i][0] = mag(addressingSlice[i]) - 1;
131 weights_[i][0] = 1;
132 }
133}
134
135
137(
138 const Foam::zero,
139 const fvMesh& procMesh,
140 const labelList& faceAddressing,
141 const labelList& cellAddressing,
142 const labelList& boundaryAddressing
143)
144:
145 procMesh_(procMesh),
146 faceAddressing_(faceAddressing),
147 cellAddressing_(cellAddressing),
148 boundaryAddressing_(boundaryAddressing),
149 // Mappers
150 patchFieldDecomposerPtrs_(),
151 processorVolPatchFieldDecomposerPtrs_(),
152 processorSurfacePatchFieldDecomposerPtrs_(),
153 faceSign_()
154{}
155
156
158(
159 const fvMesh& completeMesh,
160 const fvMesh& procMesh,
161 const labelList& faceAddressing,
162 const labelList& cellAddressing,
163 const labelList& boundaryAddressing
164)
165:
167 (
168 zero{},
169 procMesh,
170 faceAddressing,
171 cellAddressing,
172 boundaryAddressing
173 )
174{
175 reset(completeMesh);
176}
177
178
180(
181 const List<labelRange>& boundaryRanges,
182 const labelUList& faceOwner,
183 const labelUList& faceNeighbour,
184
185 const fvMesh& procMesh,
186 const labelList& faceAddressing,
187 const labelList& cellAddressing,
188 const labelList& boundaryAddressing
189)
190:
192 (
193 zero{},
194 procMesh,
195 faceAddressing,
196 cellAddressing,
197 boundaryAddressing
198 )
199{
200 reset(boundaryRanges, faceOwner, faceNeighbour);
201}
202
203
204// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
205
207{
208 return patchFieldDecomposerPtrs_.empty();
209}
210
211
213{
214 patchFieldDecomposerPtrs_.clear();
215 processorVolPatchFieldDecomposerPtrs_.clear();
216 processorSurfacePatchFieldDecomposerPtrs_.clear();
217 faceSign_.clear();
218}
219
220
222(
223 const List<labelRange>& boundaryRanges,
224 const labelUList& faceOwner,
225 const labelUList& faceNeighbour
226)
227{
228 clear();
229 const label nMappers = procMesh_.boundary().size();
230 patchFieldDecomposerPtrs_.resize(nMappers);
231 processorVolPatchFieldDecomposerPtrs_.resize(nMappers);
232 processorSurfacePatchFieldDecomposerPtrs_.resize(nMappers);
233 faceSign_.resize(nMappers);
234
235 forAll(boundaryAddressing_, patchi)
236 {
237 const label oldPatchi = boundaryAddressing_[patchi];
238 const fvPatch& fvp = procMesh_.boundary()[patchi];
239 const labelSubList localPatchSlice(fvp.patchSlice(faceAddressing_));
240
241 if
242 (
243 oldPatchi >= 0
244 && !isA<processorLduInterface>(procMesh_.boundary()[patchi])
245 )
246 {
247 patchFieldDecomposerPtrs_.set
248 (
249 patchi,
251 (
252 localPatchSlice,
253 boundaryRanges[oldPatchi].start()
254 )
255 );
256 }
257 else
258 {
259 processorVolPatchFieldDecomposerPtrs_.set
260 (
261 patchi,
263 (
264 faceOwner,
265 faceNeighbour,
266 localPatchSlice
267 )
268 );
269
270 processorSurfacePatchFieldDecomposerPtrs_.set
271 (
272 patchi,
274 (
275 static_cast<const labelUList&>(localPatchSlice)
276 )
277 );
278
279 faceSign_.set
280 (
281 patchi,
282 new scalarField(localPatchSlice.size())
283 );
284
285 {
286 scalarField& s = faceSign_[patchi];
287 forAll(s, i)
288 {
289 s[i] = sign(localPatchSlice[i]);
290 }
291 }
292 }
293 }
294}
295
296
297void Foam::fvFieldDecomposer::reset(const fvMesh& completeMesh)
298{
299 clear();
300 const label nMappers = procMesh_.boundary().size();
301 patchFieldDecomposerPtrs_.resize(nMappers);
302 processorVolPatchFieldDecomposerPtrs_.resize(nMappers);
303 processorSurfacePatchFieldDecomposerPtrs_.resize(nMappers);
304 faceSign_.resize(nMappers);
305
306 forAll(boundaryAddressing_, patchi)
307 {
308 const label oldPatchi = boundaryAddressing_[patchi];
309 const fvPatch& fvp = procMesh_.boundary()[patchi];
310 const labelSubList localPatchSlice(fvp.patchSlice(faceAddressing_));
311
312 if
313 (
314 oldPatchi >= 0
315 && !isA<processorLduInterface>(procMesh_.boundary()[patchi])
316 )
317 {
318 patchFieldDecomposerPtrs_.set
319 (
320 patchi,
322 (
323 localPatchSlice,
324 completeMesh.boundaryMesh()[oldPatchi].start()
325 )
326 );
327 }
328 else
329 {
330 processorVolPatchFieldDecomposerPtrs_.set
331 (
332 patchi,
334 (
335 completeMesh,
336 localPatchSlice
337 )
338 );
339
340 processorSurfacePatchFieldDecomposerPtrs_.set
341 (
342 patchi,
344 (
345 static_cast<const labelUList&>(localPatchSlice)
346 )
347 );
348
349 faceSign_.set
350 (
351 patchi,
352 new scalarField(localPatchSlice.size())
353 );
354
355 {
356 scalarField& s = faceSign_[patchi];
357 forAll(s, i)
358 {
359 s[i] = sign(localPatchSlice[i]);
360 }
361 }
362 }
363 }
364}
365
366
367// ************************************************************************* //
void resize(const label len)
Adjust allocated size of list.
Definition: ListI.H:139
A List obtained as a section of another List.
Definition: SubList.H:70
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:114
void reset()
Reset to defaults.
Processor patch field decomposer class. Surface field is assumed.
Processor patch field decomposer class. Maps either owner or.
Finite Volume volume and surface field decomposer.
void reset(const fvMesh &completeMesh)
Reset mappers using information from the complete mesh.
bool empty() const
True if no mappers have been allocated.
static int verbose_
Output verbosity when writing.
void clear()
Remove all mappers.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:91
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:71
const List< T >::subList patchSlice(const List< T > &l) const
Slice list to patch.
Definition: fvPatch.H:216
label start() const noexcept
The start label of boundary faces in the polyMesh face list.
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:81
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:456
A class representing the concept of 0 (zero) that can be used to avoid manipulating objects known to ...
Definition: zero.H:63
patchWriters clear()
dynamicFvMesh & mesh
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;forAll(lagrangianScalarNames, i){ word name=lagrangianScalarNames[i];IOField< scalar > s(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
dimensionedScalar sign(const dimensionedScalar &ds)
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333