lagrangianFieldDecomposerDecomposeFields.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-2015 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,
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>, 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  (
85  Type>::typeName
86  )
87  );
88 
89  lagrangianFields.set
90  (
91  cloudI,
93  (
94  lagrangianTypeObjectsA.size() + lagrangianTypeObjectsB.size()
95  )
96  );
97 
98  label lagrangianFieldi = 0;
99 
100  forAllIter(IOobjectList, lagrangianTypeObjectsA, iter)
101  {
102  lagrangianFields[cloudI].set
103  (
104  lagrangianFieldi++,
105  new CompactIOField<Field<Type>, Type>(*iter())
106  );
107  }
108 
109  forAllIter(IOobjectList, lagrangianTypeObjectsB, iter)
110  {
111  lagrangianFields[cloudI].set
112  (
113  lagrangianFieldi++,
114  new CompactIOField<Field<Type>, Type>(*iter())
115  );
116  }
117 }
118 
119 
120 template<class Type>
123 (
124  const word& cloudName,
125  const IOField<Type>& field
126 ) const
127 {
128  // Create and map the internal field values
129  Field<Type> procField(field, particleIndices_);
130 
131  // Create the field for the processor
132  return tmp<IOField<Type> >
133  (
134  new IOField<Type>
135  (
136  IOobject
137  (
138  field.name(),
139  procMesh_.time().timeName(),
140  cloud::prefix/cloudName,
141  procMesh_,
142  IOobject::NO_READ,
143  IOobject::NO_WRITE,
144  false
145  ),
146  procField
147  )
148  );
149 }
150 
151 
152 template<class Type>
155 (
156  const word& cloudName,
157  const CompactIOField<Field<Type>, Type>& field
158 ) const
159 {
160  // Create and map the internal field values
161  Field<Field<Type> > procField(field, particleIndices_);
162 
163  // Create the field for the processor
164  return tmp<CompactIOField<Field<Type>, Type> >
165  (
166  new CompactIOField<Field<Type>, Type>
167  (
168  IOobject
169  (
170  field.name(),
171  procMesh_.time().timeName(),
172  cloud::prefix/cloudName,
173  procMesh_,
174  IOobject::NO_READ,
175  IOobject::NO_WRITE,
176  false
177  ),
178  procField
179  )
180  );
181 }
182 
183 
184 template<class GeoField>
186 (
187  const word& cloudName,
189 ) const
190 {
191  if (particleIndices_.size())
192  {
193  forAll(fields, fieldI)
194  {
195  decomposeField(cloudName, fields[fieldI])().write();
196  }
197  }
198 }
199 
200 
201 template<class GeoField>
203 (
204  const word& cloudName,
206 ) const
207 {
208  if (particleIndices_.size())
209  {
210  forAll(fields, fieldI)
211  {
212  decomposeFieldField(cloudName, fields[fieldI])().write();
213  }
214  }
215 }
216 
217 
218 // ************************************************************************* //
Foam::IOobject
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:91
Foam::lagrangianFieldDecomposer::decomposeField
tmp< IOField< Type > > decomposeField(const word &cloudName, const IOField< Type > &field) const
Decompose volume field.
Foam::word
A class for handling words, derived from string.
Definition: word.H:59
forAllIter
#define forAllIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:431
Foam::IOField
A primitive field of type <T> with automated input and output.
Definition: IOField.H:50
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:406
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:118
fields
Info<< "Creating field dpdt\n"<< endl;volScalarField dpdt(IOobject("dpdt", runTime.timeName(), mesh), mesh, dimensionedScalar("dpdt", p.dimensions()/dimTime, 0));Info<< "Creating field kinetic energy K\n"<< endl;volScalarField K("K", 0.5 *magSqr(U));volScalarField p_rgh(IOobject("p_rgh", runTime.timeName(), mesh, IOobject::MUST_READ, IOobject::AUTO_WRITE), mesh);p_rgh=p - rho *gh;mesh.setFluxRequired(p_rgh.name());multivariateSurfaceInterpolationScheme< scalar >::fieldTable fields
Definition: createFields.H:127
Foam::lagrangianFieldDecomposer::decomposeFieldFields
void decomposeFieldFields(const word &cloudName, const PtrList< GeoField > &fields) const
Definition: lagrangianFieldDecomposerDecomposeFields.C:203
IOobjectList.H
Foam::CompactIOField
A Field of objects of type <T> with automated input and output using a compact storage....
Definition: CompactIOField.H:50
Foam::lagrangianFieldDecomposer::decomposeFields
void decomposeFields(const word &cloudName, const PtrList< GeoField > &fields) const
Definition: lagrangianFieldDecomposerDecomposeFields.C:186
Foam::label
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
Foam::Field< Type >
Foam::lagrangianFieldDecomposer::readFieldFields
static void readFieldFields(const label cloudI, const IOobjectList &lagrangianObjects, PtrList< PtrList< CompactIOField< Field< Type >, Type > > > &lagrangianFields)
Definition: lagrangianFieldDecomposerDecomposeFields.C:68
Foam::IOobject::name
const word & name() const
Return name.
Definition: IOobject.H:273
Foam::PtrList
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: List.H:61
Foam::IOobjectList::lookupClass
IOobjectList lookupClass(const word &className) const
Return the list for all IOobjects of a given class.
Definition: IOobjectList.C:199
Foam::IOobjectList
List of IOobjects with searching and retrieving facilities.
Definition: IOobjectList.H:50
Foam::lagrangianFieldDecomposer::readFields
static void readFields(const label cloudI, const IOobjectList &lagrangianObjects, PtrList< PtrList< IOField< Type > > > &lagrangianFields)
Definition: lagrangianFieldDecomposerDecomposeFields.C:33
cloudName
const word cloudName(propsDict.lookup("cloudName"))
write
Tcoeff write()
lagrangianFieldDecomposer.H