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