fvFieldDecomposer.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2018 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 #include "fvFieldDecomposer.H"
27 
28 
29 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
30 
32 (
33  const labelUList& addressingSlice,
34  const label addressingOffset
35 )
36 :
37  directAddressing_(addressingSlice)
38 {
39  forAll(directAddressing_, i)
40  {
41  // Subtract one to align addressing.
42  directAddressing_[i] -= addressingOffset + 1;
43  }
44 }
45 
46 
49 (
50  const fvMesh& mesh,
51  const labelUList& addressingSlice
52 )
53 :
54  directAddressing_(addressingSlice.size())
55 {
56  const labelList& own = mesh.faceOwner();
57  const labelList& neighb = mesh.faceNeighbour();
58 
59  forAll(directAddressing_, i)
60  {
61  // Subtract one to align addressing.
62  label ai = mag(addressingSlice[i]) - 1;
63 
64  if (ai < neighb.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] = neighb[ai];
75  }
76  else
77  {
78  directAddressing_[i] = own[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] = own[ai];
90  }
91  }
92 }
93 
94 
97 (
98  const labelUList& addressingSlice
99 )
100 :
101  addressing_(addressingSlice.size()),
102  weights_(addressingSlice.size())
103 {
104  forAll(addressing_, i)
105  {
106  addressing_[i].setSize(1);
107  weights_[i].setSize(1);
108 
109  addressing_[i][0] = mag(addressingSlice[i]) - 1;
110  weights_[i][0] = sign(addressingSlice[i]);
111  }
112 }
113 
114 
115 Foam::fvFieldDecomposer::fvFieldDecomposer
116 (
117  const fvMesh& completeMesh,
118  const fvMesh& procMesh,
119  const labelList& faceAddressing,
120  const labelList& cellAddressing,
121  const labelList& boundaryAddressing
122 )
123 :
124  completeMesh_(completeMesh),
125  procMesh_(procMesh),
126  faceAddressing_(faceAddressing),
127  cellAddressing_(cellAddressing),
128  boundaryAddressing_(boundaryAddressing),
129  patchFieldDecomposerPtrs_
130  (
131  procMesh_.boundary().size(),
132  static_cast<patchFieldDecomposer*>(nullptr)
133  ),
134  processorVolPatchFieldDecomposerPtrs_
135  (
136  procMesh_.boundary().size(),
137  static_cast<processorVolPatchFieldDecomposer*>(nullptr)
138  ),
139  processorSurfacePatchFieldDecomposerPtrs_
140  (
141  procMesh_.boundary().size(),
142  static_cast<processorSurfacePatchFieldDecomposer*>(nullptr)
143  )
144 {
145  forAll(boundaryAddressing_, patchi)
146  {
147  if
148  (
149  boundaryAddressing_[patchi] >= 0
150  && !isA<processorLduInterface>(procMesh.boundary()[patchi])
151  )
152  {
153  patchFieldDecomposerPtrs_[patchi] = new patchFieldDecomposer
154  (
155  procMesh_.boundary()[patchi].patchSlice(faceAddressing_),
156  completeMesh_.boundaryMesh()
157  [
158  boundaryAddressing_[patchi]
159  ].start()
160  );
161  }
162  else
163  {
164  processorVolPatchFieldDecomposerPtrs_[patchi] =
166  (
167  completeMesh_,
168  procMesh_.boundary()[patchi].patchSlice(faceAddressing_)
169  );
170 
171  processorSurfacePatchFieldDecomposerPtrs_[patchi] =
173  (
174  static_cast<const labelUList&>
175  (
176  procMesh_.boundary()[patchi].patchSlice
177  (
178  faceAddressing_
179  )
180  )
181  );
182  }
183  }
184 }
185 
186 
187 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
188 
190 {
191  forAll(patchFieldDecomposerPtrs_, patchi)
192  {
193  if (patchFieldDecomposerPtrs_[patchi])
194  {
195  delete patchFieldDecomposerPtrs_[patchi];
196  }
197  }
198 
199  forAll(processorVolPatchFieldDecomposerPtrs_, patchi)
200  {
201  if (processorVolPatchFieldDecomposerPtrs_[patchi])
202  {
203  delete processorVolPatchFieldDecomposerPtrs_[patchi];
204  }
205  }
206 
207  forAll(processorSurfacePatchFieldDecomposerPtrs_, patchi)
208  {
209  if (processorSurfacePatchFieldDecomposerPtrs_[patchi])
210  {
211  delete processorSurfacePatchFieldDecomposerPtrs_[patchi];
212  }
213  }
214 }
215 
216 // ************************************************************************* //
dimensionedScalar sign(const dimensionedScalar &ds)
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:428
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
virtual const labelList & faceNeighbour() const
Return face neighbour.
Definition: polyMesh.C:1047
processorSurfacePatchFieldDecomposer(const labelUList &addressingSlice)
Construct given addressing.
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:163
processorVolPatchFieldDecomposer(const fvMesh &mesh, const labelUList &addressingSlice)
Construct given addressing.
virtual const labelList & faceOwner() const
Return face owner.
Definition: polyMesh.C:1041
Processor patch field decomposer class. Surface field is assumed.
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:61
Processor patch field decomposer class. Maps either owner or.
patchFieldDecomposer(const labelUList &addressingSlice, const label addressingOffset)
Construct given addressing.
label patchi
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
~fvFieldDecomposer()
Destructor.
dimensioned< scalar > mag(const dimensioned< Type > &)
label size() const
Return the number of elements in the UList.
Definition: UListI.H:299
const fvBoundaryMesh & boundary() const
Return reference to boundary mesh.
Definition: fvMesh.C:545