1 /*---------------------------------------------------------------------------*\
3 \\ / F ield | foam-extend: Open Source CFD
4 \\ / O peration | Version: 3.2
5 \\ / A nd | Web: http://www.foam-extend.org
6 \\/ M anipulation | For copyright notice see file Copyright
7 -------------------------------------------------------------------------------
9 This file is part of foam-extend.
11 foam-extend is free software: you can redistribute it and/or modify it
12 under the terms of the GNU General Public License as published by the
13 Free Software Foundation, either version 3 of the License, or (at your
14 option) any later version.
16 foam-extend is distributed in the hope that it will be useful, but
17 WITHOUT ANY WARRANTY; without even the implied warranty of
18 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
19 General Public License for more details.
21 You should have received a copy of the GNU General Public License
22 along with foam-extend. If not, see <http://www.gnu.org/licenses/>.
24 \*---------------------------------------------------------------------------*/
26 #include "fvFieldDecomposer.H"
28 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
33 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
35 fvFieldDecomposer::patchFieldDecomposer::patchFieldDecomposer
37 const label sizeBeforeMapping,
38 const unallocLabelList& addressingSlice,
39 const label addressingOffset
42 sizeBeforeMapping_(sizeBeforeMapping),
43 directAddressing_(addressingSlice)
45 forAll (directAddressing_, i)
47 // Subtract one to align addressing. HJ, 5/Dec/2001
48 directAddressing_[i] -= addressingOffset + 1;
53 fvFieldDecomposer::processorVolPatchFieldDecomposer::
54 processorVolPatchFieldDecomposer
57 const unallocLabelList& addressingSlice
60 sizeBeforeMapping_(mesh.nCells()),
61 addressing_(addressingSlice.size()),
62 weights_(addressingSlice.size())
64 const scalarField& weights = mesh.weights().internalField();
65 const labelList& own = mesh.faceOwner();
66 const labelList& neighb = mesh.faceNeighbour();
68 forAll (addressing_, i)
70 // Subtract one to align addressing. HJ, 5/Dec/2001
71 label ai = mag(addressingSlice[i]) - 1;
73 if (ai < neighb.size())
75 // This is a regular face. it has been an internal face
76 // of the original mesh and now it has become a face
77 // on the parallel boundary
78 addressing_[i].setSize(2);
79 weights_[i].setSize(2);
81 addressing_[i][0] = own[ai];
82 addressing_[i][1] = neighb[ai];
84 weights_[i][0] = weights[ai];
85 weights_[i][1] = 1.0 - weights[ai];
89 // This is a face that used to be on a cyclic boundary
90 // but has now become a parallel patch face. I cannot
91 // do the interpolation properly (I would need to look
92 // up the different (face) list of data), so I will
93 // just grab the value from the owner cell
95 addressing_[i].setSize(1);
96 weights_[i].setSize(1);
98 addressing_[i][0] = own[ai];
100 weights_[i][0] = 1.0;
106 fvFieldDecomposer::processorSurfacePatchFieldDecomposer::
107 processorSurfacePatchFieldDecomposer
109 label sizeBeforeMapping,
110 const unallocLabelList& addressingSlice
113 sizeBeforeMapping_(sizeBeforeMapping),
114 addressing_(addressingSlice.size()),
115 weights_(addressingSlice.size())
117 forAll (addressing_, i)
119 addressing_[i].setSize(1);
120 weights_[i].setSize(1);
122 addressing_[i][0] = mag(addressingSlice[i]) - 1;
123 weights_[i][0] = sign(addressingSlice[i]);
128 fvFieldDecomposer::fvFieldDecomposer
130 const fvMesh& completeMesh,
131 const fvMesh& procMesh,
132 const labelList& faceAddressing,
133 const labelList& cellAddressing,
134 const labelList& boundaryAddressing
137 completeMesh_(completeMesh),
139 faceAddressing_(faceAddressing),
140 cellAddressing_(cellAddressing),
141 boundaryAddressing_(boundaryAddressing),
142 patchFieldDecomposerPtrs_
144 procMesh_.boundary().size(),
145 static_cast<patchFieldDecomposer*>(NULL)
147 processorVolPatchFieldDecomposerPtrs_
149 procMesh_.boundary().size(),
150 static_cast<processorVolPatchFieldDecomposer*>(NULL)
152 processorSurfacePatchFieldDecomposerPtrs_
154 procMesh_.boundary().size(),
155 static_cast<processorSurfacePatchFieldDecomposer*>(NULL)
158 forAll (boundaryAddressing_, patchi)
160 if (boundaryAddressing_[patchi] >= 0)
162 patchFieldDecomposerPtrs_[patchi] = new patchFieldDecomposer
164 completeMesh_.boundary()[boundaryAddressing_[patchi]].size(),
165 procMesh_.boundary()[patchi].patchSlice(faceAddressing_),
166 completeMesh_.boundaryMesh()
168 boundaryAddressing_[patchi]
174 processorVolPatchFieldDecomposerPtrs_[patchi] =
175 new processorVolPatchFieldDecomposer
178 procMesh_.boundary()[patchi].patchSlice(faceAddressing_)
181 processorSurfacePatchFieldDecomposerPtrs_[patchi] =
182 new processorSurfacePatchFieldDecomposer
184 procMesh_.boundary()[patchi].size(),
185 static_cast<const unallocLabelList&>
187 procMesh_.boundary()[patchi].patchSlice
198 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
200 fvFieldDecomposer::~fvFieldDecomposer()
202 forAll (patchFieldDecomposerPtrs_, patchi)
204 if (patchFieldDecomposerPtrs_[patchi])
206 delete patchFieldDecomposerPtrs_[patchi];
210 forAll (processorVolPatchFieldDecomposerPtrs_, patchi)
212 if (processorVolPatchFieldDecomposerPtrs_[patchi])
214 delete processorVolPatchFieldDecomposerPtrs_[patchi];
218 forAll (processorSurfacePatchFieldDecomposerPtrs_, patchi)
220 if (processorSurfacePatchFieldDecomposerPtrs_[patchi])
222 delete processorSurfacePatchFieldDecomposerPtrs_[patchi];
228 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
230 } // End namespace Foam
232 // ************************************************************************* //