Loading...
Searching...
No Matches
sizeGroup.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-2018 OpenFOAM Foundation
9 Copyright (C) 2020 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\*---------------------------------------------------------------------------*/
29#include "sizeGroup.H"
31#include "mixedFvPatchField.H"
32
33// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
34
36(
37 const word& name,
38 const dictionary& dict,
39 const phaseModel& phase,
41 const fvMesh& mesh
42)
43:
45 (
47 (
49 (
50 name,
52 (
54 velocityGroup.popBalName()
55 )
56 ),
57 mesh.time().timeName(),
58 mesh,
61 ),
62 mesh,
63 dimensionedScalar(name, dimless, dict.get<scalar>("value")),
64 velocityGroup.f().boundaryField().types()
65 ),
66 phase_(phase),
67 velocityGroup_(velocityGroup),
68 d_("d", dimLength, dict),
69 x_("x", velocityGroup.formFactor()*pow3(d_)),
70 value_(dict.get<scalar>("value"))
71{
72 // Adjust refValue at mixedFvPatchField boundaries
73 forAll(this->boundaryField(), patchi)
74 {
75 typedef mixedFvPatchField<scalar> mixedFvPatchScalarField;
76
77 if
78 (
80 )
81 {
82 mixedFvPatchScalarField& f =
84 (
85 this->boundaryFieldRef()[patchi]
86 );
87
88 f.refValue() = value_;
89 }
90 }
91}
92
93
94// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
95
97{}
98
99
100// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
101
104{
106 return nullptr;
107}
108
109
110// ************************************************************************* //
const Mesh & mesh() const noexcept
Return const reference to mesh.
const Boundary & boundaryField() const noexcept
@ READ_IF_PRESENT
Reading is optional [identical to LAZY_READ].
@ AUTO_WRITE
Automatically write from objectRegistry::writeObject().
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition IOobject.H:191
const Time & time() const noexcept
Return Time associated with the objectRegistry.
Definition IOobject.C:456
const word & name() const noexcept
Return the object name.
Definition IOobjectI.H:205
static word groupName(StringType base, const word &group)
Create dot-delimited name.group string.
bool get(const label i) const
Definition UList.H:868
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition autoPtr.H:65
sizeGroup(const word &name, const dictionary &dict, const phaseModel &phase, const velocityGroup &velocityGroup, const fvMesh &mesh)
Definition sizeGroup.C:29
virtual ~sizeGroup()
Destructor.
Definition sizeGroup.C:89
const phaseModel & phase() const
Return const-reference to the phase.
Definition sizeGroupI.H:31
autoPtr< sizeGroup > clone() const
Return clone.
Definition sizeGroup.C:96
This diameterModel is intended for use with a populationBalanceModel in order to simulate polydispers...
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
This boundary condition provides a base class for 'mixed' type boundary conditions,...
Single incompressible phase derived from the phase-fraction. Used as part of the multiPhaseMixture fo...
Definition phaseModel.H:56
A class for handling words, derived from Foam::string.
Definition word.H:66
#define NotImplemented
Issue a FatalErrorIn for a function not currently implemented.
Definition error.H:688
word timeName
Definition getTimeIndex.H:3
Type & refCast(U &obj)
A dynamic_cast (for references) to Type reference.
Definition typeInfo.H:172
const dimensionSet dimless
Dimensionless.
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
dimensionedScalar pow3(const dimensionedScalar &ds)
GeometricField< scalar, fvPatchField, volMesh > volScalarField
const Type * isA(const U &obj)
Attempt dynamic_cast to Type.
Definition typeInfo.H:87
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
labelList f(nPoints)
dictionary dict
Us boundaryFieldRef().evaluateCoupled< coupledFaPatch >()
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299