Mesh that relates to a sub-section of a Lagrangian mesh. This is used to construct fields that relate to a contiguous sub-set of the Lagrangian elements. This class only stores references and the indices defining the range of the sub-set, so it is very lightweight and can be constructed and thrown away largely without consideration of expense. More...


Public Types | |
| typedef LagrangianSubMesh | Mesh |
| Mesh type. More... | |
Public Types inherited from GeoMesh< polyMesh > | |
| typedef polyMesh | Mesh |
| Mesh type. More... | |
| typedef MESH::BoundaryMesh | BoundaryMesh |
| Boundary mesh type. More... | |
Public Member Functions | |
| ClassName ("LagrangianSubMesh") | |
| Runtime type information. More... | |
| LagrangianSubMesh (const LagrangianMesh &mesh, const LagrangianGroup group, const label size, const label start) | |
| Construct from components, except for the index which is. More... | |
| LagrangianSubMesh (const LagrangianMesh &mesh, const labelList &groupOffsets, const LagrangianGroup group) | |
| Construct from Lagrangian mesh, group offsets and group. More... | |
| ~LagrangianSubMesh () | |
| Destructor. More... | |
| const LagrangianMesh & | mesh () const |
| Return the mesh. More... | |
| LagrangianGroup | group () const |
| Return the group. More... | |
| label | size () const |
| Return size. More... | |
| label | globalSize () const |
| Return size. More... | |
| bool | empty () const |
| Return whether or not the mesh is empty. More... | |
| label | start () const |
| Return start. More... | |
| label | end () const |
| Return end. More... | |
| label | index () const |
| Return the index. More... | |
| template<class Type > | |
| SubList< Type > | sub (const List< Type > &list) const |
| Return a sub-list corresponding to this sub-mesh. More... | |
| template<class Type > | |
| SubField< Type > | sub (const Field< Type > &field) const |
| Return a sub-field corresponding to this sub-mesh. More... | |
| template<class Type , template< class > class PrimitiveField> | |
| tmp< LagrangianSubSubField< Type > > | sub (const DimensionedField< Type, LagrangianMesh, PrimitiveField > &) const |
| Return a sub-dimensioned-field corresponding to this sub-mesh. More... | |
| template<class FieldType > | |
| tmp< FieldType > | nf (const LagrangianScalarInternalDynamicField &fraction) const |
| Return the face normals at the Lagrangian locations. More... | |
| template<class FieldType > | |
| tmp< FieldType > | Uf (const LagrangianScalarInternalDynamicField &fraction) const |
| Return the face velocities at the Lagrangian locations. More... | |
| void | operator+= (const LagrangianSubMesh &) |
| Add a sub-mesh to this one. Must relate to adjacent elements. More... | |
| template<> | |
| Foam::tmp< Foam::vectorField > | nf (const LagrangianSubScalarSubField &fraction) |
| template<> | |
| Foam::tmp< Foam::vectorField > | nf (const LagrangianScalarInternalDynamicField &fraction) const |
| template<> | |
| Foam::tmp< Foam::vectorField > | Uf (const LagrangianSubScalarSubField &fraction) |
| template<> | |
| Foam::tmp< Foam::vectorField > | Uf (const LagrangianScalarInternalDynamicField &fraction) const |
| template<class Type > | |
| Foam::SubList< Type > | sub (const List< Type > &list) const |
| template<class Type > | |
| Foam::SubField< Type > | sub (const Field< Type > &field) const |
| template<class Type , template< class > class PrimitiveField> | |
| Foam::tmp< Foam::LagrangianSubSubField< Type > > | sub (const DimensionedField< Type, LagrangianMesh, PrimitiveField > &field) const |
Public Member Functions inherited from GeoMesh< polyMesh > | |
| GeoMesh (const polyMesh &mesh) | |
| Construct from MESH. More... | |
| const objectRegistry & | thisDb () const |
| Return the object registry. More... | |
| const polyMesh & | operator() () const |
| Return reference to MESH. More... | |
Static Public Member Functions | |
| static label | size (const LagrangianSubMesh &subMesh) |
| Return size. More... | |
| template<class FieldType > | |
| static tmp< FieldType > | nf (const LagrangianSubScalarSubField &fraction) |
| Return the face normals at the Lagrangian locations. More... | |
| template<class FieldType > | |
| static tmp< FieldType > | Uf (const LagrangianSubScalarSubField &fraction) |
| Return the face velocities at the Lagrangian locations. More... | |
Friends | |
| class | LagrangianMesh |
| Allow the Lagrangian mesh to construct with a specified index. More... | |
Additional Inherited Members | |
Protected Attributes inherited from GeoMesh< polyMesh > | |
| const polyMesh & | mesh_ |
| Reference to Mesh. More... | |
Mesh that relates to a sub-section of a Lagrangian mesh. This is used to construct fields that relate to a contiguous sub-set of the Lagrangian elements. This class only stores references and the indices defining the range of the sub-set, so it is very lightweight and can be constructed and thrown away largely without consideration of expense.
Definition at line 58 of file LagrangianSubMesh.H.
| typedef LagrangianSubMesh Mesh |
Mesh type.
Definition at line 108 of file LagrangianSubMesh.H.
|
explicit |
Construct from components, except for the index which is.
automatically generated from the complete mesh
Definition at line 61 of file LagrangianSubMesh.C.
|
explicit |
Construct from Lagrangian mesh, group offsets and group.
Definition at line 78 of file LagrangianSubMesh.C.
| ~LagrangianSubMesh | ( | ) |
Destructor.
Definition at line 100 of file LagrangianSubMesh.C.
| ClassName | ( | "LagrangianSubMesh" | ) |
Runtime type information.
|
inline |
Return the mesh.
Definition at line 30 of file LagrangianSubMeshI.H.
Referenced by cellPointLagrangianAccumulator::accumulate(), CrankNicolson< Type >::LagrangiancDdt(), CrankNicolson< Type >::LagrangianmInitDdt(), CrankNicolson< Type >::LagrangianmNoDdt(), patchInjection::modify(), CarrierField< scalar >::name(), LagrangianSubMesh::nf(), CarrierField< scalar >::reset(), LagrangianSubMesh::Uf(), carrierLagrangianFieldSource< Type >::value(), totalPressureVelocityLagrangianVectorFieldSource::value(), and Function1LagrangianFieldSource< Type >::value().

|
inline |
Return the group.
Definition at line 36 of file LagrangianSubMeshI.H.
Referenced by CrankNicolson< Type >::LagrangianmNoDdt(), cloudBoundaryCollisionFlux::postCrossFaces(), and cloudSurfaceDistribution::postCrossFaces().

|
inline |
Return size.
Definition at line 42 of file LagrangianSubMeshI.H.
Referenced by patchInjection::modify(), LagrangianSubMesh::nf(), cloudSurfaceDistribution::postCrossFaces(), LagrangianSubMesh::sub(), LagrangianMesh::subMeshGlobalSizes(), LagrangianSubMesh::Uf(), coneVelocityLagrangianVectorFieldSource::value(), and distributionDiameterLagrangianScalarFieldSource::value().

|
inline |
Return size.
Definition at line 48 of file LagrangianSubMeshI.H.
References Foam::returnReduce().

|
inlinestatic |
Return size.
Definition at line 54 of file LagrangianSubMeshI.H.
|
inline |
Return whether or not the mesh is empty.
Definition at line 63 of file LagrangianSubMeshI.H.
|
inline |
Return start.
Definition at line 69 of file LagrangianSubMeshI.H.
Referenced by cellPointLagrangianAccumulator::accumulate(), LagrangianMesh::crossFaces(), nonConformalCyclicLagrangianPatch::evaluate(), wedgeLagrangianPatch::evaluate(), patchInjection::modify(), LagrangianSubMesh::nf(), cloudBoundaryCollisionFlux::postCrossFaces(), LagrangianMesh::state(), LagrangianSubMesh::sub(), and LagrangianSubMesh::Uf().

|
inline |
Return end.
Definition at line 75 of file LagrangianSubMeshI.H.
Referenced by LagrangianModels::modify().

|
inline |
Return the index.
Definition at line 81 of file LagrangianSubMeshI.H.
Referenced by CloudDerivedField< Type >::ref(), and CloudStateField< Type >::ref().

Return a sub-list corresponding to this sub-mesh.
Referenced by parcel::calculate(), particle::calculate(), cloudAge::calculate(), cloudMass::calculate(), cloudSurfaceArea::calculate(), cloudVolume::calculate(), spherical::correct(), cyclicLagrangianPatch::evaluate(), nonConformalCyclicLagrangianPatch::evaluate(), wedgeLagrangianPatch::evaluate(), nonConformalProcessorCyclicLagrangianPatch::initEvaluate(), processorLagrangianPatch::initEvaluate(), CrankNicolson< Type >::LagrangiancDdt(), CrankNicolson< Type >::LagrangianmNoDdt(), patchInjection::modify(), cloudBoundaryCollisionFlux::postCrossFaces(), cloudFunctionObjectUList::preCrossFaces(), CloudStateField< Type >::ref(), cloud::solve(), internalLagrangianFieldSource< Type >::sourceValue(), Function1LagrangianFieldSource< Type >::value(), and LagrangianFieldSource< Type >::value().

Return a sub-field corresponding to this sub-mesh.
| tmp<LagrangianSubSubField<Type> > sub | ( | const DimensionedField< Type, LagrangianMesh, PrimitiveField > & | ) | const |
Return a sub-dimensioned-field corresponding to this sub-mesh.
|
static |
Return the face normals at the Lagrangian locations.
| tmp<FieldType> nf | ( | const LagrangianScalarInternalDynamicField & | fraction | ) | const |
Return the face normals at the Lagrangian locations.
|
static |
Return the face velocities at the Lagrangian locations.
| tmp<FieldType> Uf | ( | const LagrangianScalarInternalDynamicField & | fraction | ) | const |
Return the face velocities at the Lagrangian locations.
| void operator+= | ( | const LagrangianSubMesh & | subMesh | ) |
Add a sub-mesh to this one. Must relate to adjacent elements.
Definition at line 292 of file LagrangianSubMesh.C.
References Foam::exit(), Foam::FatalError, and FatalErrorInFunction.

| tmp< LagrangianSubVectorField > nf | ( | const LagrangianSubScalarSubField & | fraction | ) |
Definition at line 107 of file LagrangianSubMesh.C.
References LagrangianMesh::celli(), LagrangianMesh::coordinates(), LagrangianMesh::facei(), Foam::tracking::faceNormalAndDisplacement(), LagrangianMesh::faceTrii(), forAll, LagrangianMesh::mesh(), LagrangianSubMesh::mesh(), DimensionedField< Type, GeoMesh, PrimitiveField >::mesh(), tmp< T >::ref(), LagrangianSubMesh::size(), and LagrangianSubMesh::start().

| tmp< LagrangianSubVectorField > nf | ( | const LagrangianScalarInternalDynamicField & | fraction | ) | const |
Definition at line 155 of file LagrangianSubMesh.C.
References Foam::tracking::coordinates(), Foam::tracking::faceNormalAndDisplacement(), forAll, mesh, and tmp< T >::ref().

| tmp< LagrangianSubVectorField > Uf | ( | const LagrangianSubScalarSubField & | fraction | ) |
Definition at line 199 of file LagrangianSubMesh.C.
References LagrangianMesh::celli(), LagrangianMesh::coordinates(), TimeState::deltaTValue(), LagrangianMesh::facei(), Foam::tracking::faceNormalAndDisplacement(), LagrangianMesh::faceTrii(), forAll, LagrangianMesh::mesh(), LagrangianSubMesh::mesh(), DimensionedField< Type, GeoMesh, PrimitiveField >::mesh(), tmp< T >::ref(), LagrangianSubMesh::size(), LagrangianSubMesh::start(), and LagrangianMesh::time().

| tmp< LagrangianSubVectorField > Uf | ( | const LagrangianScalarInternalDynamicField & | fraction | ) | const |
Definition at line 247 of file LagrangianSubMesh.C.
References Foam::tracking::coordinates(), TimeState::deltaTValue(), Foam::tracking::faceNormalAndDisplacement(), forAll, mesh, tmp< T >::ref(), and fvMesh::time().

| Foam::SubList<Type> sub | ( | const List< Type > & | list | ) | const |
Definition at line 32 of file LagrangianSubMeshTemplates.C.
References LagrangianSubMesh::size(), and LagrangianSubMesh::start().

| Foam::SubField<Type> sub | ( | const Field< Type > & | field | ) | const |
Definition at line 42 of file LagrangianSubMeshTemplates.C.
| Foam::tmp<Foam::LagrangianSubSubField<Type> > sub | ( | const DimensionedField< Type, LagrangianMesh, PrimitiveField > & | field | ) | const |
Definition at line 53 of file LagrangianSubMeshTemplates.C.
References objectRegistry::cacheTemporaryObject(), IOobject::db(), DimensionedField< Type, GeoMesh, PrimitiveField >::dimensions(), Foam::constant::atomic::group, IOobject::instance(), IOobject::local(), IOobject::name(), Foam::name(), IOobject::NO_READ, and IOobject::NO_WRITE.

|
friend |
Allow the Lagrangian mesh to construct with a specified index.
Definition at line 98 of file LagrangianSubMesh.H.