Loading...
Searching...
No Matches
MapVolFields.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-2017 OpenFOAM Foundation
9 Copyright (C) 2018-2025 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
27\*---------------------------------------------------------------------------*/
28
29#ifndef Foam_MapConsistentVolFields_H
30#define Foam_MapConsistentVolFields_H
31
32#include "GeometricField.H"
33#include "meshToMesh.H"
34#include "IOobjectList.H"
36// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
37
38namespace Foam
39{
40
41template<class Type>
43(
45)
46{
47 fld.boundaryFieldRef().evaluate_if
48 (
49 [](const auto& pfld) -> bool
50 {
51 return
52 (
53 pfld.type() == pfld.patch().patch().type()
54 && polyPatch::constraintType(pfld.patch().patch().type())
55 );
56 },
58 );
59}
60
61
62template<class Type, class CombineOp>
63void MapVolFields
64(
65 const IOobjectList& objects,
66 const wordRes& selectedFields,
67 const meshToMesh& interp,
68 const CombineOp& cop
69)
70{
72
73 const fvMesh& meshSource = static_cast<const fvMesh&>(interp.srcRegion());
74 const fvMesh& meshTarget = static_cast<const fvMesh&>(interp.tgtRegion());
75
76 // Available fields, sorted order
77 for
78 (
79 const IOobject& io :
80 (
81 selectedFields.empty()
82 ? objects.csorted<fieldType>()
83 : objects.csorted<fieldType>(selectedFields)
84 )
85 )
86 {
87 const fieldType fieldSource(io, meshSource, false);
88
89 IOobject targetIO
90 (
91 io.name(),
92 meshTarget.time().timeName(),
93 meshTarget,
95 );
96
97 if (targetIO.typeHeaderOk<fieldType>(true))
98 {
99 Info<< " interpolating onto existing field "
100 << targetIO.name() << endl;
101
102 fieldType fieldTarget(targetIO, meshTarget, false);
103
104 interp.mapSrcToTgt(fieldSource, cop, fieldTarget);
105
106 evaluateConstraintTypes(fieldTarget);
107
108 fieldTarget.write();
109 }
110 else
111 {
112 Info<< " creating new field "
113 << targetIO.name() << endl;
114
115 targetIO.readOpt(IOobject::NO_READ);
116
117 tmp<fieldType> tfieldTarget
118 (
119 interp.mapSrcToTgt(fieldSource, cop)
120 );
121
122 fieldType fieldTarget(targetIO, tfieldTarget);
123
124 evaluateConstraintTypes(fieldTarget);
125
126 fieldTarget.write();
127 }
128 }
129}
130
131
132// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
133
134} // End namespace Foam
135
136// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
137
138#endif
139
140// ************************************************************************* //
Info<< nl;Info<< "Write faMesh in vtk format:"<< nl;{ vtk::uindirectPatchWriter writer(aMesh.patch(), fileName(aMesh.time().globalPath()/vtkBaseFileName));writer.writeGeometry();globalIndex procAddr(aMesh.nFaces());labelList cellIDs;if(UPstream::master()) { cellIDs.resize(procAddr.totalSize());for(const labelRange &range :procAddr.ranges()) { auto slice=cellIDs.slice(range);slice=identity(range);} } writer.beginCellData(4);writer.writeProcIDs();writer.write("cellID", cellIDs);writer.write("area", aMesh.S().field());writer.write("normal", aMesh.faceAreaNormals());writer.beginPointData(1);writer.write("normal", aMesh.pointAreaNormals());Info<< " "<< writer.output().name()<< nl;}{ vtk::lineWriter writer(aMesh.points(), aMesh.edges(), fileName(aMesh.time().globalPath()/(vtkBaseFileName+"-edges")));writer.writeGeometry();writer.beginCellData(4);writer.writeProcIDs();{ Field< scalar > fld(faMeshTools::flattenEdgeField(aMesh.magLe(), true))
Generic GeometricField class.
List of IOobjects with searching and retrieving facilities. Implemented as a HashTable,...
UPtrList< const IOobject > csorted() const
The sorted list of IOobjects with headerClassName == Type::typeName.
readOption readOpt() const noexcept
Get the read option.
@ NO_READ
Nothing to be read.
@ MUST_READ
Reading required.
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition IOobject.H:191
const word & name() const noexcept
Return the object name.
Definition IOobjectI.H:205
bool typeHeaderOk(const bool checkType=true, const bool search=true, const bool verbose=true)
Read header (respects is_globalIOobject trait) and check its info. A void type suppresses trait and t...
static word timeName(const scalar t, const int precision=precision_)
Return a time name for the given scalar time value formatted with the given precision.
Definition Time.C:714
bool empty() const noexcept
True if List is empty (ie, size() is zero).
Definition UList.H:701
static commsTypes defaultCommsType
Default commsType.
Definition UPstream.H:1045
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
const Time & time() const
Return the top-level database.
Definition fvMesh.H:360
Class to calculate the cell-addressing between two overlapping meshes.
Definition meshToMesh.H:61
const polyMesh & srcRegion() const
Return const access to the source mesh.
Definition meshToMeshI.H:26
void mapSrcToTgt(const UList< Type > &srcFld, const CombineOp &cop, List< Type > &result) const
Map field from src to tgt mesh with defined operation.
const polyMesh & tgtRegion() const
Return const access to the target mesh.
Definition meshToMeshI.H:32
static bool constraintType(const word &patchType)
Return true if the given type is a constraint type.
Definition polyPatch.C:255
A class for managing temporary objects.
Definition tmp.H:75
A List of wordRe with additional matching capabilities.
Definition wordRes.H:56
const auto & io
Namespace for OpenFOAM.
void MapVolFields(const IOobjectList &objects, const meshToMesh0 &meshToMesh0Interp, const meshToMesh0::order &mapOrder, const CombineOp &cop)
messageStream Info
Information stream (stdout output on master, null elsewhere).
void evaluateConstraintTypes(GeometricField< Type, fvPatchField, volMesh > &fld)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519