Loading...
Searching...
No Matches
MapGeometricFields.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) 2017-2023 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
27Class
28 Foam::MapInternalField
29
30Description
31 Generic internal field mapper. For "real" mapping, add template
32 specialisations for mapping of internal fields depending on mesh
33 type.
34
35\*---------------------------------------------------------------------------*/
36
37#ifndef MapGeometricFields_H
38#define MapGeometricFields_H
39
40#include "polyMesh.H"
42// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
43
44namespace Foam
46
47template<class Type, class MeshMapper, class GeoMesh>
49{
50public:
51
52 MapInternalField() = default;
53
54 void operator()
55 (
57 const MeshMapper& mapper
58 ) const;
59};
60
61
62//- Generic Geometric field mapper.
63// For "real" mapping, add template specialisations
64// for mapping of internal fields depending on mesh type.
65template
66<
67 class Type,
68 template<class> class PatchField,
69 class MeshMapper,
70 class GeoMesh
71>
73(
74 const MeshMapper& mapper
75)
76{
78
80 (
81 mapper.thisDb().objectRegistry::template csorted<FieldType>()
82 );
83
84 // It is necessary to enforce that all old-time fields are stored
85 // before the mapping is performed. Otherwise, if the
86 // old-time-level field is mapped before the field itself, sizes
87 // will not match.
88
89 for (const auto& field : fields)
90 {
91 //Note: check can be removed once pointFields are actually stored on
92 // the pointMesh instead of now on the polyMesh!
93 if (&field.mesh() == &mapper.mesh())
94 {
95 field.storeOldTimes();
96 }
97 }
98
99 for (const auto& field : fields)
100 {
101 auto& fld = field.constCast();
102
103 if (&field.mesh() == &mapper.mesh())
104 {
105 if (polyMesh::debug)
106 {
107 Info<< "Mapping "
108 << FieldType::typeName << ' ' << field.name() << endl;
109 }
110
111 // Map the internal field
113 (
114 fld.internalFieldRef(),
115 mapper
116 );
117
118 // Map the patch fields
119 auto& bfield = fld.boundaryFieldRef();
120
121 forAll(bfield, patchi)
122 {
123 // Cannot check sizes for patch fields because of
124 // empty fields in FV and because point fields get their size
125 // from the patch which has already been resized
126 //
127
128 bfield[patchi].autoMap(mapper.boundaryMap()[patchi]);
129 }
130
131 fld.instance() = fld.time().timeName();
132 }
133 else if (polyMesh::debug)
134 {
135 Info<< "Not mapping "
136 << FieldType::typeName << ' ' << field.name()
137 << " since originating mesh differs from that of mapper."
138 << endl;
139 }
140 }
141}
142
143
144// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
145
146} // End namespace Foam
147
148// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
149
150#endif
151
152// ************************************************************************* //
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))
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Generic mesh wrapper used by volMesh, surfaceMesh, pointMesh etc.
Definition GeoMesh.H:46
Generic GeometricField class.
Generic internal field mapper. For "real" mapping, add template specialisations for mapping of intern...
A list of pointers to objects of type <T>, without allocation/deallocation management of the pointers...
Definition UPtrList.H:101
rDeltaTY field()
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:519
void MapGeometricFields(const MeshMapper &mapper)
Generic Geometric field mapper.
multivariateSurfaceInterpolationScheme< scalar >::fieldTable fields
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299