mapFields.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 | www.openfoam.com
6 \\/ M anipulation |
7-------------------------------------------------------------------------------
8 Copyright (C) 2016-2020 OpenCFD Ltd.
9-------------------------------------------------------------------------------
10License
11 This file is part of OpenFOAM.
12
13 OpenFOAM is free software: you can redistribute it and/or modify it
14 under the terms of the GNU General Public License as published by
15 the Free Software Foundation, either version 3 of the License, or
16 (at your option) any later version.
17
18 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21 for more details.
22
23 You should have received a copy of the GNU General Public License
24 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25
26\*---------------------------------------------------------------------------*/
27
28#include "mapFields.H"
29#include "meshToMesh.H"
31
32// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
33
34namespace Foam
35{
36namespace functionObjects
37{
40}
41}
42
43
44// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
45
46void Foam::functionObjects::mapFields::createInterpolation
47(
48 const dictionary& dict
49)
50{
51 const fvMesh& meshTarget = mesh_;
52 const word mapRegionName(dict.get<word>("mapRegion"));
53
54 Info<< name() << ':' << nl
55 << " Reading mesh " << mapRegionName << endl;
56
57 mapRegionPtr_.reset
58 (
59 new fvMesh
60 (
61 IOobject
62 (
63 mapRegionName,
64 meshTarget.time().constant(),
65 meshTarget.time(),
67 )
68 )
69 );
70
71 const fvMesh& mapRegion = mapRegionPtr_();
72 const word mapMethodName(dict.get<word>("mapMethod"));
74 {
76 << type() << " " << name() << ": unknown map method "
77 << mapMethodName << nl
78 << "Available methods include: "
80 << exit(FatalError);
81 }
82
84 (
86 );
87
88 // Lookup corresponding AMI method
89 word patchMapMethodName = meshToMesh::interpolationMethodAMI(mapMethod);
90
91 // Optionally override
92 if (dict.readIfPresent("patchMapMethod", patchMapMethodName))
93 {
94 Info<< " Patch mapping method: " << patchMapMethodName << endl;
95 }
96
97 Info<< " Creating mesh to mesh interpolation" << endl;
98
99 if (dict.get<bool>("consistent"))
100 {
101 interpPtr_.reset
102 (
103 new meshToMesh
104 (
105 mapRegion,
106 meshTarget,
107 mapMethodName,
108 patchMapMethodName
109 )
110 );
111 }
112 else
113 {
114 HashTable<word> patchMap;
115 wordList cuttingPatches;
116
117 dict.readEntry("patchMap", patchMap);
118 dict.readEntry("cuttingPatches", cuttingPatches);
119
120 interpPtr_.reset
121 (
122 new meshToMesh
123 (
124 mapRegion,
125 meshTarget,
126 mapMethodName,
127 patchMapMethodName,
128 patchMap,
129 cuttingPatches
130 )
131 );
132 }
133}
134
135
136// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
137
139(
140 const word& name,
141 const Time& runTime,
142 const dictionary& dict
143)
144:
146 mapRegionPtr_(),
147 interpPtr_(),
148 fieldNames_()
149{
150 read(dict);
151}
152
153
154// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
155
157{
159 {
160 dict.readEntry("fields", fieldNames_);
161 createInterpolation(dict);
162
163 return true;
164 }
165
166 return false;
167}
168
169
171{
172 Log << type() << " " << name() << " execute:" << nl;
173
174 bool ok = false;
175
176 ok = mapFieldType<scalar>() || ok;
177 ok = mapFieldType<vector>() || ok;
178 ok = mapFieldType<sphericalTensor>() || ok;
179 ok = mapFieldType<symmTensor>() || ok;
180 ok = mapFieldType<tensor>() || ok;
181
182 if (log)
183 {
184 if (!ok)
185 {
186 Info<< " none" << nl;
187 }
188
189 Info<< endl;
190 }
191 return true;
192}
193
194
196{
197 Log << type() << " " << name() << " write:" << nl;
198
199 bool ok = false;
200
201 ok = writeFieldType<scalar>() || ok;
202 ok = writeFieldType<vector>() || ok;
203 ok = writeFieldType<sphericalTensor>() || ok;
204 ok = writeFieldType<symmTensor>() || ok;
205 ok = writeFieldType<tensor>() || ok;
206
207 if (log)
208 {
209 if (!ok)
210 {
211 Info<< " none" << nl;
212 }
213
214 Info<< endl;
215 }
216
217 return true;
218}
219
220
221// ************************************************************************* //
#define Log
Definition: PDRblock.C:35
bool found
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
virtual bool read()
Re-read model coefficients if they have changed.
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:80
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
T get(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
bool readIfPresent(const word &keyword, T &val, enum keyType::option matchOpt=keyType::REGEX) const
bool readEntry(const word &keyword, T &val, enum keyType::option matchOpt=keyType::REGEX, bool mandatory=true) const
Abstract base-class for Time/database function objects.
const word & name() const noexcept
Return the name of this functionObject.
virtual const word & type() const =0
Runtime type information.
Specialization of Foam::functionObject for an Foam::fvMesh, providing a reference to the Foam::fvMesh...
const fvMesh & mesh_
Reference to the fvMesh.
Computes the natural logarithm of an input volScalarField.
Definition: log.H:230
Maps input fields from local mesh to secondary mesh at runtime.
Definition: mapFields.H:225
virtual bool execute()
Execute, currently does nothing.
Definition: mapFields.C:170
virtual bool write()
Calculate the mapFields and write.
Definition: mapFields.C:195
virtual bool read(const dictionary &)
Read the mapFields data.
Definition: mapFields.C:156
interpolationMethod
Enumeration specifying interpolation method.
Definition: meshToMesh.H:72
static word interpolationMethodAMI(const interpolationMethod method)
Conversion between mesh and patch interpolation methods.
Definition: meshToMesh.C:644
static const Enum< interpolationMethod > interpolationMethodNames_
Definition: meshToMesh.H:79
A class for handling words, derived from Foam::string.
Definition: word.H:68
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition: className.H:121
engineTime & runTime
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
Namespace for OpenFOAM.
List< word > wordList
A List of words.
Definition: fileName.H:63
messageStream Info
Information stream (stdout output on master, null elsewhere)
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: MSwindows.C:598
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
error FatalError
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
constexpr char nl
The newline '\n' character (0x0a)
Definition: Ostream.H:53
dictionary dict