patchCellsSource Class Reference

Source defined by a boundary condition applied to cells next to patches. This fvOption needs to be used with a boundarySourcePatch type of boundary condition (e.g. speciesSorption and enthalpySorption.) More...

Inheritance diagram for patchCellsSource:
[legend]
Collaboration diagram for patchCellsSource:
[legend]

Public Member Functions

 TypeName ("patchCellsSource")
 Runtime type information. More...
 
 patchCellsSource (const word &sourceName, const word &modelType, const dictionary &dict, const fvMesh &mesh)
 Construct from explicit source name and mesh. More...
 
 patchCellsSource (const patchCellsSource &)=delete
 No copy construct. More...
 
void operator= (const patchCellsSource &)=delete
 No copy assignment. More...
 
virtual ~patchCellsSource ()=default
 Destructor. More...
 
virtual void addSup (const volScalarField &rho, fvMatrix< scalar > &eqn, const label fieldi)
 
virtual bool read (const dictionary &dict)
 Read source dictionary. More...
 
- Public Member Functions inherited from option
 TypeName ("option")
 Runtime type information. More...
 
 declareRunTimeSelectionTable (autoPtr, option, dictionary,(const word &name, const word &modelType, const dictionary &dict, const fvMesh &mesh),(name, modelType, dict, mesh))
 
 option (const word &name, const word &modelType, const dictionary &dict, const fvMesh &mesh)
 Construct from components. More...
 
autoPtr< optionclone () const
 Return clone. More...
 
virtual ~option ()=default
 Destructor. More...
 
const wordname () const noexcept
 Return const access to the source name. More...
 
const fvMeshmesh () const noexcept
 Return const access to the mesh database. More...
 
const dictionarycoeffs () const noexcept
 Return dictionary. More...
 
bool active () const noexcept
 True if source is active. More...
 
void setApplied (const label fieldi)
 Set the applied flag to true for field index fieldi. More...
 
bool active (const bool on) noexcept
 Change source active flag, return previous value. More...
 
virtual bool isActive ()
 Is the source active? More...
 
virtual label applyToField (const word &fieldName) const
 Return index of field name if found in fieldNames list. More...
 
virtual void checkApplied () const
 Check that the source has been applied. More...
 
virtual void addSup (fvMatrix< scalar > &eqn, const label fieldi)
 
virtual void addSup (fvMatrix< vector > &eqn, const label fieldi)
 
virtual void addSup (fvMatrix< symmTensor > &eqn, const label fieldi)
 
virtual void addSup (fvMatrix< sphericalTensor > &eqn, const label fieldi)
 
virtual void addSup (fvMatrix< tensor > &eqn, const label fieldi)
 
virtual void addSup (const volScalarField &rho, fvMatrix< scalar > &eqn, const label fieldi)
 
virtual void addSup (const volScalarField &rho, fvMatrix< vector > &eqn, const label fieldi)
 
virtual void addSup (const volScalarField &rho, fvMatrix< symmTensor > &eqn, const label fieldi)
 
virtual void addSup (const volScalarField &rho, fvMatrix< sphericalTensor > &eqn, const label fieldi)
 
virtual void addSup (const volScalarField &rho, fvMatrix< tensor > &eqn, const label fieldi)
 
virtual void addSup (const volScalarField &alpha, const volScalarField &rho, fvMatrix< scalar > &eqn, const label fieldi)
 
virtual void addSup (const volScalarField &alpha, const volScalarField &rho, fvMatrix< vector > &eqn, const label fieldi)
 
virtual void addSup (const volScalarField &alpha, const volScalarField &rho, fvMatrix< symmTensor > &eqn, const label fieldi)
 
virtual void addSup (const volScalarField &alpha, const volScalarField &rho, fvMatrix< sphericalTensor > &eqn, const label fieldi)
 
virtual void addSup (const volScalarField &alpha, const volScalarField &rho, fvMatrix< tensor > &eqn, const label fieldi)
 
virtual void constrain (fvMatrix< scalar > &eqn, const label fieldi)
 
virtual void constrain (fvMatrix< vector > &eqn, const label fieldi)
 
virtual void constrain (fvMatrix< sphericalTensor > &eqn, const label fieldi)
 
virtual void constrain (fvMatrix< symmTensor > &eqn, const label fieldi)
 
virtual void constrain (fvMatrix< tensor > &eqn, const label fieldi)
 
virtual void correct (volScalarField &field)
 
virtual void correct (volVectorField &field)
 
virtual void correct (volSphericalTensorField &field)
 
virtual void correct (volSymmTensorField &field)
 
virtual void correct (volTensorField &field)
 
virtual void postProcessSens (scalarField &sensField, const word &fieldName=word::null, const word &designVariablesName=word::null)
 
virtual void postProcessSens (vectorField &sensField, const word &fieldName=word::null, const word &designVariablesName=word::null)
 
virtual void postProcessSens (tensorField &sensField, const word &fieldName=word::null, const word &designVariablesName=word::null)
 
virtual void writeHeader (Ostream &) const
 Write the source header information. More...
 
virtual void writeFooter (Ostream &) const
 Write the source footer information. More...
 
virtual void writeData (Ostream &) const
 Write the source properties. More...
 
virtual bool read (const dictionary &dict)
 Read source dictionary. More...
 

Additional Inherited Members

- Static Public Member Functions inherited from option
static autoPtr< optionNew (const word &name, const dictionary &dict, const fvMesh &mesh)
 Return a reference to the selected fvOption model. More...
 
- Public Attributes inherited from option
bool log
 Switch write log to Info. More...
 
- Protected Member Functions inherited from option
void resetApplied ()
 Resize/reset applied flag list for all fieldNames_ entries. More...
 
- Protected Attributes inherited from option
const word name_
 Source name. More...
 
const word modelType_
 Model type. More...
 
const fvMeshmesh_
 Reference to the mesh database. More...
 
dictionary dict_
 Top level source dictionary. More...
 
dictionary coeffs_
 Dictionary containing source coefficients. More...
 
wordList fieldNames_
 Field names to apply source to - populated by derived models. More...
 
List< boolapplied_
 Applied flag list - corresponds to each fieldNames_ entry. More...
 
bool active_
 Source active flag. More...
 

Detailed Description

Source defined by a boundary condition applied to cells next to patches. This fvOption needs to be used with a boundarySourcePatch type of boundary condition (e.g. speciesSorption and enthalpySorption.)

Usage
Minimal example by using constant/fvOptions:
<fvOptionName>
{
    // Mandatory entries
    type        patchCellsSource;

    // Optional entries (one only!)
    U           <word>;
    he          <word>;
    species     <word>;

    // Inherited entries
    ...
}

where the entries mean:

Property Description Type Reqd Deflt
type Type name: patchCellsSource word yes -
U Name of operand velocity field word choice -
he Name of operand energy field word choice -
species Name of operand species field word choice -

The inherited entries are elaborated in:

Source files

Definition at line 119 of file patchCellsSource.H.

Constructor & Destructor Documentation

◆ patchCellsSource() [1/2]

patchCellsSource ( const word sourceName,
const word modelType,
const dictionary dict,
const fvMesh mesh 
)

Construct from explicit source name and mesh.

Definition at line 47 of file patchCellsSource.C.

References option::coeffs_, Foam::endl(), Foam::exit(), Foam::FatalIOError, FatalIOErrorInFunction, option::fieldNames_, dictionary::readIfPresent(), option::resetApplied(), and List< T >::resize().

Here is the call graph for this function:

◆ patchCellsSource() [2/2]

patchCellsSource ( const patchCellsSource )
delete

No copy construct.

◆ ~patchCellsSource()

virtual ~patchCellsSource ( )
virtualdefault

Destructor.

Member Function Documentation

◆ TypeName()

TypeName ( "patchCellsSource"  )

Runtime type information.

◆ operator=()

void operator= ( const patchCellsSource )
delete

No copy assignment.

◆ addSup()

void addSup ( const volScalarField rho,
fvMatrix< scalar > &  eqn,
const label  fieldi 
)
virtual

Add explicit contribution to compressible (momentum, enthalpy, species) equation

Reimplemented from option.

Definition at line 105 of file patchCellsSource.C.

References tmp< T >::cref(), fvMatrix< Type >::dimensions(), Foam::dimVolume, Foam::endl(), forAll, Foam::gMinMax(), Foam::Info, Time::New(), psi, fvMatrix< Type >::psi(), Foam::type(), and Foam::Zero.

Here is the call graph for this function:

◆ read()

bool read ( const dictionary dict)
virtual

Read source dictionary.

Reimplemented from option.

Definition at line 165 of file patchCellsSource.C.

References dict, and kEpsilonLopesdaCosta< BasicTurbulenceModel >::read().

Here is the call graph for this function:

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