1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171 172 173 174 175 176 177 178 179 180 181 182 183 184
|
/*---------------------------------------------------------------------------*\
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\ / A nd | Copyright (C) 2011-2016 OpenFOAM Foundation
\\/ M anipulation |
-------------------------------------------------------------------------------
License
This file is part of OpenFOAM.
OpenFOAM is free software: you can redistribute it and/or modify it
under the terms of the GNU General Public License as published by
the Free Software Foundation, either version 3 of the License, or
(at your option) any later version.
OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
for more details.
You should have received a copy of the GNU General Public License
along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
InNamespace
Foam
Description
Gets the indices of (source)particles that have been appended to the
target cloud and maps the lagrangian fields accordingly.
\*---------------------------------------------------------------------------*/
#ifndef MapLagrangianFields_H
#define MapLagrangianFields_H
#include "cloud.H"
#include "GeometricField.H"
#include "meshToMesh0.H"
#include "IOobjectList.H"
#include "CompactIOField.H"
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
namespace Foam
{
//- Gets the indices of (source)particles that have been appended to the
// target cloud and maps the lagrangian fields accordingly.
template<class Type>
void MapLagrangianFields
(
const string& cloudName,
const IOobjectList& objects,
const meshToMesh0& meshToMesh0Interp,
const labelList& addParticles
)
{
const fvMesh& meshTarget = meshToMesh0Interp.toMesh();
{
IOobjectList fields = objects.lookupClass(IOField<Type>::typeName);
forAllIter(IOobjectList, fields, fieldIter)
{
Info<< " mapping lagrangian field "
<< fieldIter()->name() << endl;
// Read field (does not need mesh)
IOField<Type> fieldSource(*fieldIter());
// Map
IOField<Type> fieldTarget
(
IOobject
(
fieldIter()->name(),
meshTarget.time().timeName(),
cloud::prefix/cloudName,
meshTarget,
IOobject::NO_READ,
IOobject::NO_WRITE,
false
),
addParticles.size()
);
forAll(addParticles, i)
{
fieldTarget[i] = fieldSource[addParticles[i]];
}
// Write field
fieldTarget.write();
}
}
{
IOobjectList fieldFields =
objects.lookupClass(IOField<Field<Type>>::typeName);
forAllIter(IOobjectList, fieldFields, fieldIter)
{
Info<< " mapping lagrangian fieldField "
<< fieldIter()->name() << endl;
// Read field (does not need mesh)
IOField<Field<Type>> fieldSource(*fieldIter());
// Map - use CompactIOField to automatically write in
// compact form for binary format.
CompactIOField<Field<Type>, Type> fieldTarget
(
IOobject
(
fieldIter()->name(),
meshTarget.time().timeName(),
cloud::prefix/cloudName,
meshTarget,
IOobject::NO_READ,
IOobject::NO_WRITE,
false
),
addParticles.size()
);
forAll(addParticles, i)
{
fieldTarget[i] = fieldSource[addParticles[i]];
}
// Write field
fieldTarget.write();
}
}
{
IOobjectList fieldFields =
objects.lookupClass(CompactIOField<Field<Type>, Type>::typeName);
forAllIter(IOobjectList, fieldFields, fieldIter)
{
Info<< " mapping lagrangian fieldField "
<< fieldIter()->name() << endl;
// Read field (does not need mesh)
CompactIOField<Field<Type>, Type> fieldSource(*fieldIter());
// Map
CompactIOField<Field<Type>, Type> fieldTarget
(
IOobject
(
fieldIter()->name(),
meshTarget.time().timeName(),
cloud::prefix/cloudName,
meshTarget,
IOobject::NO_READ,
IOobject::NO_WRITE,
false
),
addParticles.size()
);
forAll(addParticles, i)
{
fieldTarget[i] = fieldSource[addParticles[i]];
}
// Write field
fieldTarget.write();
}
}
}
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
} // End namespace Foam
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
#endif
// ************************************************************************* //
|