initial commit for version 1.6.x patch release
[OpenFOAM-1.6.x.git] / applications / utilities / parallelProcessing / decomposePar / fvFieldDecomposer.C
blobf5a6451c1266b6a0e02c165c3501c037c7f85b64
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 1991-2009 OpenCFD Ltd.
6      \\/     M anipulation  |
7 -------------------------------------------------------------------------------
8 License
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 the
13     Free Software Foundation; either version 2 of the License, or (at your
14     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
19     for more details.
21     You should have received a copy of the GNU General Public License
22     along with OpenFOAM; if not, write to the Free Software Foundation,
23     Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301 USA
25 \*---------------------------------------------------------------------------*/
27 #include "fvFieldDecomposer.H"
29 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
31 namespace Foam
34 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
36 fvFieldDecomposer::patchFieldDecomposer::patchFieldDecomposer
38     const unallocLabelList& addressingSlice,
39     const label addressingOffset
42     directAddressing_(addressingSlice)
44     forAll (directAddressing_, i)
45     {
46         // Subtract one to align addressing.  
47         directAddressing_[i] -= addressingOffset + 1;
48     }
52 fvFieldDecomposer::processorVolPatchFieldDecomposer::
53 processorVolPatchFieldDecomposer
55     const fvMesh& mesh,
56     const unallocLabelList& addressingSlice
59     directAddressing_(addressingSlice.size())
61     const labelList& own = mesh.faceOwner();
62     const labelList& neighb = mesh.faceNeighbour();
64     forAll (directAddressing_, i)
65     {
66         // Subtract one to align addressing.  
67         label ai = mag(addressingSlice[i]) - 1;
69         if (ai < neighb.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.
76             if (addressingSlice[i] >= 0)
77             {
78                 // I have the owner so use the neighbour value
79                 directAddressing_[i] = neighb[ai];
80             }
81             else
82             {
83                 directAddressing_[i] = own[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
94             directAddressing_[i] = own[ai];
95         }
96     }
100 fvFieldDecomposer::processorSurfacePatchFieldDecomposer::
101 processorSurfacePatchFieldDecomposer
103     const unallocLabelList& addressingSlice
106     addressing_(addressingSlice.size()),
107     weights_(addressingSlice.size())
109     forAll (addressing_, i)
110     {
111         addressing_[i].setSize(1);
112         weights_[i].setSize(1);
114         addressing_[i][0] = mag(addressingSlice[i]) - 1;
115         weights_[i][0] = sign(addressingSlice[i]);
116     }
120 fvFieldDecomposer::fvFieldDecomposer
122     const fvMesh& completeMesh,
123     const fvMesh& procMesh,
124     const labelList& faceAddressing,
125     const labelList& cellAddressing,
126     const labelList& boundaryAddressing
129     completeMesh_(completeMesh),
130     procMesh_(procMesh),
131     faceAddressing_(faceAddressing),
132     cellAddressing_(cellAddressing),
133     boundaryAddressing_(boundaryAddressing),
134     patchFieldDecomposerPtrs_
135     (
136         procMesh_.boundary().size(),
137         static_cast<patchFieldDecomposer*>(NULL)
138     ),
139     processorVolPatchFieldDecomposerPtrs_
140     (
141         procMesh_.boundary().size(),
142         static_cast<processorVolPatchFieldDecomposer*>(NULL)
143     ),
144     processorSurfacePatchFieldDecomposerPtrs_
145     (
146         procMesh_.boundary().size(),
147         static_cast<processorSurfacePatchFieldDecomposer*>(NULL)
148     )
150     forAll (boundaryAddressing_, patchi)
151     {
152         if (boundaryAddressing_[patchi] >= 0)
153         {
154             patchFieldDecomposerPtrs_[patchi] = new patchFieldDecomposer
155             (
156                 procMesh_.boundary()[patchi].patchSlice(faceAddressing_),
157                 completeMesh_.boundaryMesh()
158                 [
159                     boundaryAddressing_[patchi]
160                 ].start()
161             );
162         }
163         else
164         {
165             processorVolPatchFieldDecomposerPtrs_[patchi] = 
166                 new processorVolPatchFieldDecomposer
167                 (
168                     completeMesh_,
169                     procMesh_.boundary()[patchi].patchSlice(faceAddressing_)
170                 );
172             processorSurfacePatchFieldDecomposerPtrs_[patchi] = 
173                 new processorSurfacePatchFieldDecomposer
174                 (
175                     static_cast<const unallocLabelList&>
176                     (
177                         procMesh_.boundary()[patchi].patchSlice
178                         (
179                             faceAddressing_
180                         )
181                     )
182                 );
183         }
184     }
188 // * * * * * * * * * * * * * * * * Destructor  * * * * * * * * * * * * * * * //
190 fvFieldDecomposer::~fvFieldDecomposer()
192     forAll (patchFieldDecomposerPtrs_, patchi)
193     {
194         if (patchFieldDecomposerPtrs_[patchi])
195         {
196             delete patchFieldDecomposerPtrs_[patchi];
197         }
198     }
200     forAll (processorVolPatchFieldDecomposerPtrs_, patchi)
201     {
202         if (processorVolPatchFieldDecomposerPtrs_[patchi])
203         {
204             delete processorVolPatchFieldDecomposerPtrs_[patchi];
205         }
206     }
208     forAll (processorSurfacePatchFieldDecomposerPtrs_, patchi)
209     {
210         if (processorSurfacePatchFieldDecomposerPtrs_[patchi])
211         {
212             delete processorSurfacePatchFieldDecomposerPtrs_[patchi];
213         }
214     }
218 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
220 } // End namespace Foam
222 // ************************************************************************* //