waveForcing Class Reference

This fvModel applies forcing to the liquid phase-fraction field and all components of the vector field to relax the fields towards those calculated from the current wave distribution. More...

Inheritance diagram for waveForcing:
Collaboration diagram for waveForcing:

Public Member Functions

 TypeName ("waveForcing")
 Runtime type information. More...
 
 waveForcing (const word &name, const word &modelType, const fvMesh &mesh, const dictionary &dict)
 Construct from components. More...
 
virtual ~waveForcing ()
 Destructor. More...
 
virtual wordList addSupFields () const
 Return the list of fields for which the fvModel adds source term. More...
 
virtual void addSup (fvMatrix< scalar > &eqn, const word &fieldName) const
 Source term to VoF phase-fraction equation. More...
 
virtual void addSup (const volScalarField &rho, fvMatrix< vector > &eqn, const word &fieldName) const
 Source term to momentum equation. More...
 
virtual bool movePoints ()
 Update for mesh motion. More...
 
virtual void topoChange (const polyTopoChangeMap &)
 Update topology using the given map. More...
 
virtual void mapMesh (const polyMeshMap &)
 Update from another mesh using the given map. More...
 
virtual void distribute (const polyDistributionMap &)
 Redistribute or update using the given distribution map. More...
 
virtual void correct ()
 Correct the wave forcing coefficients. More...
 
- Public Member Functions inherited from forcing
 TypeName ("forcing")
 Runtime type information. More...
 
 forcing (const word &name, const word &modelType, const fvMesh &mesh, const dictionary &dict)
 Construct from components. More...
 
virtual ~forcing ()
 Destructor. More...
 
virtual bool read (const dictionary &dict)
 Read dictionary. More...
 
- Public Member Functions inherited from fvModel
 TypeName ("fvModel")
 Runtime type information. More...
 
 declareRunTimeSelectionTable (autoPtr, fvModel, dictionary,(const word &name, const word &modelType, const fvMesh &mesh, const dictionary &dict),(name, modelType, mesh, dict))
 
 fvModel (const word &name, const word &modelType, const fvMesh &mesh, const dictionary &dict)
 Construct from components. More...
 
autoPtr< fvModelclone () const
 Return clone. More...
 
virtual ~fvModel ()
 Destructor. More...
 
const wordname () const
 Return const access to the source name. More...
 
const fvMeshmesh () const
 Return const access to the mesh database. More...
 
const dictionarycoeffs () const
 Return dictionary. More...
 
virtual bool addsSupToField (const word &fieldName) const
 Return true if the fvModel adds a source term to the given. More...
 
virtual scalar maxDeltaT () const
 Return the maximum time-step for stable operation. More...
 
 FOR_ALL_FIELD_TYPES (DEFINE_FV_MODEL_ADD_SUP)
 Add a source term to an equation. More...
 
 FOR_ALL_FIELD_TYPES (DEFINE_FV_MODEL_ADD_RHO_SUP)
 Add a source term to a compressible equation. More...
 
 FOR_ALL_FIELD_TYPES (DEFINE_FV_MODEL_ADD_ALPHA_RHO_SUP)
 Add a source term to a phase equation. More...
 
template<class Type >
tmp< fvMatrix< Type > > source (const VolField< Type > &field) const
 Return source for an equation. More...
 
template<class Type >
tmp< fvMatrix< Type > > source (const VolField< Type > &field, const word &fieldName) const
 Return source for an equation with a specified name. More...
 
template<class Type >
tmp< fvMatrix< Type > > source (const volScalarField &rho, const VolField< Type > &field) const
 Return source for a compressible equation. More...
 
template<class Type >
tmp< fvMatrix< Type > > source (const volScalarField &rho, const VolField< Type > &field, const word &fieldName) const
 Return source for a compressible equation with a specified name. More...
 
template<class Type >
tmp< fvMatrix< Type > > source (const volScalarField &alpha, const volScalarField &rho, const VolField< Type > &field) const
 Return source for a phase equation. More...
 
template<class Type >
tmp< fvMatrix< Type > > source (const volScalarField &alpha, const volScalarField &rho, const VolField< Type > &field, const word &fieldName) const
 Return source for a phase equation with a specified name. More...
 
template<class Type >
tmp< fvMatrix< Type > > source (const volScalarField &alpha, const geometricOneField &rho, const VolField< Type > &field) const
 Return source for a phase equation. More...
 
template<class Type >
tmp< fvMatrix< Type > > source (const geometricOneField &alpha, const volScalarField &rho, const VolField< Type > &field) const
 Return source for a phase equation. More...
 
template<class Type >
tmp< fvMatrix< Type > > source (const geometricOneField &alpha, const geometricOneField &rho, const VolField< Type > &field) const
 Return source for a phase equation. More...
 
template<class Type >
tmp< fvMatrix< Type > > d2dt2 (const VolField< Type > &field) const
 Return source for an equation with a second time derivative. More...
 
template<class Type >
tmp< fvMatrix< Type > > d2dt2 (const VolField< Type > &field, const word &fieldName) const
 Return source for an equation with a second time derivative. More...
 
virtual void preUpdateMesh ()
 Prepare for mesh update. More...
 
template<class Type , class AlphaRhoFieldType , class ... AlphaRhoFieldTypes>
Foam::dimensionSet sourceDims (const VolField< Type > &field, const dimensionSet &ds, const AlphaRhoFieldType &alphaRho, const AlphaRhoFieldTypes &... alphaRhos)
 
template<class Type >
Foam::dimensionSet sourceDims (const VolField< Type > &field, const dimensionSet &ds)
 
template<class Type , class ... AlphaRhoFieldTypes>
Foam::tmp< Foam::fvMatrix< Type > > source (const VolField< Type > &field, const word &fieldName, const dimensionSet &ds, const AlphaRhoFieldTypes &... alphaRhos) const
 
template<class Type >
Foam::tmp< Foam::fvMatrix< Type > > source (const VolField< Type > &field) const
 
template<class Type >
Foam::tmp< Foam::fvMatrix< Type > > source (const VolField< Type > &field, const word &fieldName) const
 
template<class Type >
Foam::tmp< Foam::fvMatrix< Type > > source (const volScalarField &rho, const VolField< Type > &field) const
 
template<class Type >
Foam::tmp< Foam::fvMatrix< Type > > source (const volScalarField &rho, const VolField< Type > &field, const word &fieldName) const
 
template<class Type >
Foam::tmp< Foam::fvMatrix< Type > > source (const volScalarField &alpha, const volScalarField &rho, const VolField< Type > &field) const
 
template<class Type >
Foam::tmp< Foam::fvMatrix< Type > > source (const volScalarField &alpha, const volScalarField &rho, const VolField< Type > &field, const word &fieldName) const
 
template<class Type >
Foam::tmp< Foam::fvMatrix< Type > > source (const geometricOneField &alpha, const geometricOneField &rho, const VolField< Type > &field) const
 
template<class Type >
Foam::tmp< Foam::fvMatrix< Type > > source (const volScalarField &alpha, const geometricOneField &rho, const VolField< Type > &field) const
 
template<class Type >
Foam::tmp< Foam::fvMatrix< Type > > source (const geometricOneField &alpha, const volScalarField &rho, const VolField< Type > &field) const
 
template<class Type >
Foam::tmp< Foam::fvMatrix< Type > > d2dt2 (const VolField< Type > &field) const
 
template<class Type >
Foam::tmp< Foam::fvMatrix< Type > > d2dt2 (const VolField< Type > &field, const word &fieldName) const
 

Additional Inherited Members

- Static Public Member Functions inherited from fvModel
template<class Type , class AlphaRhoFieldType , class ... AlphaRhoFieldTypes>
static dimensionSet sourceDims (const VolField< Type > &field, const dimensionSet &ds, const AlphaRhoFieldType &alphaRho, const AlphaRhoFieldTypes &... alphaRhos)
 Return the dimensions of the matrix of a source term. More...
 
template<class Type >
static dimensionSet sourceDims (const VolField< Type > &field, const dimensionSet &ds)
 Return the dimensions of the matrix of a source term (base. More...
 
static autoPtr< fvModelNew (const word &name, const fvMesh &mesh, const dictionary &dict)
 Return a reference to the selected fvModel. More...
 
- Protected Member Functions inherited from forcing
void readCoeffs ()
 Non-virtual read. More...
 
tmp< volScalarField::Internalscale () const
 Return the scale distribution. More...
 
tmp< volScalarField::InternalforceCoeff (const volScalarField::Internal &scale) const
 Return the force coefficient given the scale distribution. More...
 
tmp< volScalarField::InternalforceCoeff () const
 Return the force coefficient. More...
 
- Protected Member Functions inherited from fvModel
template<class Type >
void addSupType (fvMatrix< Type > &eqn, const word &fieldName) const
 Add a source term to an equation. More...
 
template<class Type >
void addSupType (const volScalarField &rho, fvMatrix< Type > &eqn, const word &fieldName) const
 Add a source term to a compressible equation. More...
 
template<class Type >
void addSupType (const volScalarField &alpha, const volScalarField &rho, fvMatrix< Type > &eqn, const word &fieldName) const
 Add a source term to a phase equation. More...
 
template<class Type , class ... AlphaRhoFieldTypes>
tmp< fvMatrix< Type > > source (const VolField< Type > &field, const word &fieldName, const dimensionSet &ds, const AlphaRhoFieldTypes &... alphaRhos) const
 Return source for equation with specified name and dimensions. More...
 
- Protected Attributes inherited from forcing
dimensionedScalar lambda_
 Damping coefficient [1/s]. More...
 
dimensionedScalar lambdaBoundary_
 Optional boundary damping coefficient [1/s]. More...
 
autoPtr< Function1< scalar > > scale_
 The scaling function. More...
 
vectorField origins_
 Origins of the scaling coordinate. More...
 
vectorField directions_
 Directions of increasing scaling coordinate. More...
 

Detailed Description

This fvModel applies forcing to the liquid phase-fraction field and all components of the vector field to relax the fields towards those calculated from the current wave distribution.

The forcing force coefficient $\lambda$ should be set based on the desired level of forcing and the residence time the waves through the forcing zone. For example, if waves moving at 2 [m/s] are travelling through a forcing zone 8 [m] in length, then the residence time is 4 [s]. If it is deemed necessary to force for 5 time-scales, then $\lambda$ should be set to equal 5/(4 [s]) = 1.2 [1/s].

Usage
Example usage:
waveForcing1
{
    type            waveForcing;

    libs            ("libwaves.so");

    liquidPhase     water;

    // Define the line along which to apply the graduation
    origin          (600 0 0);
    direction       (-1 0 0);

    // // Or, define multiple lines
    // origins         ((600 0 0) (600 -300 0) (600 300 0));
    // directions      ((-1 0 0) (0 1 0) (0 -1 0));

    scale
    {
        type        halfCosineRamp;
        start       0;
        duration    300;
    }

    lambda          0.5; // Forcing coefficient

 // lambda          2;   // Optional boundary forcing coefficient
                         // Useful when wave BCs are specified
                         // without mean flow
}
See also
Foam::fv::forcing
Source files

Definition at line 97 of file waveForcing.H.

Constructor & Destructor Documentation

◆ waveForcing()

waveForcing ( const word name,
const word modelType,
const fvMesh mesh,
const dictionary dict 
)

Construct from components.

Definition at line 49 of file waveForcing.C.

◆ ~waveForcing()

virtual ~waveForcing ( )
inlinevirtual

Destructor.

Definition at line 144 of file waveForcing.H.

Member Function Documentation

◆ TypeName()

TypeName ( "waveForcing"  )

Runtime type information.

◆ addSupFields()

Foam::wordList addSupFields ( ) const
virtual

Return the list of fields for which the fvModel adds source term.

to the transport equation

Reimplemented from fvModel.

Definition at line 68 of file waveForcing.C.

◆ addSup() [1/2]

void addSup ( fvMatrix< scalar > &  eqn,
const word fieldName 
) const
virtual

Source term to VoF phase-fraction equation.

Definition at line 74 of file waveForcing.C.

References fvMatrix< Type >::psi(), and Foam::fvm::Sp().

Here is the call graph for this function:

◆ addSup() [2/2]

void addSup ( const volScalarField rho,
fvMatrix< vector > &  eqn,
const word fieldName 
) const
virtual

Source term to momentum equation.

Definition at line 90 of file waveForcing.C.

References Foam::fvc::ddt(), Foam::fvc::div(), fvMatrix< Type >::psi(), rho, and Foam::fvm::Sp().

Here is the call graph for this function:

◆ movePoints()

bool movePoints ( )
virtual

Update for mesh motion.

Implements fvModel.

Definition at line 116 of file waveForcing.C.

◆ topoChange()

void topoChange ( const polyTopoChangeMap )
virtual

Update topology using the given map.

Implements fvModel.

Definition at line 123 of file waveForcing.C.

◆ mapMesh()

void mapMesh ( const polyMeshMap map)
virtual

Update from another mesh using the given map.

Implements fvModel.

Definition at line 129 of file waveForcing.C.

◆ distribute()

void distribute ( const polyDistributionMap )
virtual

Redistribute or update using the given distribution map.

Implements fvModel.

Definition at line 135 of file waveForcing.C.

◆ correct()

void correct ( )
virtual

Correct the wave forcing coefficients.

Reimplemented from fvModel.

Definition at line 141 of file waveForcing.C.

References Foam::dimless, Foam::dimVelocity, Foam::constant::universal::h, Foam::levelSetAverage(), Foam::levelSetFraction(), and DimensionedField< Type, GeoMesh >::New().

Here is the call graph for this function:

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