Forward compatibility: flex
[foam-extend-3.2.git] / applications / utilities / parallelProcessing / decomposePar / fvFieldDecomposerDecomposeFields.C
blob1a8009441c00a8c199686460e60c3685d33d31d2
1 /*---------------------------------------------------------------------------*\
2   =========                 |
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 -------------------------------------------------------------------------------
8 License
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"
27 #include "processorFvPatchField.H"
28 #include "processorFvsPatchField.H"
30 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
32 namespace Foam
35 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
37 template<class Type>
38 tmp<GeometricField<Type, fvPatchField, volMesh> >
39 fvFieldDecomposer::decomposeField
41     const GeometricField<Type, fvPatchField, volMesh>& field
42 ) const
44     // Create and map the internal field values
45     Field<Type> internalField(field.internalField(), cellAddressing_);
47     // Create and map the patch field values
48     PtrList<fvPatchField<Type> > patchFields(boundaryAddressing_.size());
50     forAll (boundaryAddressing_, patchi)
51     {
52         if (boundaryAddressing_[patchi] >= 0)
53         {
54             patchFields.set
55             (
56                 patchi,
57                 fvPatchField<Type>::New
58                 (
59                     field.boundaryField()[boundaryAddressing_[patchi]],
60                     procMesh_.boundary()[patchi],
61                     DimensionedField<Type, volMesh>::null(),
62                     *patchFieldDecomposerPtrs_[patchi]
63                 )
64             );
65         }
66         else
67         {
68             patchFields.set
69             (
70                 patchi,
71                 new processorFvPatchField<Type>
72                 (
73                     procMesh_.boundary()[patchi],
74                     DimensionedField<Type, volMesh>::null(),
75                     Field<Type>
76                     (
77                         field.internalField(),
78                         *processorVolPatchFieldDecomposerPtrs_[patchi]
79                     )
80                 )
81             );
82         }
83     }
85     // Create the field for the processor
86     return tmp<GeometricField<Type, fvPatchField, volMesh> >
87     (
88         new GeometricField<Type, fvPatchField, volMesh>
89         (
90             IOobject
91             (
92                 field.name(),
93                 procMesh_.time().timeName(),
94                 procMesh_,
95                 IOobject::NO_READ,
96                 IOobject::NO_WRITE
97             ),
98             procMesh_,
99             field.dimensions(),
100             internalField,
101             patchFields
102         )
103     );
107 template<class Type>
108 tmp<GeometricField<Type, fvsPatchField, surfaceMesh> >
109 fvFieldDecomposer::decomposeField
111     const GeometricField<Type, fvsPatchField, surfaceMesh>& field
112 ) const
114     labelList mapAddr
115     (
116         labelList::subList
117         (
118             faceAddressing_,
119             procMesh_.nInternalFaces()
120         )
121     );
122     forAll (mapAddr, i)
123     {
124         mapAddr[i] -= 1;
125     }
127     // Create and map the internal field values
128     Field<Type> internalField
129     (
130         field.internalField(),
131         mapAddr
132     );
134     // Problem with addressing when a processor patch picks up both internal
135     // faces and faces from cyclic boundaries. This is a bit of a hack, but
136     // I cannot find a better solution without making the internal storage
137     // mechanism for surfaceFields correspond to the one of faces in polyMesh
138     // (i.e. using slices)
139     Field<Type> allFaceField(field.mesh().nFaces());
141     forAll (field.internalField(), i)
142     {
143         allFaceField[i] = field.internalField()[i];
144     }
146     forAll (field.boundaryField(), patchi)
147     {
148         const Field<Type> & p = field.boundaryField()[patchi];
150         const label patchStart = field.mesh().boundaryMesh()[patchi].start();
152         forAll (p, i)
153         {
154             allFaceField[patchStart + i] = p[i];
155         }
156     }
158     // Create and map the patch field values
159     PtrList<fvsPatchField<Type> > patchFields(boundaryAddressing_.size());
161     forAll (boundaryAddressing_, patchi)
162     {
163         if (boundaryAddressing_[patchi] >= 0)
164         {
165             patchFields.set
166             (
167                 patchi,
168                 fvsPatchField<Type>::New
169                 (
170                     field.boundaryField()[boundaryAddressing_[patchi]],
171                     procMesh_.boundary()[patchi],
172                     DimensionedField<Type, surfaceMesh>::null(),
173                     *patchFieldDecomposerPtrs_[patchi]
174                 )
175             );
176         }
177         else
178         {
179             patchFields.set
180             (
181                 patchi,
182                 new processorFvsPatchField<Type>
183                 (
184                     procMesh_.boundary()[patchi],
185                     DimensionedField<Type, surfaceMesh>::null(),
186                     Field<Type>
187                     (
188                         allFaceField,
189                         *processorSurfacePatchFieldDecomposerPtrs_[patchi]
190                     )
191                 )
192             );
193         }
194     }
196     // Create the field for the processor
197     return tmp<GeometricField<Type, fvsPatchField, surfaceMesh> >
198     (
199         new GeometricField<Type, fvsPatchField, surfaceMesh>
200         (
201             IOobject
202             (
203                 field.name(),
204                 procMesh_.time().timeName(),
205                 procMesh_,
206                 IOobject::NO_READ,
207                 IOobject::NO_WRITE
208             ),
209             procMesh_,
210             field.dimensions(),
211             internalField,
212             patchFields
213         )
214     );
218 template<class GeoField>
219 void fvFieldDecomposer::decomposeFields
221     const PtrList<GeoField>& fields
222 ) const
224     forAll (fields, fieldI)
225     {
226         decomposeField(fields[fieldI])().write();
227     }
231 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
233 } // End namespace Foam
235 // ************************************************************************* //