initial commit for version 1.6.x patch release
[OpenFOAM-1.6.x.git] / applications / utilities / parallelProcessing / decomposePar / fvFieldDecomposerDecomposeFields.C
blob512ea77d1d7473017a06395001b8413f1ef79de4
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"
28 #include "processorFvPatchField.H"
29 #include "processorFvsPatchField.H"
31 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
33 namespace Foam
36 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
38 template<class Type>
39 tmp<GeometricField<Type, fvPatchField, volMesh> >
40 fvFieldDecomposer::decomposeField
42     const GeometricField<Type, fvPatchField, volMesh>& field
43 ) const
45     // Create and map the internal field values
46     Field<Type> internalField(field.internalField(), cellAddressing_);
48     // Create and map the patch field values
49     PtrList<fvPatchField<Type> > patchFields(boundaryAddressing_.size());
51     forAll (boundaryAddressing_, patchi)
52     {
53         if (boundaryAddressing_[patchi] >= 0)
54         {
55             patchFields.set
56             (
57                 patchi,
58                 fvPatchField<Type>::New
59                 (
60                     field.boundaryField()[boundaryAddressing_[patchi]],
61                     procMesh_.boundary()[patchi],
62                     DimensionedField<Type, volMesh>::null(),
63                     *patchFieldDecomposerPtrs_[patchi]
64                 )
65             );
66         }
67         else
68         {
69             patchFields.set
70             (
71                 patchi,
72                 new processorFvPatchField<Type>
73                 (
74                     procMesh_.boundary()[patchi],
75                     DimensionedField<Type, volMesh>::null(),
76                     Field<Type>
77                     (
78                         field.internalField(),
79                         *processorVolPatchFieldDecomposerPtrs_[patchi]
80                     )
81                 )
82             );
83         }
84     }
86     // Create the field for the processor
87     return tmp<GeometricField<Type, fvPatchField, volMesh> >
88     (
89         new GeometricField<Type, fvPatchField, volMesh>
90         (
91             IOobject
92             (
93                 field.name(),
94                 procMesh_.time().timeName(),
95                 procMesh_,
96                 IOobject::NO_READ,
97                 IOobject::NO_WRITE
98             ),
99             procMesh_,
100             field.dimensions(),
101             internalField,
102             patchFields
103         )
104     );
108 template<class Type>
109 tmp<GeometricField<Type, fvsPatchField, surfaceMesh> >
110 fvFieldDecomposer::decomposeField
112     const GeometricField<Type, fvsPatchField, surfaceMesh>& field
113 ) const
115     labelList mapAddr
116     (
117         labelList::subList
118         (
119             faceAddressing_,
120             procMesh_.nInternalFaces()
121         )
122     );
123     forAll (mapAddr, i)
124     {
125         mapAddr[i] -= 1;
126     }
128     // Create and map the internal field values
129     Field<Type> internalField
130     (
131         field.internalField(),
132         mapAddr
133     );
135     // Problem with addressing when a processor patch picks up both internal
136     // faces and faces from cyclic boundaries. This is a bit of a hack, but
137     // I cannot find a better solution without making the internal storage
138     // mechanism for surfaceFields correspond to the one of faces in polyMesh
139     // (i.e. using slices)
140     Field<Type> allFaceField(field.mesh().nFaces());
142     forAll (field.internalField(), i)
143     {
144         allFaceField[i] = field.internalField()[i];
145     }
147     forAll (field.boundaryField(), patchi)
148     {
149         const Field<Type> & p = field.boundaryField()[patchi];
151         const label patchStart = field.mesh().boundaryMesh()[patchi].start();
153         forAll (p, i)
154         {
155             allFaceField[patchStart + i] = p[i];
156         }
157     }
159     // Create and map the patch field values
160     PtrList<fvsPatchField<Type> > patchFields(boundaryAddressing_.size());
162     forAll (boundaryAddressing_, patchi)
163     {
164         if (boundaryAddressing_[patchi] >= 0)
165         {
166             patchFields.set
167             (
168                 patchi,
169                 fvsPatchField<Type>::New
170                 (
171                     field.boundaryField()[boundaryAddressing_[patchi]],
172                     procMesh_.boundary()[patchi],
173                     DimensionedField<Type, surfaceMesh>::null(),
174                     *patchFieldDecomposerPtrs_[patchi]
175                 )
176             );
177         }
178         else
179         {
180             patchFields.set
181             (
182                 patchi,
183                 new processorFvsPatchField<Type>
184                 (
185                     procMesh_.boundary()[patchi],
186                     DimensionedField<Type, surfaceMesh>::null(),
187                     Field<Type>
188                     (
189                         allFaceField,
190                         *processorSurfacePatchFieldDecomposerPtrs_[patchi]
191                     )
192                 )
193             );
194         }
195     }
197     // Create the field for the processor
198     return tmp<GeometricField<Type, fvsPatchField, surfaceMesh> >
199     (
200         new GeometricField<Type, fvsPatchField, surfaceMesh>
201         (
202             IOobject
203             (
204                 field.name(),
205                 procMesh_.time().timeName(),
206                 procMesh_,
207                 IOobject::NO_READ,
208                 IOobject::NO_WRITE
209             ),
210             procMesh_,
211             field.dimensions(),
212             internalField,
213             patchFields
214         )
215     );
219 template<class GeoField>
220 void fvFieldDecomposer::decomposeFields
222     const PtrList<GeoField>& fields
223 ) const
225     forAll (fields, fieldI)
226     {
227         decomposeField(fields[fieldI])().write();
228     }
232 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
234 } // End namespace Foam
236 // ************************************************************************* //