Loading...
Searching...
No Matches
PhaseLimitStabilization.C
Go to the documentation of this file.
1/*---------------------------------------------------------------------------*\
2 ========= |
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4 \\ / O peration |
5 \\ / A nd | www.openfoam.com
6 \\/ M anipulation |
7-------------------------------------------------------------------------------
8 Copyright (C) 2017 OpenFOAM Foundation
9 Copyright (C) 2020-2021 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
12 This file is part of OpenFOAM.
13
14 OpenFOAM is free software: you can redistribute it and/or modify it
15 under the terms of the GNU General Public License as published by
16 the Free Software Foundation, either version 3 of the License, or
17 (at your option) any later version.
18
19 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22 for more details.
23
24 You should have received a copy of the GNU General Public License
25 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26
27\*---------------------------------------------------------------------------*/
28
30#include "fvMatrices.H"
31#include "fvmSup.H"
33
34// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
35
36template<class Type>
38(
39 const word& name,
40 const word& modelType,
41 const dictionary& dict,
42 const fvMesh& mesh
43)
44:
45 fv::option(name, modelType, dict, mesh),
46 fieldName_(coeffs_.get<word>("field")),
47 rateName_(coeffs_.get<word>("rate")),
48 residualAlpha_(coeffs_.get<scalar>("residualAlpha"))
49{
50 fieldNames_.resize(1, fieldName_);
52}
53
54
55// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
56
57template<class Type>
59(
60 const volScalarField& alpha,
61 const volScalarField& rho,
62 fvMatrix<Type>& eqn,
63 const label fieldi
64)
65{
67
68 auto& rate =
69 mesh_.lookupObjectRef<uniformDimensionedScalarField>(rateName_);
70
71 eqn -= fvm::Sp(max(residualAlpha_ - alpha, scalar(0))*rho*rate, psi);
72}
73
74
75template<class Type>
77{
79 {
80 coeffs_.readEntry("residualAlpha", residualAlpha_);
81
82 return true;
83 }
84
85 return false;
86}
87
88
89// ************************************************************************* //
Generic GeometricField class.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition fvMatrix.H:118
const GeometricField< Type, fvPatchField, volMesh > & psi(const label i=0) const
Return psi.
Definition fvMatrix.H:487
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
virtual void addSup(const volScalarField &alpha, const volScalarField &rho, fvMatrix< Type > &eqn, const label fieldi)
Source term to compressible phase equation.
PhaseLimitStabilization(const word &name, const word &modelType, const dictionary &dict, const fvMesh &mesh)
Construct from components.
virtual bool read(const dictionary &dict)
Read dictionary.
const word & name() const noexcept
Return const access to the source name.
Definition fvOptionI.H:24
const fvMesh & mesh_
Reference to the mesh database.
Definition fvOption.H:142
wordList fieldNames_
Field names to apply source to - populated by derived models.
Definition fvOption.H:157
option(const word &name, const word &modelType, const dictionary &dict, const fvMesh &mesh)
Construct from components.
Definition fvOption.C:51
virtual bool read(const dictionary &dict)
Read source dictionary.
Definition fvOptionIO.C:48
dictionary coeffs_
Dictionary containing source coefficients.
Definition fvOption.H:152
void resetApplied()
Resize/reset applied flag list for all fieldNames_ entries.
Definition fvOption.C:41
const fvMesh & mesh() const noexcept
Return const access to the mesh database.
Definition fvOptionI.H:30
A class for handling words, derived from Foam::string.
Definition word.H:66
const volScalarField & psi
A special matrix type and solver, designed for finite volume solutions of scalar equations.
Calculate the finiteVolume matrix for implicit and explicit sources.
Namespace for finite-volume.
zeroField Sp(const Foam::zero, const GeometricField< Type, fvPatchField, volMesh > &)
A no-op source.
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:40
GeometricField< scalar, fvPatchField, volMesh > volScalarField
UniformDimensionedField< scalar > uniformDimensionedScalarField
volScalarField & alpha
dictionary dict
Various UniformDimensionedField types.