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 | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2016 OpenFOAM Foundation
9  Copyright (C) 2019 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 InNamespace
28  Foam
29 
30 Description
31  Gets the indices of (source)particles that have been appended to the
32  target cloud and maps the lagrangian fields accordingly.
33 
34 \*---------------------------------------------------------------------------*/
35 
36 #ifndef MapLagrangianFields_H
37 #define MapLagrangianFields_H
38 
39 #include "cloud.H"
40 #include "GeometricField.H"
41 #include "meshToMesh.H"
42 #include "IOobjectList.H"
43 #include "CompactIOField.H"
44 
45 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
46 
47 namespace Foam
48 {
49 
50 //- Gets the indices of (source)particles that have been appended to the
51 // target cloud and maps the lagrangian fields accordingly.
52 template<class Type>
54 (
55  const string& cloudName,
56  const IOobjectList& objects,
57  const polyMesh& meshTarget,
58  const labelList& addParticles
59 )
60 {
61  {
63 
64  forAllConstIters(fields, fieldIter)
65  {
66  const word& fieldName = (*fieldIter)->name();
67 
68  Info<< " mapping lagrangian field " << fieldName << endl;
69 
70  // Read field (does not need mesh)
71  IOField<Type> fieldSource(*fieldIter());
72 
73  // Map
74  IOField<Type> fieldTarget
75  (
76  IOobject
77  (
78  fieldName,
79  meshTarget.time().timeName(),
81  meshTarget,
84  false
85  ),
86  addParticles.size()
87  );
88 
89  forAll(addParticles, i)
90  {
91  fieldTarget[i] = fieldSource[addParticles[i]];
92  }
93 
94  // Write field
95  fieldTarget.write();
96  }
97  }
98 
99  {
100  IOobjectList fieldFields =
101  objects.lookupClass(IOField<Field<Type>>::typeName);
102 
103  forAllConstIters(fieldFields, fieldIter)
104  {
105  const word& fieldName = (*fieldIter)->name();
106 
107  Info<< " mapping lagrangian fieldField " << fieldName << endl;
108 
109  // Read field (does not need mesh)
110  // Note: some fieldFields are 0 size (e.g. collision records) if
111  // not used
112  IOField<Field<Type>> fieldSource(*fieldIter());
113 
114  // Map - use CompactIOField to automatically write in
115  // compact form for binary format.
116  CompactIOField<Field<Type>, Type> fieldTarget
117  (
118  IOobject
119  (
120  fieldName,
121  meshTarget.time().timeName(),
123  meshTarget,
126  false
127  ),
128  min(fieldSource.size(), addParticles.size()) // handle 0 size
129  );
130 
131  if (fieldSource.size())
132  {
133  forAll(addParticles, i)
134  {
135  fieldTarget[i] = fieldSource[addParticles[i]];
136  }
137  }
138  else if (cloud::debug)
139  {
140  Pout<< "Not mapping " << fieldName << " since source size "
141  << fieldSource.size() << " different to"
142  << " cloud size " << addParticles.size()
143  << endl;
144  }
145 
146  // Write field
147  fieldTarget.write();
148  }
149  }
150 
151  {
152  IOobjectList fieldFields =
153  objects.lookupClass(CompactIOField<Field<Type>, Type>::typeName);
154 
155  forAllConstIters(fieldFields, fieldIter)
156  {
157  const word& fieldName = (*fieldIter)->name();
158 
159  Info<< " mapping lagrangian fieldField " << fieldName << endl;
160 
161  // Read field (does not need mesh)
162  CompactIOField<Field<Type>, Type> fieldSource(*fieldIter());
163 
164  // Map
165  CompactIOField<Field<Type>, Type> fieldTarget
166  (
167  IOobject
168  (
169  fieldName,
170  meshTarget.time().timeName(),
172  meshTarget,
175  false
176  ),
177  min(fieldSource.size(), addParticles.size()) // handle 0 size
178  );
179 
180  if (fieldSource.size())
181  {
182  forAll(addParticles, i)
183  {
184  fieldTarget[i] = fieldSource[addParticles[i]];
185  }
186  }
187  else if (cloud::debug)
188  {
189  Pout<< "Not mapping " << fieldName << " since source size "
190  << fieldSource.size() << " different to"
191  << " cloud size " << addParticles.size()
192  << endl;
193  }
194 
195  // Write field
196  fieldTarget.write();
197  }
198  }
199 }
200 
201 
202 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
203 
204 } // End namespace Foam
205 
206 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
207 
208 #endif
209 
210 // ************************************************************************* //
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:118
List of IOobjects with searching and retrieving facilities. Implemented as a HashTable, so the various sorted methods should be used if traversing in parallel.
Definition: IOobjectList.H:55
IOobjectList lookupClass(const char *clsName) const
The list of IOobjects with the given headerClassName.
Definition: IOobjectList.C:231
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:487
Ignore writing from objectRegistry::writeObject()
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.
multivariateSurfaceInterpolationScheme< scalar >::fieldTable fields
Definition: createFields.H:97
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:413
Generic templated field type.
Definition: Field.H:61
const word cloudName(propsDict.get< word >("cloud"))
A class for handling words, derived from Foam::string.
Definition: word.H:63
const Time & time() const noexcept
Return time registry.
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:26
static word timeName(const scalar t, const int precision=precision_)
Return time name of given scalar time formatted with the given precision.
Definition: Time.C:760
int debug
Static debugging option.
A Field of objects of type <T> with automated input and output using a compact storage. Behaves like IOField except when binary output in case it writes a CompactListList.
Nothing to be read.
messageStream Info
Information stream (stdout output on master, null elsewhere)
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:73
virtual bool write(const bool valid=true) const
Write using setting from DB.
Defines the attributes of an object for which implicit objectRegistry management is supported...
Definition: IOobject.H:166
A primitive field of type <T> with automated input and output.
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
Namespace for OpenFOAM.
forAllConstIters(mixture.phases(), phase)
Definition: pEqn.H:28
static const word prefix
The prefix to local: lagrangian.
Definition: cloud.H:93