lagrangianFieldDecomposerDecomposeFields.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2022 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 
27 #include "IOobjectList.H"
28 
29 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
30 
31 template<class Type>
33 (
34  const label cloudI,
35  const IOobjectList& lagrangianObjects,
36  PtrList<PtrList<IOField<Type>>>& lagrangianFields
37 )
38 {
39  // Search list of objects for lagrangian fields
40  IOobjectList lagrangianTypeObjects
41  (
42  lagrangianObjects.lookupClass(IOField<Type>::typeName)
43  );
44 
45  lagrangianFields.set
46  (
47  cloudI,
48  new PtrList<IOField<Type>>
49  (
50  lagrangianTypeObjects.size()
51  )
52  );
53 
54  label lagrangianFieldi = 0;
55  forAllIter(IOobjectList, lagrangianTypeObjects, iter)
56  {
57  lagrangianFields[cloudI].set
58  (
59  lagrangianFieldi++,
60  new IOField<Type>(*iter())
61  );
62  }
63 }
64 
65 
66 template<class Type>
68 (
69  const label cloudI,
70  const IOobjectList& lagrangianObjects,
71  PtrList<PtrList<CompactIOField<Field<Type>>>>& lagrangianFields
72 )
73 {
74  // Search list of objects for lagrangian fields
75  IOobjectList lagrangianTypeObjectsA
76  (
77  lagrangianObjects.lookupClass(IOField<Field<Type>>::typeName)
78  );
79 
80  IOobjectList lagrangianTypeObjectsB
81  (
82  lagrangianObjects.lookupClass
83  (
84  CompactIOField<Field<Type>>::typeName
85  )
86  );
87 
88  lagrangianFields.set
89  (
90  cloudI,
91  new PtrList<CompactIOField<Field<Type>>>
92  (
93  lagrangianTypeObjectsA.size() + lagrangianTypeObjectsB.size()
94  )
95  );
96 
97  label lagrangianFieldi = 0;
98 
99  forAllIter(IOobjectList, lagrangianTypeObjectsA, iter)
100  {
101  lagrangianFields[cloudI].set
102  (
103  lagrangianFieldi++,
104  new CompactIOField<Field<Type>>(*iter())
105  );
106  }
107 
108  forAllIter(IOobjectList, lagrangianTypeObjectsB, iter)
109  {
110  lagrangianFields[cloudI].set
111  (
112  lagrangianFieldi++,
113  new CompactIOField<Field<Type>>(*iter())
114  );
115  }
116 }
117 
118 
119 template<class Type>
122 (
123  const word& cloudName,
124  const IOField<Type>& field
125 ) const
126 {
127  // Create and map the internal field values
128  Field<Type> procField(field, particleIndices_);
129 
130  // Create the field for the processor
131  return tmp<IOField<Type>>
132  (
133  new IOField<Type>
134  (
135  IOobject
136  (
137  field.name(),
138  procMesh_.time().name(),
140  procMesh_,
143  false
144  ),
145  procField
146  )
147  );
148 }
149 
150 
151 template<class Type>
154 (
155  const word& cloudName,
156  const CompactIOField<Field<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>>>
164  (
165  new CompactIOField<Field<Type>>
166  (
167  IOobject
168  (
169  field.name(),
170  procMesh_.time().name(),
172  procMesh_,
175  false
176  ),
177  procField
178  )
179  );
180 }
181 
182 
183 template<class GeoField>
185 (
186  const word& cloudName,
187  const PtrList<GeoField>& fields
188 ) const
189 {
190  const bool write = particleIndices_.size() > 0;
191  forAll(fields, fieldi)
192  {
193  decomposeField(cloudName, fields[fieldi])().write(write);
194  }
195 }
196 
197 
198 template<class GeoField>
200 (
201  const word& cloudName,
202  const PtrList<GeoField>& fields
203 ) const
204 {
205  const bool write = particleIndices_.size() > 0;
206  forAll(fields, fieldi)
207  {
208  decomposeFieldField(cloudName, fields[fieldi])().write(write);
209  }
210 }
211 
212 
213 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
#define forAllIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:459
static const char *const typeName
Definition: Field.H:105
static const word prefix
The prefix to local: lagrangian.
Definition: cloud.H:63
static void readFields(const label cloudI, const IOobjectList &lagrangianObjects, PtrList< PtrList< IOField< Type >>> &lagrangianFields)
tmp< IOField< Type > > decomposeField(const word &cloudName, const IOField< Type > &field) const
Decompose volume field.
static void readFieldFields(const label cloudI, const IOobjectList &lagrangianObjects, PtrList< PtrList< CompactIOField< Field< Type >>>> &lagrangianFields)
tmp< CompactIOField< Field< Type > > > decomposeFieldField(const word &cloudName, const CompactIOField< Field< Type >> &field) const
void decomposeFields(const word &cloudName, const PtrList< GeoField > &fields) const
void decomposeFieldFields(const word &cloudName, const PtrList< GeoField > &fields) const
A class for managing temporary objects.
Definition: tmp.H:55
Info<< "Calculating turbulent flame speed field St\n"<< endl;volScalarField St(IOobject("St", runTime.name(), mesh, IOobject::NO_READ, IOobject::AUTO_WRITE), flameWrinkling->Xi() *Su);multivariateSurfaceInterpolationScheme< scalar >::fieldTable fields
Definition: createFields.H:230
void write(std::ostream &os, const bool binary, List< floatScalar > &fField)
Write floats ascii or binary.
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
const word cloudName(propsDict.lookup("cloudName"))