46void Foam::functionObjects::mapFields::createInterpolation
57 <<
type() <<
" " <<
name() <<
": unknown map method "
58 << mapMethodName <<
nl
59 <<
"Available methods include: "
73 if (
dict.readIfPresent(
"patchMapMethod", patchMapMethodName))
76 <<
" Patch mapping method: " << patchMapMethodName <<
endl;
80 <<
" Creating mesh to mesh interpolation" <<
endl;
82 if (
dict.get<
bool>(
"consistent"))
97 HashTable<word> patchMap;
99 bool createPatchMap =
dict.getOrDefault<
bool>(
"createPatchMap",
false);
101 const entry* ePtr =
dict.findEntry(
"patchMap");
103 if (createPatchMap && ePtr)
106 <<
"Requested 'createPatchMap' but 'patchMap' entry provided. "
107 <<
"Please remove one of these entries"
111 if (!createPatchMap && !ePtr)
114 <<
"Either the 'createPatchMap' or 'patchMap' entry must be "
115 <<
"provided when not using the 'consistent' mapping option"
119 const polyBoundaryMesh&
pbm =
mesh_.boundaryMesh();
120 const polyBoundaryMesh& pbmT = mapRegion.
boundaryMesh();
121 DynamicList<label> cuttingIndices;
125 Info<<
" Creating patchMap and cuttingPatches" <<
endl;
127 if (
dict.found(
"cuttingPatches"))
130 <<
"Ignoring user supplied cuttingPatches when "
131 <<
"createPatchMap option is active"
137 const polyPatch& ppT = pbmT[patchiT];
141 const word& patchNameT = ppT.name();
142 const label patchi =
pbm.findPatchID(patchNameT);
146 Info<<
" Adding to cuttingPatches: "
147 << ppT.name() <<
endl;
149 cuttingIndices.push_back(ppT.index());
155 Info<<
" Adding to patchMap: " << ppT.name()
158 patchMap.set(patchNameT, patchNameT);
167 patchMap = HashTable<word>(ePtr->stream());
170 const wordRes cuttingPatchRes =
dict.get<wordRes>(
"cuttingPatches");
171 cuttingIndices.push_back(pbmT.indices(cuttingPatchRes));
174 const wordList cuttingPatches(pbmT.names(), cuttingIndices);
181 for (
const auto& ppT : pbmT)
186 && !patchMap.found(ppT.name())
190 unknownPatchNames.insert(ppT.name());
194 for (
const label patchiT : cuttingIndices)
196 const word& patchName = pbmT[patchiT].name();
198 unknownPatchNames.erase(patchName);
200 if (patchMap.found(patchName))
202 Info<<
" Removing cutting patch from patchMap: "
203 << patchName <<
endl;
205 patchMap.erase(patchName);
209 if (unknownPatchNames.size())
212 <<
"Patches not present in source mesh. "
213 <<
"Add to cuttingPatches? Patches: "
235const Foam::fvMesh& Foam::functionObjects::mapFields::lookupMapRegion()
const
238 const word mapRegionName(dict_.get<word>(
"mapRegion"));
239 const auto* mapRegionPtr = mesh_.time().findObject<fvMesh>(mapRegionName);
242 return *mapRegionPtr;
247 <<
" Reading mesh " << mapRegionName <<
endl;
249 fvMesh* ptr =
new fvMesh
254 mesh_.time().constant(),
261 ptr->objectRegistry::store();
291 dict_.readEntry(
"fields", fieldNames_);
293 createInterpolation(dict_, mesh_, lookupMapRegion());
308 const word mapRegionName(dict_.get<
word>(
"mapRegion"));
309 const fvMesh& mapRegion = lookupMapRegion();
311 if (mesh_.changing() || mapRegion.
changing())
313 createInterpolation(dict_, mesh_, mapRegion);
316 ok = mapFieldType<scalar>() || ok;
317 ok = mapFieldType<vector>() || ok;
318 ok = mapFieldType<sphericalTensor>() || ok;
319 ok = mapFieldType<symmTensor>() || ok;
320 ok = mapFieldType<tensor>() || ok;
341 ok = writeFieldType<scalar>() || ok;
342 ok = writeFieldType<vector>() || ok;
343 ok = writeFieldType<sphericalTensor>() || ok;
344 ok = writeFieldType<symmTensor>() || ok;
345 ok = writeFieldType<tensor>() || ok;
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
const polyBoundaryMesh & pbm
@ MUST_READ
Reading required.
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Abstract base-class for Time/database function objects.
const word & name() const noexcept
Return the name of this functionObject.
virtual bool read(const dictionary &dict)
Read and set the function object if its data have changed.
virtual const word & type() const =0
Runtime type information.
const fvMesh & mesh_
Reference to the fvMesh.
fvMeshFunctionObject(const fvMeshFunctionObject &)=delete
No copy construct.
Computes the natural logarithm of an input volScalarField.
Maps input fields from local mesh to secondary mesh at runtime.
mapFields(const word &name, const Time &runTime, const dictionary &dict)
Construct from name, Time and dictionary.
virtual bool read(const dictionary &dict)
Read the function-object dictionary.
virtual bool execute()
Execute the function-object operations.
virtual bool write()
Write the function-object results.
Mesh data needed to do the Finite Volume discretisation.
interpolationMethod
Enumeration specifying interpolation method.
static word interpolationMethodAMI(const interpolationMethod method)
Conversion between mesh and patch interpolation methods.
static const Enum< interpolationMethod > interpolationMethodNames_
const polyBoundaryMesh & boundaryMesh() const noexcept
Return boundary mesh.
bool changing() const noexcept
Is mesh changing (topology changing and/or moving).
static bool constraintType(const word &patchType)
Return true if the given type is a constraint type.
A class for handling words, derived from Foam::string.
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
#define WarningInFunction
Report a warning using Foam::Warning.
Function objects are OpenFOAM utilities to ease workflow configurations and enhance workflows.
List< word > wordList
List of word.
HashSet< word, Hash< word > > wordHashSet
A HashSet of words, uses string hasher.
bool read(const char *buf, int32_t &val)
Same as readInt32.
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.
dimensionedScalar log(const dimensionedScalar &ds)
T returnReduce(const T &value, BinaryOp bop, const int tag=UPstream::msgType(), const int communicator=UPstream::worldComm)
Perform reduction on a copy, using specified binary operation.
const word GlobalIOList< Tuple2< scalar, vector > >::typeName("scalarVectorTable")
Ostream & endl(Ostream &os)
Add newline and flush stream.
Ostream & indent(Ostream &os)
Indent stream.
IOerror FatalIOError
Error stream (stdout output on all processes), with additional 'FOAM FATAL IO ERROR' header text and ...
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
errorManipArg< error, int > exit(error &err, const int errNo=1)
constexpr char nl
The newline '\n' character (0x0a).
#define forAll(list, i)
Loop across all elements in list.