ENH: decomposePar.C: add shortcircuit to avoid allocating point mappers
[OpenFOAM-2.0.x.git] / applications / utilities / parallelProcessing / decomposePar / fvFieldDecomposerDecomposeFields.C
blob0bce0023cc8b3c23e94e2213a79464351921ed08
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 2011 OpenFOAM Foundation
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
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
19     for more details.
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"
27 #include "processorFvPatchField.H"
28 #include "processorFvsPatchField.H"
29 #include "processorCyclicFvPatchField.H"
30 #include "processorCyclicFvsPatchField.H"
32 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
34 template<class Type>
35 Foam::tmp<Foam::GeometricField<Type, Foam::fvPatchField, Foam::volMesh> >
36 Foam::fvFieldDecomposer::decomposeField
38     const GeometricField<Type, fvPatchField, volMesh>& field
39 ) const
41     // Create and map the internal field values
42     Field<Type> internalField(field.internalField(), cellAddressing_);
44     // Create and map the patch field values
45     PtrList<fvPatchField<Type> > patchFields(boundaryAddressing_.size());
47     forAll(boundaryAddressing_, patchi)
48     {
49         if (patchFieldDecomposerPtrs_[patchi])
50         {
51             patchFields.set
52             (
53                 patchi,
54                 fvPatchField<Type>::New
55                 (
56                     field.boundaryField()[boundaryAddressing_[patchi]],
57                     procMesh_.boundary()[patchi],
58                     DimensionedField<Type, volMesh>::null(),
59                     *patchFieldDecomposerPtrs_[patchi]
60                 )
61             );
62         }
63         else if (isA<processorCyclicFvPatch>(procMesh_.boundary()[patchi]))
64         {
65             patchFields.set
66             (
67                 patchi,
68                 new processorCyclicFvPatchField<Type>
69                 (
70                     procMesh_.boundary()[patchi],
71                     DimensionedField<Type, volMesh>::null(),
72                     Field<Type>
73                     (
74                         field.internalField(),
75                         *processorVolPatchFieldDecomposerPtrs_[patchi]
76                     )
77                 )
78             );
79         }
80         else if (isA<processorFvPatch>(procMesh_.boundary()[patchi]))
81         {
82             patchFields.set
83             (
84                 patchi,
85                 new processorFvPatchField<Type>
86                 (
87                     procMesh_.boundary()[patchi],
88                     DimensionedField<Type, volMesh>::null(),
89                     Field<Type>
90                     (
91                         field.internalField(),
92                         *processorVolPatchFieldDecomposerPtrs_[patchi]
93                     )
94                 )
95             );
96         }
97         else
98         {
99             FatalErrorIn("fvFieldDecomposer::decomposeField()")
100                 << "Unknown type." << abort(FatalError);
101         }
102     }
104     // Create the field for the processor
105     return tmp<GeometricField<Type, fvPatchField, volMesh> >
106     (
107         new GeometricField<Type, fvPatchField, volMesh>
108         (
109             IOobject
110             (
111                 field.name(),
112                 procMesh_.time().timeName(),
113                 procMesh_,
114                 IOobject::NO_READ,
115                 IOobject::NO_WRITE
116             ),
117             procMesh_,
118             field.dimensions(),
119             internalField,
120             patchFields
121         )
122     );
126 template<class Type>
127 Foam::tmp<Foam::GeometricField<Type, Foam::fvsPatchField, Foam::surfaceMesh> >
128 Foam::fvFieldDecomposer::decomposeField
130     const GeometricField<Type, fvsPatchField, surfaceMesh>& field
131 ) const
133     labelList mapAddr
134     (
135         labelList::subList
136         (
137             faceAddressing_,
138             procMesh_.nInternalFaces()
139         )
140     );
141     forAll(mapAddr, i)
142     {
143         mapAddr[i] -= 1;
144     }
146     // Create and map the internal field values
147     Field<Type> internalField
148     (
149         field.internalField(),
150         mapAddr
151     );
153     // Problem with addressing when a processor patch picks up both internal
154     // faces and faces from cyclic boundaries. This is a bit of a hack, but
155     // I cannot find a better solution without making the internal storage
156     // mechanism for surfaceFields correspond to the one of faces in polyMesh
157     // (i.e. using slices)
158     Field<Type> allFaceField(field.mesh().nFaces());
160     forAll(field.internalField(), i)
161     {
162         allFaceField[i] = field.internalField()[i];
163     }
165     forAll(field.boundaryField(), patchi)
166     {
167         const Field<Type> & p = field.boundaryField()[patchi];
169         const label patchStart = field.mesh().boundaryMesh()[patchi].start();
171         forAll(p, i)
172         {
173             allFaceField[patchStart + i] = p[i];
174         }
175     }
177     // Create and map the patch field values
178     PtrList<fvsPatchField<Type> > patchFields(boundaryAddressing_.size());
180     forAll(boundaryAddressing_, patchi)
181     {
182         if (patchFieldDecomposerPtrs_[patchi])
183         {
184             patchFields.set
185             (
186                 patchi,
187                 fvsPatchField<Type>::New
188                 (
189                     field.boundaryField()[boundaryAddressing_[patchi]],
190                     procMesh_.boundary()[patchi],
191                     DimensionedField<Type, surfaceMesh>::null(),
192                     *patchFieldDecomposerPtrs_[patchi]
193                 )
194             );
195         }
196         else if (isA<processorCyclicFvPatch>(procMesh_.boundary()[patchi]))
197         {
198             patchFields.set
199             (
200                 patchi,
201                 new processorCyclicFvsPatchField<Type>
202                 (
203                     procMesh_.boundary()[patchi],
204                     DimensionedField<Type, surfaceMesh>::null(),
205                     Field<Type>
206                     (
207                         allFaceField,
208                         *processorSurfacePatchFieldDecomposerPtrs_[patchi]
209                     )
210                 )
211             );
212         }
213         else if (isA<processorFvPatch>(procMesh_.boundary()[patchi]))
214         {
215             patchFields.set
216             (
217                 patchi,
218                 new processorFvsPatchField<Type>
219                 (
220                     procMesh_.boundary()[patchi],
221                     DimensionedField<Type, surfaceMesh>::null(),
222                     Field<Type>
223                     (
224                         allFaceField,
225                         *processorSurfacePatchFieldDecomposerPtrs_[patchi]
226                     )
227                 )
228             );
229         }
230         else
231         {
232             FatalErrorIn("fvFieldDecomposer::decomposeField()")
233                 << "Unknown type." << abort(FatalError);
234         }
235     }
237     // Create the field for the processor
238     return tmp<GeometricField<Type, fvsPatchField, surfaceMesh> >
239     (
240         new GeometricField<Type, fvsPatchField, surfaceMesh>
241         (
242             IOobject
243             (
244                 field.name(),
245                 procMesh_.time().timeName(),
246                 procMesh_,
247                 IOobject::NO_READ,
248                 IOobject::NO_WRITE
249             ),
250             procMesh_,
251             field.dimensions(),
252             internalField,
253             patchFields
254         )
255     );
259 template<class GeoField>
260 void Foam::fvFieldDecomposer::decomposeFields
262     const PtrList<GeoField>& fields
263 ) const
265     forAll(fields, fieldI)
266     {
267         decomposeField(fields[fieldI])().write();
268     }
272 // ************************************************************************* //