1 /*---------------------------------------------------------------------------*\
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
5 \\ / A nd | Copyright (C) 2011 OpenFOAM Foundation
7 -------------------------------------------------------------------------------
9 This file is part of OpenFOAM.
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.
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
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/>.
24 \*---------------------------------------------------------------------------*/
26 #include "fvFieldDecomposer.H"
29 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
31 Foam::fvFieldDecomposer::patchFieldDecomposer::patchFieldDecomposer
33 const labelUList& addressingSlice,
34 const label addressingOffset
37 directAddressing_(addressingSlice)
39 forAll(directAddressing_, i)
41 // Subtract one to align addressing.
42 directAddressing_[i] -= addressingOffset + 1;
47 Foam::fvFieldDecomposer::processorVolPatchFieldDecomposer::
48 processorVolPatchFieldDecomposer
51 const labelUList& addressingSlice
54 directAddressing_(addressingSlice.size())
56 const labelList& own = mesh.faceOwner();
57 const labelList& neighb = mesh.faceNeighbour();
59 forAll(directAddressing_, i)
61 // Subtract one to align addressing.
62 label ai = mag(addressingSlice[i]) - 1;
64 if (ai < neighb.size())
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.
71 if (addressingSlice[i] >= 0)
73 // I have the owner so use the neighbour value
74 directAddressing_[i] = neighb[ai];
78 directAddressing_[i] = own[ai];
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
89 directAddressing_[i] = own[ai];
95 Foam::fvFieldDecomposer::processorSurfacePatchFieldDecomposer::
96 processorSurfacePatchFieldDecomposer
98 const labelUList& addressingSlice
101 addressing_(addressingSlice.size()),
102 weights_(addressingSlice.size())
104 forAll(addressing_, i)
106 addressing_[i].setSize(1);
107 weights_[i].setSize(1);
109 addressing_[i][0] = mag(addressingSlice[i]) - 1;
110 weights_[i][0] = sign(addressingSlice[i]);
115 Foam::fvFieldDecomposer::fvFieldDecomposer
117 const fvMesh& completeMesh,
118 const fvMesh& procMesh,
119 const labelList& faceAddressing,
120 const labelList& cellAddressing,
121 const labelList& boundaryAddressing
124 completeMesh_(completeMesh),
126 faceAddressing_(faceAddressing),
127 cellAddressing_(cellAddressing),
128 boundaryAddressing_(boundaryAddressing),
129 patchFieldDecomposerPtrs_
131 procMesh_.boundary().size(),
132 static_cast<patchFieldDecomposer*>(NULL)
134 processorVolPatchFieldDecomposerPtrs_
136 procMesh_.boundary().size(),
137 static_cast<processorVolPatchFieldDecomposer*>(NULL)
139 processorSurfacePatchFieldDecomposerPtrs_
141 procMesh_.boundary().size(),
142 static_cast<processorSurfacePatchFieldDecomposer*>(NULL)
145 forAll(boundaryAddressing_, patchi)
149 boundaryAddressing_[patchi] >= 0
150 && !isA<processorLduInterface>(procMesh.boundary()[patchi])
153 patchFieldDecomposerPtrs_[patchi] = new patchFieldDecomposer
155 procMesh_.boundary()[patchi].patchSlice(faceAddressing_),
156 completeMesh_.boundaryMesh()
158 boundaryAddressing_[patchi]
164 processorVolPatchFieldDecomposerPtrs_[patchi] =
165 new processorVolPatchFieldDecomposer
168 procMesh_.boundary()[patchi].patchSlice(faceAddressing_)
171 processorSurfacePatchFieldDecomposerPtrs_[patchi] =
172 new processorSurfacePatchFieldDecomposer
174 static_cast<const labelUList&>
176 procMesh_.boundary()[patchi].patchSlice
187 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
189 Foam::fvFieldDecomposer::~fvFieldDecomposer()
191 forAll(patchFieldDecomposerPtrs_, patchi)
193 if (patchFieldDecomposerPtrs_[patchi])
195 delete patchFieldDecomposerPtrs_[patchi];
199 forAll(processorVolPatchFieldDecomposerPtrs_, patchi)
201 if (processorVolPatchFieldDecomposerPtrs_[patchi])
203 delete processorVolPatchFieldDecomposerPtrs_[patchi];
207 forAll(processorSurfacePatchFieldDecomposerPtrs_, patchi)
209 if (processorSurfacePatchFieldDecomposerPtrs_[patchi])
211 delete processorSurfacePatchFieldDecomposerPtrs_[patchi];
216 // ************************************************************************* //