linear< Type > Class Template Reference

Centred interpolation interpolation scheme class. More...

Inheritance diagram for linear< Type >:
Collaboration diagram for linear< Type >:

Public Member Functions

 TypeName ("linear")
 Runtime type information. More...
 
 linear (const fvMesh &mesh)
 Construct from mesh. More...
 
 linear (const fvMesh &mesh, Istream &)
 Construct from Istream. More...
 
 linear (const fvMesh &mesh, const surfaceScalarField &, Istream &)
 Construct from faceFlux and Istream. More...
 
tmp< surfaceScalarFieldweights (const VolField< Type > &) const
 Return the interpolation weighting factors. More...
 
void operator= (const linear &)=delete
 Disallow default bitwise assignment. More...
 
 linear (const Specie &sp, const scalar psi, const scalar rho0)
 Construct from components. More...
 
 linear (const word &name, const dictionary &dict)
 Construct from name and dictionary. More...
 
 linear (const word &name, const linear &)
 Construct as named copy. More...
 
autoPtr< linearclone () const
 Construct and return a clone. More...
 
scalar rho (scalar p, scalar T) const
 Return density [kg/m^3]. More...
 
scalar h (const scalar p, const scalar T) const
 Return enthalpy contribution [J/kg]. More...
 
scalar Cp (scalar p, scalar T) const
 Return Cp contribution [J/(kg K]. More...
 
scalar e (const scalar p, const scalar T) const
 Return internal energy contribution [J/kg]. More...
 
scalar Cv (scalar p, scalar T) const
 Return Cv contribution [J/(kg K]. More...
 
scalar sp (const scalar p, const scalar T) const
 Return entropy contribution to the integral of Cp/T [J/kg/K]. More...
 
scalar sv (const scalar p, const scalar T) const
 Return entropy contribution to the integral of Cv/T [J/kg/K]. More...
 
scalar psi (scalar p, scalar T) const
 Return compressibility [s^2/m^2]. More...
 
scalar Z (scalar p, scalar T) const
 Return compression factor []. More...
 
scalar CpMCv (scalar p, scalar T) const
 Return (Cp - Cv) [J/(kg K]. More...
 
scalar alphav (const scalar p, const scalar T) const
 Return volumetric coefficient of thermal expansion [1/T]. More...
 
void write (Ostream &os) const
 Write to Ostream. More...
 
void operator+= (const linear &)
 
void operator*= (const scalar)
 
- Public Member Functions inherited from surfaceInterpolationScheme< Type >
 TypeName ("surfaceInterpolationScheme")
 Runtime type information. More...
 
 declareRunTimeSelectionTable (tmp, surfaceInterpolationScheme, Mesh,(const fvMesh &mesh, Istream &schemeData),(mesh, schemeData))
 
 declareRunTimeSelectionTable (tmp, surfaceInterpolationScheme, MeshFlux,(const fvMesh &mesh, const surfaceScalarField &faceFlux, Istream &schemeData),(mesh, faceFlux, schemeData))
 
 surfaceInterpolationScheme (const fvMesh &mesh)
 Construct from mesh. More...
 
 surfaceInterpolationScheme (const surfaceInterpolationScheme &)=delete
 Disallow default bitwise copy construction. More...
 
virtual ~surfaceInterpolationScheme ()
 Destructor. More...
 
const fvMeshmesh () const
 Return mesh reference. More...
 
virtual bool corrected () const
 Return true if this scheme uses an explicit correction. More...
 
virtual tmp< SurfaceField< Type > > correction (const VolField< Type > &) const
 Return the explicit correction to the face-interpolate. More...
 
virtual tmp< SurfaceField< typename innerProduct< vector, Type >::type > > dotInterpolate (const surfaceVectorField &Sf, const VolField< Type > &vf) const
 Return the face-interpolate of the given cell field. More...
 
tmp< SurfaceField< typename innerProduct< vector, Type >::type > > dotInterpolate (const surfaceVectorField &Sf, const tmp< VolField< Type >> &) const
 Return the face-interpolate of the given tmp cell field. More...
 
virtual tmp< SurfaceField< Type > > interpolate (const VolField< Type > &) const
 Return the face-interpolate of the given cell field. More...
 
tmp< SurfaceField< Type > > interpolate (const tmp< VolField< Type >> &) const
 Return the face-interpolate of the given tmp cell field. More...
 
void operator= (const surfaceInterpolationScheme &)=delete
 Disallow default bitwise assignment. More...
 
tmp< SurfaceField< typename innerProduct< vector, scalar >::type > > dotInterpolate (const surfaceVectorField &Sf, const VolField< scalar > &) const
 
- Public Member Functions inherited from refCount
int count () const
 Return the current reference count. More...
 
bool unique () const
 Return true if the reference count is zero. More...
 
void operator++ ()
 Increment the reference count. More...
 
void operator++ (int)
 Increment the reference count. More...
 
void operator-- ()
 Decrement the reference count. More...
 
void operator-- (int)
 Decrement the reference count. More...
 
void operator= (const refCount &)=delete
 Disallow bitwise assignment. More...
 

Static Public Member Functions

static word typeName ()
 Return the instantiated type name. More...
 
- Static Public Member Functions inherited from surfaceInterpolationScheme< Type >
static tmp< surfaceInterpolationScheme< Type > > New (const fvMesh &mesh, Istream &schemeData)
 Return new tmp interpolation scheme. More...
 
static tmp< surfaceInterpolationScheme< Type > > New (const fvMesh &mesh, const surfaceScalarField &faceFlux, Istream &schemeData)
 Return new tmp interpolation scheme. More...
 
static tmp< SurfaceField< Type > > interpolate (const VolField< Type > &, const tmp< surfaceScalarField > &, const tmp< surfaceScalarField > &)
 Return the face-interpolate of the given cell field. More...
 
template<class SFType >
static tmp< SurfaceField< typename innerProduct< typename SFType::value_type, Type >::type > > dotInterpolate (const SFType &Sf, const VolField< Type > &vf, const tmp< surfaceScalarField > &tlambdas)
 Return the face-interpolate of the given cell field. More...
 
static tmp< SurfaceField< Type > > interpolate (const VolField< Type > &, const tmp< surfaceScalarField > &)
 Return the face-interpolate of the given cell field. More...
 

Static Public Attributes

static const bool incompressible = false
 Is the equation of state is incompressible i.e. rho != f(p) More...
 
static const bool isochoric = false
 Is the equation of state is isochoric i.e. rho = const. More...
 

Friends

linear operator+ (const linear &, const linear &)
 
linear operator* (const scalar s, const linear &)
 
linear operator== (const linear &, const linear &)
 
Ostreamoperator (Ostream &, const linear &)
 

Additional Inherited Members

- Protected Member Functions inherited from refCount
 refCount ()
 Construct null initialising count to 0. More...
 
 refCount (const refCount &)=delete
 Disallow copy. More...
 

Detailed Description

template<class Type>
class Foam::linear< Type >

Centred interpolation interpolation scheme class.

Linear equation of state with constant compressibility:

Source files

    rho = rho0 + psi*p

Coefficient mixing is very inaccurate and not supported, so this equation of state is not applicable to mixtures.

Usage
Property Description
rho0 Reference density
psi Constant compressibility

Example specification of the linear equation of state:

    equationOfState
    {
        rho0        1000;
        psi         1e-5;
    }
Source files

Definition at line 50 of file linear.H.

Constructor & Destructor Documentation

◆ linear() [1/6]

linear ( const fvMesh mesh)
inline

Construct from mesh.

Definition at line 64 of file linear.H.

Referenced by Foam::operator<<().

Here is the caller graph for this function:

◆ linear() [2/6]

linear ( const fvMesh mesh,
Istream  
)
inline

Construct from Istream.

Definition at line 70 of file linear.H.

◆ linear() [3/6]

linear ( const fvMesh mesh,
const surfaceScalarField ,
Istream  
)
inline

Construct from faceFlux and Istream.

Definition at line 76 of file linear.H.

◆ linear() [4/6]

linear ( const Specie &  sp,
const scalar  psi,
const scalar  rho0 
)
inline

Construct from components.

Definition at line 31 of file linearI.H.

◆ linear() [5/6]

linear ( const word name,
const dictionary dict 
)

Construct from name and dictionary.

Definition at line 32 of file linear.C.

◆ linear() [6/6]

linear ( const word name,
const linear< Type > &   
)
inline

Construct as named copy.

Member Function Documentation

◆ TypeName()

TypeName ( "linear< Type >"  )

Runtime type information.

◆ weights()

tmp<surfaceScalarField> weights ( const VolField< Type > &  ) const
inlinevirtual

Return the interpolation weighting factors.

Implements surfaceInterpolationScheme< Type >.

Definition at line 90 of file linear.H.

References surfaceInterpolationScheme< Type >::mesh().

Here is the call graph for this function:

◆ operator=()

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

Disallow default bitwise assignment.

◆ clone()

Foam::autoPtr< Foam::linear< Specie > > clone
inline

Construct and return a clone.

Definition at line 59 of file linearI.H.

◆ typeName()

static word typeName ( )
inlinestatic

Return the instantiated type name.

Definition at line 153 of file linear.H.

◆ rho()

Foam::scalar rho ( scalar  p,
scalar  T 
) const
inline

Return density [kg/m^3].

Definition at line 68 of file linearI.H.

References p.

◆ h()

Foam::scalar h ( const scalar  p,
const scalar  T 
) const
inline

Return enthalpy contribution [J/kg].

Definition at line 75 of file linearI.H.

References Foam::log(), p, rho, and Foam::T().

Here is the call graph for this function:

◆ Cp()

Foam::scalar Cp ( scalar  p,
scalar  T 
) const
inline

Return Cp contribution [J/(kg K].

Definition at line 82 of file linearI.H.

◆ e()

Foam::scalar e ( const scalar  p,
const scalar  T 
) const
inline

Return internal energy contribution [J/kg].

Definition at line 89 of file linearI.H.

References Foam::log(), p, rho, and Foam::T().

Here is the call graph for this function:

◆ Cv()

Foam::scalar Cv ( scalar  p,
scalar  T 
) const
inline

Return Cv contribution [J/(kg K].

Definition at line 98 of file linearI.H.

◆ sp()

Foam::scalar sp ( const scalar  p,
const scalar  T 
) const
inline

Return entropy contribution to the integral of Cp/T [J/kg/K].

Definition at line 105 of file linearI.H.

References Foam::log(), p, Foam::constant::standard::Pstd, and Foam::T().

Referenced by Foam::operator<<().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ sv()

Foam::scalar sv ( const scalar  p,
const scalar  T 
) const
inline

Return entropy contribution to the integral of Cv/T [J/kg/K].

Definition at line 112 of file linearI.H.

References NotImplemented.

◆ psi()

Foam::scalar psi ( scalar  p,
scalar  T 
) const
inline

Return compressibility [s^2/m^2].

Definition at line 120 of file linearI.H.

Referenced by Foam::operator<<().

Here is the caller graph for this function:

◆ Z()

Foam::scalar Z ( scalar  p,
scalar  T 
) const
inline

Return compression factor [].

Definition at line 127 of file linearI.H.

References p, Foam::R(), rho, and Foam::T().

Here is the call graph for this function:

◆ CpMCv()

Foam::scalar CpMCv ( scalar  p,
scalar  T 
) const
inline

Return (Cp - Cv) [J/(kg K].

Definition at line 134 of file linearI.H.

◆ alphav()

Foam::scalar alphav ( const scalar  p,
const scalar  T 
) const
inline

Return volumetric coefficient of thermal expansion [1/T].

Definition at line 141 of file linearI.H.

◆ write()

void write ( Ostream os) const

Write to Ostream.

Definition at line 43 of file linear.C.

References dictionary::add(), dict, dictionaryName::dictName(), Foam::indent(), and Foam::vtkWriteOps::write().

Here is the call graph for this function:

◆ operator+=()

void operator+= ( const linear< Type > &  )
inline

Definition at line 150 of file linearI.H.

References noCoefficientMixing.

◆ operator*=()

void operator*= ( const scalar  s)
inline

Definition at line 160 of file linearI.H.

References s().

Here is the call graph for this function:

Friends And Related Function Documentation

◆ operator+

linear operator+ ( const linear< Type > &  ,
const linear< Type > &   
)
friend

◆ operator*

linear operator* ( const scalar  s,
const linear< Type > &   
)
friend

◆ operator==

linear operator== ( const linear< Type > &  ,
const linear< Type > &   
)
friend

◆ operator

Ostream& operator ( Ostream ,
const linear< Type > &   
)
friend

Member Data Documentation

◆ incompressible

const bool incompressible = false
static

Is the equation of state is incompressible i.e. rho != f(p)

Definition at line 162 of file linear.H.

◆ isochoric

const bool isochoric = false
static

Is the equation of state is isochoric i.e. rho = const.

Definition at line 165 of file linear.H.


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