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"
27 #include "processorFvPatchField.H"
28 #include "processorFvsPatchField.H"
29 #include "processorCyclicFvPatchField.H"
30 #include "processorCyclicFvsPatchField.H"
32 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
35 Foam::tmp<Foam::GeometricField<Type, Foam::fvPatchField, Foam::volMesh> >
36 Foam::fvFieldDecomposer::decomposeField
38 const GeometricField<Type, fvPatchField, volMesh>& field
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)
49 if (patchFieldDecomposerPtrs_[patchi])
54 fvPatchField<Type>::New
56 field.boundaryField()[boundaryAddressing_[patchi]],
57 procMesh_.boundary()[patchi],
58 DimensionedField<Type, volMesh>::null(),
59 *patchFieldDecomposerPtrs_[patchi]
63 else if (isA<processorCyclicFvPatch>(procMesh_.boundary()[patchi]))
68 new processorCyclicFvPatchField<Type>
70 procMesh_.boundary()[patchi],
71 DimensionedField<Type, volMesh>::null(),
74 field.internalField(),
75 *processorVolPatchFieldDecomposerPtrs_[patchi]
80 else if (isA<processorFvPatch>(procMesh_.boundary()[patchi]))
85 new processorFvPatchField<Type>
87 procMesh_.boundary()[patchi],
88 DimensionedField<Type, volMesh>::null(),
91 field.internalField(),
92 *processorVolPatchFieldDecomposerPtrs_[patchi]
99 FatalErrorIn("fvFieldDecomposer::decomposeField()")
100 << "Unknown type." << abort(FatalError);
104 // Create the field for the processor
105 return tmp<GeometricField<Type, fvPatchField, volMesh> >
107 new GeometricField<Type, fvPatchField, volMesh>
112 procMesh_.time().timeName(),
127 Foam::tmp<Foam::GeometricField<Type, Foam::fvsPatchField, Foam::surfaceMesh> >
128 Foam::fvFieldDecomposer::decomposeField
130 const GeometricField<Type, fvsPatchField, surfaceMesh>& field
138 procMesh_.nInternalFaces()
146 // Create and map the internal field values
147 Field<Type> internalField
149 field.internalField(),
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)
162 allFaceField[i] = field.internalField()[i];
165 forAll(field.boundaryField(), patchi)
167 const Field<Type> & p = field.boundaryField()[patchi];
169 const label patchStart = field.mesh().boundaryMesh()[patchi].start();
173 allFaceField[patchStart + i] = p[i];
177 // Create and map the patch field values
178 PtrList<fvsPatchField<Type> > patchFields(boundaryAddressing_.size());
180 forAll(boundaryAddressing_, patchi)
182 if (patchFieldDecomposerPtrs_[patchi])
187 fvsPatchField<Type>::New
189 field.boundaryField()[boundaryAddressing_[patchi]],
190 procMesh_.boundary()[patchi],
191 DimensionedField<Type, surfaceMesh>::null(),
192 *patchFieldDecomposerPtrs_[patchi]
196 else if (isA<processorCyclicFvPatch>(procMesh_.boundary()[patchi]))
201 new processorCyclicFvsPatchField<Type>
203 procMesh_.boundary()[patchi],
204 DimensionedField<Type, surfaceMesh>::null(),
208 *processorSurfacePatchFieldDecomposerPtrs_[patchi]
213 else if (isA<processorFvPatch>(procMesh_.boundary()[patchi]))
218 new processorFvsPatchField<Type>
220 procMesh_.boundary()[patchi],
221 DimensionedField<Type, surfaceMesh>::null(),
225 *processorSurfacePatchFieldDecomposerPtrs_[patchi]
232 FatalErrorIn("fvFieldDecomposer::decomposeField()")
233 << "Unknown type." << abort(FatalError);
237 // Create the field for the processor
238 return tmp<GeometricField<Type, fvsPatchField, surfaceMesh> >
240 new GeometricField<Type, fvsPatchField, surfaceMesh>
245 procMesh_.time().timeName(),
259 template<class GeoField>
260 void Foam::fvFieldDecomposer::decomposeFields
262 const PtrList<GeoField>& fields
265 forAll(fields, fieldI)
267 decomposeField(fields[fieldI])().write();
272 // ************************************************************************* //