LagrangianAverage< Type > Class Template Referenceabstract

Basic Lagrangian averaging process in which values are averaged in the cells and then interpolated assuming a constant value within the cell. More...

Inheritance diagram for LagrangianAverage< Type >:

Public Member Functions

 TypeName ("LagrangianAverage")
 Runtime type information. More...
 
 declareRunTimeSelectionTable (autoPtr, LagrangianAverage, dictionary,(const word &name, const LagrangianMesh &mesh, const dimensionSet &dimensions, const Field< scalar > &cellWeightSum),(name, mesh, dimensions, cellWeightSum))
 Declare run-time constructor selection tables. More...
 
 LagrangianAverage (const word &name, const LagrangianMesh &mesh, const dimensionSet &dimensions)
 Construct with a name, for a mesh and with given dimensions. More...
 
virtual ~LagrangianAverage ()
 Destructor. More...
 
virtual void remove (const LagrangianSubSubField< scalar > &weight, const LagrangianSubSubField< Type > &psi)=0
 Remove weighted values from the average. More...
 
void remove (const LagrangianSubSubField< Type > &weightPsi)
 Remove values from the average. More...
 
virtual void add (const LagrangianSubSubField< scalar > &weight, const LagrangianSubSubField< Type > &psi, const bool cache)=0
 Add weighted values to the average. More...
 
void add (const LagrangianSubSubField< Type > &weightPsi, const bool cache)
 Add values to the average. More...
 
virtual void correct (const LagrangianSubSubField< scalar > &weight, const LagrangianSubSubField< Type > &psi, const bool cache)=0
 Correct weighted values in the average. More...
 
void correct (const LagrangianSubSubField< Type > &weightPsi, const bool cache)
 Correct values in the average. More...
 
tmp< LagrangianSubField< Type > > interpolate (const LagrangianSubMesh &subMesh) const
 Interpolate to a sub-mesh. More...
 
void operator= (const LagrangianAverage< Type > &)=delete
 Disallow default bitwise assignment. More...
 
template<template< class > class WeightPF, template< class > class PsiPF>
Foam::autoPtr< Foam::LagrangianAverage< Type > > New (const word &type, const word &name, const DimensionedField< scalar, LagrangianMesh, WeightPF > &weight, const DimensionedField< Type, LagrangianMesh, PsiPF > &psi)
 
template<class CellMesh , template< class > class PsiPF>
Foam::LagrangianAverageNewReturnType< CellMesh, Type >::type New (const word &type, const word &name, const DimensionedField< scalar, CellMesh > &cellWeightSum, const DimensionedField< Type, LagrangianMesh, PsiPF > &weightPsi)
 

Static Public Member Functions

template<template< class > class WeightPF, template< class > class PsiPF>
static autoPtr< LagrangianAverage< Type > > New (const word &type, const word &name, const DimensionedField< scalar, LagrangianMesh, WeightPF > &weight, const DimensionedField< Type, LagrangianMesh, PsiPF > &psi)
 Select to average a field. More...
 
template<class CellMesh , template< class > class WeightPsiPF>
static LagrangianAverageNewReturnType< CellMesh, Type >::type New (const word &type, const word &name, const DimensionedField< scalar, CellMesh > &cellWeightSum, const DimensionedField< Type, LagrangianMesh, WeightPsiPF > &weightPsi)
 Select to average a field with a given cell weight sum. More...
 

Detailed Description

template<class Type>
class Foam::LagrangianAverage< Type >

Basic Lagrangian averaging process in which values are averaged in the cells and then interpolated assuming a constant value within the cell.

Base class for methods of local-averaging of Lagrangian properties.

Lagrangian averaging process in which values are averaged onto the mesh cell-centres and points using the basis functions associated with the tetrahedral decomposition. Interpolation back to the Lagrangian points is then done with the same basis functions.

Source files

Source files

Source files

Definition at line 67 of file LagrangianAverage.H.

Constructor & Destructor Documentation

◆ LagrangianAverage()

LagrangianAverage ( const word name,
const LagrangianMesh mesh,
const dimensionSet dimensions 
)

Construct with a name, for a mesh and with given dimensions.

Definition at line 32 of file LagrangianAverage.C.

◆ ~LagrangianAverage()

Destructor.

Definition at line 138 of file LagrangianAverage.C.

Member Function Documentation

◆ TypeName()

TypeName ( "LagrangianAverage< Type >"  )

Runtime type information.

◆ declareRunTimeSelectionTable()

declareRunTimeSelectionTable ( autoPtr  ,
LagrangianAverage< Type >  ,
dictionary  ,
(const word &name, const LagrangianMesh &mesh, const dimensionSet &dimensions, const Field< scalar > &cellWeightSum)  ,
(name, mesh, dimensions, cellWeightSum)   
)

Declare run-time constructor selection tables.

◆ New() [1/4]

static autoPtr<LagrangianAverage<Type> > New ( const word type,
const word name,
const DimensionedField< scalar, LagrangianMesh, WeightPF > &  weight,
const DimensionedField< Type, LagrangianMesh, PsiPF > &  psi 
)
static

Select to average a field.

◆ New() [2/4]

static LagrangianAverageNewReturnType<CellMesh, Type>::type New ( const word type,
const word name,
const DimensionedField< scalar, CellMesh > &  cellWeightSum,
const DimensionedField< Type, LagrangianMesh, WeightPsiPF > &  weightPsi 
)
static

Select to average a field with a given cell weight sum.

◆ remove() [1/2]

virtual void remove ( const LagrangianSubSubField< scalar > &  weight,
const LagrangianSubSubField< Type > &  psi 
)
pure virtual

Remove weighted values from the average.

Implemented in cellPoint< Type >, and cell< Type >.

◆ remove() [2/2]

void remove ( const LagrangianSubSubField< Type > &  weightPsi)

Remove values from the average.

Definition at line 145 of file LagrangianAverage.C.

References Foam::NullObjectRef().

Here is the call graph for this function:

◆ add() [1/2]

virtual void add ( const LagrangianSubSubField< scalar > &  weight,
const LagrangianSubSubField< Type > &  psi,
const bool  cache 
)
pure virtual

Add weighted values to the average.

Implemented in cellPoint< Type >, and cell< Type >.

◆ add() [2/2]

void add ( const LagrangianSubSubField< Type > &  weightPsi,
const bool  cache 
)

Add values to the average.

Definition at line 155 of file LagrangianAverage.C.

References Foam::add(), and Foam::NullObjectRef().

Here is the call graph for this function:

◆ correct() [1/2]

virtual void correct ( const LagrangianSubSubField< scalar > &  weight,
const LagrangianSubSubField< Type > &  psi,
const bool  cache 
)
pure virtual

Correct weighted values in the average.

Implemented in cellPoint< Type >, and cell< Type >.

◆ correct() [2/2]

void correct ( const LagrangianSubSubField< Type > &  weightPsi,
const bool  cache 
)

Correct values in the average.

Definition at line 166 of file LagrangianAverage.C.

References Foam::MULES::correct(), and Foam::NullObjectRef().

Here is the call graph for this function:

◆ interpolate()

Foam::tmp< Foam::LagrangianSubField< Type > > interpolate ( const LagrangianSubMesh subMesh) const

Interpolate to a sub-mesh.

Definition at line 178 of file LagrangianAverage.C.

References Foam::fvc::interpolate(), DimensionedField< Type, GeoMesh, PrimitiveField >::New(), and tmp< T >::ref().

Here is the call graph for this function:

◆ operator=()

void operator= ( const LagrangianAverage< Type > &  )
delete

Disallow default bitwise assignment.

◆ New() [3/4]

◆ New() [4/4]

Foam::LagrangianAverageNewReturnType<CellMesh, Type>::type New ( const word type,
const word name,
const DimensionedField< scalar, CellMesh > &  cellWeightSum,
const DimensionedField< Type, LagrangianMesh, PsiPF > &  weightPsi 
)

The documentation for this class was generated from the following files: