parLagrangianRedistributorRedistributeFields.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) 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 
26 #include "parLagrangianRedistributor.H"
27 #include "Time.H"
28 #include "IOobjectList.H"
29 #include "mapDistributePolyMesh.H"
30 #include "cloud.H"
31 #include "CompactIOField.H"
32 #include "passiveParticleCloud.H"
33 
34 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
35 
36 template<class Container>
38 (
39  const IOobjectList& objects,
40  const HashSet<word>& selectedFields
41 )
42 {
43  const word fieldClassName(Container::typeName);
44 
45  // Parallel synchronise
46  wordList fieldNames(objects.names(fieldClassName));
47  Pstream::combineGather(fieldNames, ListUniqueEqOp<word>());
49 
50  if (!selectedFields.empty())
51  {
52  DynamicList<word> selectedNames(fieldNames.size());
53  forAll(fieldNames, i)
54  {
55  if (selectedFields.found(fieldNames[i]))
56  {
57  selectedNames.append(fieldNames[i]);
58  }
59  }
60  fieldNames.transfer(selectedNames);
61  }
62  return fieldNames;
63 }
64 
65 
66 template<class Type>
68 (
69  const mapDistributeBase& map,
70  const word& cloudName,
71  const IOobjectList& objects,
72  const HashSet<word>& selectedFields
73 ) const
74 {
75  const wordList objectNames
76  (
77  filterObjects<IOField<Type> >
78  (
79  objects,
80  selectedFields
81  )
82  );
83 
84  if (objectNames.size())
85  {
86  const word fieldClassName(IOField<Type>::typeName);
87 
88  Info<< " Redistributing lagrangian "
89  << fieldClassName << "s\n" << endl;
90 
91  forAll(objectNames, i)
92  {
93  Info<< " " << objectNames[i] << endl;
94 
95  // Read if present
96  IOField<Type> field
97  (
98  IOobject
99  (
100  objectNames[i],
101  srcMesh_.time().timeName(),
103  srcMesh_,
106  false
107  ),
108  0
109  );
110 
111  map.distribute(field);
112 
113 
114  if (field.size())
115  {
116  IOField<Type>
117  (
118  IOobject
119  (
120  objectNames[i],
121  tgtMesh_.time().timeName(),
123  tgtMesh_,
126  false
127  ),
128  xferMove<Field<Type> >(field)
129  ).write();
130  }
131  }
132 
133  Info<< endl;
134  }
135 }
136 
137 
138 template<class Type>
140 (
141  const mapDistributeBase& map,
142  const word& cloudName,
143  const IOobjectList& objects,
144  const HashSet<word>& selectedFields
145 ) const
146 {
147  wordList objectNames
148  (
149  filterObjects<CompactIOField<Field<Type>, Type> >
150  (
151  objects,
152  selectedFields
153  )
154  );
155 
156  // Append IOField names
157  {
158  const wordList ioFieldNames
159  (
160  filterObjects<IOField<Field<Type> > >
161  (
162  objects,
163  selectedFields
164  )
165  );
166  objectNames.append(ioFieldNames);
167  }
168 
169 
170  if (objectNames.size())
171  {
172  const word fieldClassName(CompactIOField<Field<Type>, Type>::typeName);
173 
174  Info<< " Redistributing lagrangian "
175  << fieldClassName << "s\n" << endl;
176 
177  forAll(objectNames, i)
178  {
179  Info<< " " << objectNames[i] << endl;
180 
181  // Read if present
182  CompactIOField<Field<Type>, Type > field
183  (
184  IOobject
185  (
186  objectNames[i],
187  srcMesh_.time().timeName(),
189  srcMesh_,
192  false
193  ),
194  0
195  );
196 
197  // Distribute
198  map.distribute(field);
199 
200  // Write
201  if (field.size())
202  {
203  CompactIOField<Field<Type>, Type>
204  (
205  IOobject
206  (
207  objectNames[i],
208  tgtMesh_.time().timeName(),
210  tgtMesh_,
213  false
214  ),
215  xferMove<Field<Field<Type> > >(field)
216  ).write();
217  }
218  }
219  }
220 }
221 
222 
223 template<class Container>
225 (
226  const passiveParticleCloud& cloud,
227  const IOobjectList& objects,
228  const HashSet<word>& selectedFields
229 )
230 {
231  const wordList objectNames
232  (
233  filterObjects<Container>
234  (
235  objects,
236  selectedFields
237  )
238  );
239 
240  if (objectNames.size())
241  {
242  const word fieldClassName(Container::typeName);
243 
244  Info<< " Reading lagrangian "
245  << fieldClassName << "s\n" << endl;
246 
247  forAll(objectNames, i)
248  {
249  Info<< " " << objectNames[i] << endl;
250 
251  // Read if present
252  Container* fieldPtr = new Container
253  (
254  IOobject
255  (
256  objectNames[i],
257  cloud.time().timeName(),
258  cloud,
261  ),
262  0
263  );
264 
265  fieldPtr->store();
266  }
267  }
268 }
269 
270 
271 template<class Container>
273 (
274  const mapDistributeBase& map,
275  passiveParticleCloud& cloud
276 ) const
277 {
278  HashTable<Container*> fields
279  (
280  cloud.lookupClass<Container >()
281  );
282 
283  if (fields.size())
284  {
285  const word fieldClassName(Container::typeName);
286 
287  Info<< " Redistributing lagrangian "
288  << fieldClassName << "s\n" << endl;
289 
290  forAllIter(typename HashTable<Container*>, fields, iter)
291  {
292  Container& field = *iter();
293 
294  Info<< " " << field.name() << endl;
295 
296  map.distribute(field);
297 
298  if (field.size())
299  {
300  Container
301  (
302  IOobject
303  (
304  field.name(),
305  tgtMesh_.time().timeName(),
306  cloud::prefix/cloud.name(),
307  tgtMesh_,
310  false
311  ),
312  xferMove<Field<typename Container::value_type> >(field)
313  ).write();
314  }
315  }
316  }
317 }
318 
319 
320 // ************************************************************************* //
Foam::parLagrangianRedistributor::redistributeStoredLagrangianFields
void redistributeStoredLagrangianFields(const mapDistributeBase &map, passiveParticleCloud &cloud) const
Redistribute and write stored lagrangian fields.
Definition: parLagrangianRedistributorRedistributeFields.C:273
Foam::parLagrangianRedistributor::readLagrangianFields
static void readLagrangianFields(const passiveParticleCloud &cloud, const IOobjectList &objects, const HashSet< word > &selectedFields)
Read and store all fields of a cloud.
Definition: parLagrangianRedistributorRedistributeFields.C:225
forAllIter
#define forAllIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:431
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:406
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::Field< Type >::typeName
static const char *const typeName
Definition: Field.H:94
Foam::List::transfer
void transfer(List< T > &)
Transfer the contents of the argument List into this list.
cloud.H
IOobjectList.H
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
Foam::IOobject::NO_WRITE
@ NO_WRITE
Definition: IOobject.H:118
Foam::parLagrangianRedistributor::redistributeLagrangianFields
void redistributeLagrangianFields(const mapDistributeBase &map, const word &cloudName, const IOobjectList &objects, const HashSet< word > &selectedFields) const
Read, redistribute and write all/selected lagrangian fields.
Definition: parLagrangianRedistributorRedistributeFields.C:68
Foam::wordList
List< word > wordList
A List of words.
Definition: fileName.H:54
mapDistributePolyMesh.H
Foam::Ostream::write
virtual Ostream & write(const token &)=0
Write next token to stream.
Foam::IOobject::NO_READ
@ NO_READ
Definition: IOobject.H:111
Foam::Info
messageStream Info
Foam::parLagrangianRedistributor::filterObjects
static wordList filterObjects(const IOobjectList &objects, const HashSet< word > &selectedFields)
Pick up any fields of a given type.
Foam::Pstream::combineGather
static void combineGather(const List< commsStruct > &comms, T &Value, const CombineOp &cop, const int tag, const label comm)
Definition: combineGatherScatter.C:50
Foam::cloud::prefix
static const word prefix
The prefix to local: lagrangian.
Definition: cloud.H:71
Foam::xferMove
Xfer< T > xferMove(T &)
Construct by transferring the contents of the arg.
passiveParticleCloud.H
Foam::parLagrangianRedistributor::redistributeLagrangianFieldFields
void redistributeLagrangianFieldFields(const mapDistributeBase &map, const word &cloudName, const IOobjectList &objects, const HashSet< word > &selectedFields) const
Read, redistribute and write all/selected lagrangian fieldFields.
Definition: parLagrangianRedistributorRedistributeFields.C:140
Foam::fieldNames
wordList fieldNames(const IOobjectList &objects, const bool syncPar)
Get sorted names of fields of type. If syncPar and running in parallel.
Definition: ReadFields.C:36
Foam::List
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: HashTable.H:59
cloudName
const word cloudName(propsDict.lookup("cloudName"))
Foam::List::size
void size(const label)
Override size to be inconsistent with allocated storage.
CompactIOField.H
Foam::IOobject::READ_IF_PRESENT
@ READ_IF_PRESENT
Definition: IOobject.H:110
Foam::Pstream::combineScatter
static void combineScatter(const List< commsStruct > &comms, T &Value, const int tag, const label comm)
Scatter data. Reverse of combineGather.
Definition: combineGatherScatter.C:178