lagrangianFieldDecomposerFields.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 | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2017 OpenFOAM Foundation
9 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 \*---------------------------------------------------------------------------*/
27 
29 #include "IOobjectList.H"
30 
31 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
32 
33 template<class Type>
35 (
36  const label cloudI,
37  const IOobjectList& lagrangianObjects,
38  PtrList<PtrList<IOField<Type>>>& lagrangianFields
39 )
40 {
41  // Search list of objects for lagrangian fields
42  IOobjectList lagrangianTypeObjects
43  (
44  lagrangianObjects.lookupClass(IOField<Type>::typeName)
45  );
46 
47  lagrangianFields.set
48  (
49  cloudI,
50  new PtrList<IOField<Type>>
51  (
52  lagrangianTypeObjects.size()
53  )
54  );
55 
56  label lagrangianFieldi = 0;
57  forAllConstIters(lagrangianTypeObjects, iter)
58  {
59  lagrangianFields[cloudI].set
60  (
61  lagrangianFieldi++,
62  new IOField<Type>(*iter())
63  );
64  }
65 }
66 
67 
68 template<class Type>
70 (
71  const label cloudI,
72  const IOobjectList& lagrangianObjects,
73  PtrList<PtrList<CompactIOField<Field<Type>, Type>>>& lagrangianFields
74 )
75 {
76  // Search list of objects for lagrangian fields
77  IOobjectList lagrangianTypeObjectsA
78  (
79  lagrangianObjects.lookupClass(IOField<Field<Type>>::typeName)
80  );
81 
82  IOobjectList lagrangianTypeObjectsB
83  (
84  lagrangianObjects.lookupClass
85  (
86  CompactIOField<Field<Type>,
87  Type>::typeName
88  )
89  );
90 
91  lagrangianFields.set
92  (
93  cloudI,
94  new PtrList<CompactIOField<Field<Type>, Type>>
95  (
96  lagrangianTypeObjectsA.size() + lagrangianTypeObjectsB.size()
97  )
98  );
99 
100  label lagrangianFieldi = 0;
101 
102  forAllConstIters(lagrangianTypeObjectsA, iter)
103  {
104  lagrangianFields[cloudI].set
105  (
106  lagrangianFieldi++,
107  new CompactIOField<Field<Type>, Type>(*iter())
108  );
109  }
110 
111  forAllConstIters(lagrangianTypeObjectsB, iter)
112  {
113  lagrangianFields[cloudI].set
114  (
115  lagrangianFieldi++,
116  new CompactIOField<Field<Type>, Type>(*iter())
117  );
118  }
119 }
120 
121 
122 template<class Type>
125 (
126  const word& cloudName,
127  const IOField<Type>& field
128 ) const
129 {
130  // Create and map the internal field values
131  Field<Type> procField(field, particleIndices_);
132 
133  // Create the field for the processor
134  return tmp<IOField<Type>>::New
135  (
136  IOobject
137  (
138  field.name(),
139  procMesh_.time().timeName(),
141  procMesh_,
144  false
145  ),
146  procField
147  );
148 }
149 
150 
151 template<class Type>
154 (
155  const word& cloudName,
156  const CompactIOField<Field<Type>, Type>& field
157 ) const
158 {
159  // Create and map the internal field values
160  Field<Field<Type>> procField(field, particleIndices_);
161 
162  // Create the field for the processor
163  return tmp<CompactIOField<Field<Type>, Type>>::New
164  (
165  IOobject
166  (
167  field.name(),
168  procMesh_.time().timeName(),
170  procMesh_,
173  false
174  ),
175  procField
176  );
177 }
178 
179 
180 template<class GeoField>
182 (
183  const word& cloudName,
184  const PtrList<GeoField>& fields
185 ) const
186 {
187  //if (particleIndices_.size())
188  {
189  bool valid = particleIndices_.size() > 0;
190  forAll(fields, fieldi)
191  {
192  decomposeField(cloudName, fields[fieldi])().write(valid);
193  }
194  }
195 }
196 
197 
198 template<class GeoField>
200 (
201  const word& cloudName,
202  const PtrList<GeoField>& fields
203 ) const
204 {
205  //if (particleIndices_.size())
206  {
207  bool valid = particleIndices_.size() > 0;
208  forAll(fields, fieldi)
209  {
210  decomposeFieldField(cloudName, fields[fieldi])().write(valid);
211  }
212  }
213 }
214 
215 
216 // ************************************************************************* //
Foam::IOobject::NO_WRITE
@ NO_WRITE
Definition: IOobject.H:191
Foam::cloud::prefix
static const word prefix
Definition: cloud.H:83
Foam::lagrangianFieldDecomposer::decomposeField
tmp< IOField< Type > > decomposeField(const word &cloudName, const IOField< Type > &field) const
cloudName
const word cloudName(propsDict.get< word >("cloud"))
Foam::lagrangianFieldDecomposer::decomposeFieldField
tmp< CompactIOField< Field< Type >, Type > > decomposeFieldField(const word &cloudName, const CompactIOField< Field< Type >, Type > &field) const
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:57
Foam::lagrangianFieldDecomposer::decomposeFieldFields
void decomposeFieldFields(const word &cloudName, const PtrList< GeoField > &fields) const
IOobjectList.H
Foam::lagrangianFieldDecomposer::readFields
static void readFields(const label cloudI, const IOobjectList &lagrangianObjects, PtrList< PtrList< IOField< Type >>> &lagrangianFields)
Foam::lagrangianFieldDecomposer::decomposeFields
void decomposeFields(const word &cloudName, const PtrList< GeoField > &fields) const
forAll
#define forAll(list, i)
Definition: stdFoam.H:349
field
rDeltaTY field()
Foam::New
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
Definition: DimensionedFieldReuseFunctions.H:100
forAllConstIters
forAllConstIters(mixture.phases(), phase)
Definition: pEqn.H:28
Foam::lagrangianFieldDecomposer::readFieldFields
static void readFieldFields(const label cloudI, const IOobjectList &lagrangianObjects, PtrList< PtrList< CompactIOField< Field< Type >, Type >> > &lagrangianFields)
Foam::vtk::write
void write(vtk::formatter &fmt, const Type &val, const label n=1)
Definition: foamVtkOutputTemplates.C:29
Foam::IOobject::NO_READ
@ NO_READ
Definition: IOobject.H:184
lagrangianFieldDecomposer.H
fields
multivariateSurfaceInterpolationScheme< scalar >::fieldTable fields
Definition: createFields.H:97