Files
openfoam/applications/utilities/postProcessing/dataConversion/foamToEnsight/meshSubsetHelper.H
Mark Olesen e57ca15bda ENH: complete reworking of foamToEnsight to make into a library (issue #241)
- eliminate ensightAsciiStream, ensightBinaryStream, ensightStream in
  favour of using ensightFile and ensightGeoFile classes throughout.

- encapsulate mesh-parts sorting with the ensightCells, ensightFaces
  class.

- handle of patches/faceZones entirely within ensightMesh for a lighter
  interaction with field output. Both faceZones and point fields need
  more testing to see if they behave properly for all cases.

- move some output functionality into its own namespace
  'ensightOutput', move into a library.

- use the ensightCase class to open new ensight output streams
  in the proper sub-directory locations.
2016-10-06 10:43:22 +02:00

194 lines
4.9 KiB
C++

/*---------------------------------------------------------------------------*\
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\ / A nd | Copyright (C) 2011-2016 OpenFOAM Foundation
\\/ M anipulation | Copyright (C) 2016 OpenCFD Ltd.
-------------------------------------------------------------------------------
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/>.
Class
Foam::meshSubsetHelper
Description
Simple helper to hold a mesh or mesh-subset and provide uniform access.
SourceFiles
meshSubsetHelper.C
meshSubsetHelperTemplates.C
\*---------------------------------------------------------------------------*/
#ifndef meshSubsetHelper_H
#define meshSubsetHelper_H
#include "fvMeshSubset.H"
#include "zeroGradientFvPatchField.H"
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
namespace Foam
{
// Forward declaration of classes
class Time;
/*---------------------------------------------------------------------------*\
Class meshSubsetHelper Declaration
\*---------------------------------------------------------------------------*/
class meshSubsetHelper
{
// Private data
//- Reference to mesh
fvMesh& baseMesh_;
//- Subsetting engine + sub-fvMesh
fvMeshSubset subsetter_;
//- Name of current cellSet/cellZone (or empty)
const word name_;
//- Internal book-keeping. 0 = unused, 1 = set, 2 = zone
const int type_;
// Private Member Functions
//- Disallow default bitwise copy construct
meshSubsetHelper(const meshSubsetHelper&) = delete;
//- Disallow default bitwise assignment
void operator=(const meshSubsetHelper&) = delete;
public:
// Constructors
//- Construct from components
meshSubsetHelper
(
fvMesh& baseMesh,
const word& name = word::null,
const bool isCellSet = false
);
// Member Functions
// Access
//- The entire base mesh
inline const fvMesh& baseMesh() const
{
return baseMesh_;
}
//- The mesh subsetter
inline const fvMeshSubset& subsetter() const
{
return subsetter_;
}
//- Check if running a sub-mesh is being used
inline bool useSubMesh() const
{
return type_;
}
//- Access either mesh or submesh
inline const fvMesh& mesh() const
{
if (useSubMesh())
{
return subsetter_.subMesh();
}
else
{
return baseMesh_;
}
}
// Edit
//- Update mesh subset
void correct(bool verbose = false);
//- Read mesh
polyMesh::readUpdateState readUpdate();
//- Construct volField (with zeroGradient) from an internal field
template<class Type>
static tmp<GeometricField<Type, fvPatchField, volMesh>>
zeroGradientField
(
const typename GeometricField
<
Type,
fvPatchField,
volMesh
>::Internal& df
);
//- Wrapper for field or the subsetted field.
// Map volume field (does in fact do very little interpolation;
// just copied from fvMeshSubset)
template<class Type>
tmp<GeometricField<Type, fvPatchField, volMesh>>
interpolate
(
const GeometricField<Type, fvPatchField, volMesh>&
) const;
//- Convert an internal field to a volume field
template<class Type>
tmp<GeometricField<Type, fvPatchField, volMesh>>
interpolate
(
const typename GeometricField
<
Type,
fvPatchField,
volMesh
>::Internal&
) const;
};
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
} // End namespace Foam
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
#ifdef NoRepository
#include "meshSubsetHelperTemplates.C"
#endif
#endif
// ************************************************************************* //