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-------------------------------------------------------------------------------
11License
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
27InNamespace
28 Foam
29
30Description
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 "meshToMesh0.H"
42#include "IOobjectList.H"
43#include "CompactIOField.H"
44
45// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
46
47namespace 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.
52template<class Type>
54(
55 const string& cloudName,
56 const IOobjectList& objects,
57 const meshToMesh0& meshToMesh0Interp,
58 const labelList& addParticles
59)
60{
61 const fvMesh& meshTarget = meshToMesh0Interp.toMesh();
62
63 {
65
66 forAllConstIters(fields, fieldIter)
67 {
68 Info<< " mapping lagrangian field "
69 << (*fieldIter)->name() << endl;
70
71 // Read field (does not need mesh)
72 IOField<Type> fieldSource(*fieldIter());
73
74 // Map
75 IOField<Type> fieldTarget
76 (
78 (
79 (*fieldIter)->name(),
80 meshTarget.time().timeName(),
82 meshTarget,
85 false
86 ),
87 addParticles.size()
88 );
89
90 forAll(addParticles, i)
91 {
92 fieldTarget[i] = fieldSource[addParticles[i]];
93 }
94
95 // Write field
96 fieldTarget.write();
97 }
98 }
99
100 {
101 IOobjectList fieldFields =
102 objects.lookupClass(IOField<Field<Type>>::typeName);
103
104 forAllConstIters(fieldFields, fieldIter)
105 {
106 Info<< " mapping lagrangian fieldField "
107 << (*fieldIter)->name() << endl;
108
109 // Read field (does not need mesh)
110 IOField<Field<Type>> fieldSource(*fieldIter());
111
112 // Map - use CompactIOField to automatically write in
113 // compact form for binary format.
114 CompactIOField<Field<Type>, Type> fieldTarget
115 (
117 (
118 (*fieldIter)->name(),
119 meshTarget.time().timeName(),
121 meshTarget,
124 false
125 ),
126 addParticles.size()
127 );
128
129 forAll(addParticles, i)
130 {
131 fieldTarget[i] = fieldSource[addParticles[i]];
132 }
133
134 // Write field
135 fieldTarget.write();
136 }
137 }
138
139 {
140 IOobjectList fieldFields =
141 objects.lookupClass(CompactIOField<Field<Type>, Type>::typeName);
142
143 forAllConstIters(fieldFields, fieldIter)
144 {
145 Info<< " mapping lagrangian fieldField "
146 << (*fieldIter)->name() << endl;
147
148 // Read field (does not need mesh)
149 CompactIOField<Field<Type>, Type> fieldSource(*fieldIter());
150
151 // Map
152 CompactIOField<Field<Type>, Type> fieldTarget
153 (
155 (
156 (*fieldIter)->name(),
157 meshTarget.time().timeName(),
159 meshTarget,
162 false
163 ),
164 addParticles.size()
165 );
166
167 forAll(addParticles, i)
168 {
169 fieldTarget[i] = fieldSource[addParticles[i]];
170 }
171
172 // Write field
173 fieldTarget.write();
174 }
175 }
176}
177
178
179// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
180
181} // End namespace Foam
182
183// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
184
185#endif
186
187// ************************************************************************* //
A Field of objects of type <T> with automated input and output using a compact storage....
Generic templated field type.
Definition: Field.H:82
A primitive field of type <T> with automated input and output.
Definition: IOField.H:58
List of IOobjects with searching and retrieving facilities.
Definition: IOobjectList.H:59
IOobjectList lookupClass(const char *clsName) const
The list of IOobjects with the given headerClassName.
Definition: IOobjectList.C:313
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:170
static word timeName(const scalar t, const int precision=precision_)
Definition: Time.C:780
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:114
static const word prefix
The prefix to local: lagrangian.
Definition: cloud.H:87
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:91
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:290
Serial mesh to mesh interpolation class.
Definition: meshToMesh0.H:66
const fvMesh & toMesh() const
Definition: meshToMesh0.H:238
virtual bool write(const bool valid=true) const
Write using setting from DB.
Namespace for OpenFOAM.
messageStream Info
Information stream (stdout output on master, null elsewhere)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
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:333
#define forAllConstIters(container, iter)
Iterate across all elements of the container object with const access.
Definition: stdFoam.H:278
const word cloudName(propsDict.get< word >("cloud"))