pointFieldDecomposerDecomposeFields.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | Copyright (C) 2011-2016 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
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.
15 
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.
20 
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/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 #include "pointFieldDecomposer.H"
28 
29 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
30 
31 template<class Type>
34 (
35  const GeometricField<Type, pointPatchField, pointMesh>& field
36 ) const
37 {
38  // Create and map the internal field values
39  Field<Type> internalField(field.primitiveField(), pointAddressing_);
40 
41  // Create a list of pointers for the patchFields
42  PtrList<pointPatchField<Type>> patchFields(boundaryAddressing_.size());
43 
44  // Create and map the patch field values
45  forAll(boundaryAddressing_, patchi)
46  {
47  if (patchFieldDecomposerPtrs_[patchi])
48  {
49  patchFields.set
50  (
51  patchi,
53  (
54  field.boundaryField()[boundaryAddressing_[patchi]],
55  procMesh_.boundary()[patchi],
57  *patchFieldDecomposerPtrs_[patchi]
58  )
59  );
60  }
61  else
62  {
63  patchFields.set
64  (
65  patchi,
66  new processorPointPatchField<Type>
67  (
68  procMesh_.boundary()[patchi],
70  )
71  );
72  }
73  }
74 
75  // Create the field for the processor
76  return tmp<GeometricField<Type, pointPatchField, pointMesh>>
77  (
78  new GeometricField<Type, pointPatchField, pointMesh>
79  (
80  IOobject
81  (
82  field.name(),
83  procMesh_().time().timeName(),
84  procMesh_(),
87  false
88  ),
89  procMesh_,
90  field.dimensions(),
91  internalField,
92  patchFields
93  )
94  );
95 }
96 
97 
98 template<class GeoField>
100 (
101  const PtrList<GeoField>& fields
102 ) const
103 {
104  forAll(fields, fieldi)
105  {
106  decomposeField(fields[fieldi])().write();
107  }
108 }
109 
110 
111 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:428
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:76
static word timeName(const scalar, const int precision=precision_)
Return time name of given scalar time.
Definition: Time.C:715
const pointBoundaryMesh & boundary() const
Return reference to boundary mesh.
Definition: pointMesh.H:106
static const DimensionedField< Type, GeoMesh > & null()
Return a null DimensionedField.
const Time & time() const
Return time.
Definition: IOobject.C:227
tmp< GeometricField< Type, pointPatchField, pointMesh > > decomposeField(const GeometricField< Type, pointPatchField, pointMesh > &) const
Decompose point field.
label patchi
static autoPtr< pointPatchField< Type > > New(const word &, const pointPatch &, const DimensionedField< Type, pointMesh > &)
Return a pointer to a new patchField created on freestore given.
A class for managing temporary objects.
Definition: PtrList.H:54
void decomposeFields(const PtrList< GeoField > &fields) const
runTime write()