1 /*---------------------------------------------------------------------------*\
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
5 \\ / A nd | Copyright held by original author
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 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
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 "faFieldDecomposer.H"
28 #include "processorFaPatchField.H"
29 #include "processorFaePatchField.H"
31 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
36 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
39 tmp<GeometricField<Type, faPatchField, areaMesh> >
40 faFieldDecomposer::decomposeField
42 const GeometricField<Type, faPatchField, areaMesh>& field
45 // Create and map the internal field values
46 Field<Type> internalField(field.internalField(), faceAddressing_);
48 // Create and map the patch field values
49 PtrList<faPatchField<Type> > patchFields(boundaryAddressing_.size());
51 forAll (boundaryAddressing_, patchi)
53 if (boundaryAddressing_[patchi] >= 0)
58 faPatchField<Type>::New
60 field.boundaryField()[boundaryAddressing_[patchi]],
61 procMesh_.boundary()[patchi],
62 DimensionedField<Type, areaMesh>::null(),
63 *patchFieldDecomposerPtrs_[patchi]
72 new processorFaPatchField<Type>
74 procMesh_.boundary()[patchi],
75 DimensionedField<Type, areaMesh>::null(),
78 field.internalField(),
79 *processorAreaPatchFieldDecomposerPtrs_[patchi]
86 // Create the field for the processor
87 return tmp<GeometricField<Type, faPatchField, areaMesh> >
89 new GeometricField<Type, faPatchField, areaMesh>
94 procMesh_.time().timeName(),
109 tmp<GeometricField<Type, faePatchField, edgeMesh> >
110 faFieldDecomposer::decomposeField
112 const GeometricField<Type, faePatchField, edgeMesh>& field
120 procMesh_.nInternalEdges()
128 // Create and map the internal field values
129 Field<Type> internalField
131 field.internalField(),
135 // Problem with addressing when a processor patch picks up both internal
136 // edges and edges 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 edgeFields correspond to the one of edges in polyMesh
139 // (i.e. using slices)
140 Field<Type> allEdgeField(field.mesh().nEdges());
142 forAll (field.internalField(), i)
144 allEdgeField[i] = field.internalField()[i];
147 forAll (field.boundaryField(), patchi)
149 const Field<Type> & p = field.boundaryField()[patchi];
151 const label patchStart = field.mesh().boundary()[patchi].start();
155 allEdgeField[patchStart + i] = p[i];
159 // Create and map the patch field values
160 PtrList<faePatchField<Type> > patchFields(boundaryAddressing_.size());
162 forAll (boundaryAddressing_, patchi)
164 if (boundaryAddressing_[patchi] >= 0)
169 faePatchField<Type>::New
171 field.boundaryField()[boundaryAddressing_[patchi]],
172 procMesh_.boundary()[patchi],
173 DimensionedField<Type, edgeMesh>::null(),
174 *patchFieldDecomposerPtrs_[patchi]
183 new processorFaePatchField<Type>
185 procMesh_.boundary()[patchi],
186 DimensionedField<Type, edgeMesh>::null(),
190 *processorEdgePatchFieldDecomposerPtrs_[patchi]
197 // Create the field for the processor
198 return tmp<GeometricField<Type, faePatchField, edgeMesh> >
200 new GeometricField<Type, faePatchField, edgeMesh>
205 procMesh_.time().timeName(),
219 template<class GeoField>
220 void faFieldDecomposer::decomposeFields
222 const PtrList<GeoField>& fields
225 forAll (fields, fieldI)
227 decomposeField(fields[fieldI])().write();
232 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
234 } // End namespace Foam
236 // ************************************************************************* //