MapLagrangianFields.H
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 InNamespace
25  Foam
26 
27 Description
28  Gets the indices of (source)particles that have been appended to the
29  target cloud and maps the lagrangian fields accordingly.
30 
31 \*---------------------------------------------------------------------------*/
32 
33 #ifndef MapLagrangianFields_H
34 #define MapLagrangianFields_H
35 
36 #include "cloud.H"
37 #include "GeometricField.H"
38 #include "meshToMesh.H"
39 #include "IOobjectList.H"
40 #include "CompactIOField.H"
41 
42 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
43 
44 namespace Foam
45 {
46 
47 //- Gets the indices of (source)particles that have been appended to the
48 // target cloud and maps the lagrangian fields accordingly.
49 template<class Type>
51 (
52  const string& cloudName,
53  const IOobjectList& objects,
54  const polyMesh& meshTarget,
55  const labelList& addParticles
56 )
57 {
58  {
60 
61  forAllIter(IOobjectList, fields, fieldIter)
62  {
63  const word& fieldName = fieldIter()->name();
64 
65  Info<< " mapping lagrangian field " << fieldName << endl;
66 
67  // Read field (does not need mesh)
68  IOField<Type> fieldSource(*fieldIter());
69 
70  // Map
71  IOField<Type> fieldTarget
72  (
73  IOobject
74  (
75  fieldName,
76  meshTarget.time().timeName(),
78  meshTarget,
81  false
82  ),
83  addParticles.size()
84  );
85 
86  forAll(addParticles, i)
87  {
88  fieldTarget[i] = fieldSource[addParticles[i]];
89  }
90 
91  // Write field
92  fieldTarget.write();
93  }
94  }
95 
96  {
97  IOobjectList fieldFields =
98  objects.lookupClass(IOField<Field<Type> >::typeName);
99 
100  forAllIter(IOobjectList, fieldFields, fieldIter)
101  {
102  const word& fieldName = fieldIter()->name();
103 
104  Info<< " mapping lagrangian fieldField " << fieldName << endl;
105 
106  // Read field (does not need mesh)
107  // Note: some fieldFields are 0 size (e.g. collision records) if
108  // not used
109  IOField<Field<Type> > fieldSource(*fieldIter());
110 
111  // Map - use CompactIOField to automatically write in
112  // compact form for binary format.
113  CompactIOField<Field<Type>, Type> fieldTarget
114  (
115  IOobject
116  (
117  fieldName,
118  meshTarget.time().timeName(),
120  meshTarget,
123  false
124  ),
125  min(fieldSource.size(), addParticles.size()) // handle 0 size
126  );
127 
128  if (fieldSource.size())
129  {
130  forAll(addParticles, i)
131  {
132  fieldTarget[i] = fieldSource[addParticles[i]];
133  }
134  }
135  else if (cloud::debug)
136  {
137  Pout<< "Not mapping " << fieldName << " since source size "
138  << fieldSource.size() << " different to"
139  << " cloud size " << addParticles.size()
140  << endl;
141  }
142 
143  // Write field
144  fieldTarget.write();
145  }
146  }
147 
148  {
149  IOobjectList fieldFields =
150  objects.lookupClass(CompactIOField<Field<Type>, Type>::typeName);
151 
152  forAllIter(IOobjectList, fieldFields, fieldIter)
153  {
154  const word& fieldName = fieldIter()->name();
155 
156  Info<< " mapping lagrangian fieldField " << fieldName << endl;
157 
158  // Read field (does not need mesh)
159  CompactIOField<Field<Type>, Type> fieldSource(*fieldIter());
160 
161  // Map
162  CompactIOField<Field<Type>, Type> fieldTarget
163  (
164  IOobject
165  (
166  fieldName,
167  meshTarget.time().timeName(),
169  meshTarget,
172  false
173  ),
174  min(fieldSource.size(), addParticles.size()) // handle 0 size
175  );
176 
177  if (fieldSource.size())
178  {
179  forAll(addParticles, i)
180  {
181  fieldTarget[i] = fieldSource[addParticles[i]];
182  }
183  }
184  else if (cloud::debug)
185  {
186  Pout<< "Not mapping " << fieldName << " since source size "
187  << fieldSource.size() << " different to"
188  << " cloud size " << addParticles.size()
189  << endl;
190  }
191 
192  // Write field
193  fieldTarget.write();
194  }
195  }
196 }
197 
198 
199 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
200 
201 } // End namespace Foam
202 
203 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
204 
205 #endif
206 
207 // ************************************************************************* //
Foam::IOobject
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:91
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
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::objectRegistry::time
const Time & time() const
Return time.
Definition: objectRegistry.H:117
meshToMesh.H
cloud.H
Foam::regIOobject::write
virtual bool write() const
Write using setting from DB.
Definition: regIOobjectWrite.C:126
IOobjectList.H
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
Foam::CompactIOField
A Field of objects of type <T> with automated input and output using a compact storage....
Definition: CompactIOField.H:50
Foam::polyMesh
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
Foam::IOobject::NO_WRITE
@ NO_WRITE
Definition: IOobject.H:118
Foam::Field< Type >
Foam::IOobject::NO_READ
@ NO_READ
Definition: IOobject.H:111
Foam::Info
messageStream Info
Foam::MapLagrangianFields
void MapLagrangianFields(const string &cloudName, const IOobjectList &objects, const meshToMesh0 &meshToMesh0Interp, const labelList &addParticles)
Gets the indices of (source)particles that have been appended to the.
Definition: MapLagrangianFields.H:51
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
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
GeometricField.H
Foam::cloud::prefix
static const word prefix
The prefix to local: lagrangian.
Definition: cloud.H:71
Foam::Pout
prefixOSstream Pout(cout, "Pout")
Definition: IOstreams.H:53
Foam::Time::timeName
static word timeName(const scalar, const int precision=precision_)
Return time name of given scalar time.
Definition: Time.C:741
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::min
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)