accelerationSource.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2018-2021 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 #include "fvMesh.H"
27 #include "fvMatrix.H"
28 #include "geometricOneField.H"
29 #include "accelerationSource.H"
31 
32 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
33 
34 namespace Foam
35 {
36 namespace fv
37 {
38  defineTypeNameAndDebug(accelerationSource, 0);
39  addToRunTimeSelectionTable(fvModel, accelerationSource, dictionary);
40 }
41 }
42 
43 
44 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
45 
46 void Foam::fv::accelerationSource::readCoeffs()
47 {
48  UName_ = coeffs().lookupOrDefault<word>("U", "U");
49 
50  velocity_ = Function1<vector>::New("velocity", coeffs());
51 }
52 
53 
54 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
55 
57 (
58  const word& name,
59  const word& modelType,
60  const dictionary& dict,
61  const fvMesh& mesh
62 )
63 :
64  fvModel(name, modelType, dict, mesh),
65  set_(coeffs(), mesh),
66  UName_(word::null),
67  velocity_(nullptr)
68 {
69  readCoeffs();
70 }
71 
72 
73 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
74 
76 {
77  return wordList(1, UName_);
78 }
79 
80 
82 (
83  fvMatrix<vector>& eqn,
84  const word& fieldName
85 ) const
86 {
87  add(geometricOneField(), eqn, fieldName);
88 }
89 
90 
92 (
93  const volScalarField& rho,
94  fvMatrix<vector>& eqn,
95  const word& fieldName
96 ) const
97 {
98  add(rho, eqn, fieldName);
99 }
100 
101 
103 (
104  const volScalarField& alpha,
105  const volScalarField& rho,
106  fvMatrix<vector>& eqn,
107  const word& fieldName
108 ) const
109 {
110  add((alpha*rho)(), eqn, fieldName);
111 }
112 
113 
115 {
116  set_.updateMesh(mpm);
117 }
118 
119 
121 {
122  if (fvModel::read(dict))
123  {
124  set_.read(coeffs());
125  readCoeffs();
126  return true;
127  }
128  else
129  {
130  return false;
131  }
132 }
133 
134 
135 // ************************************************************************* //
defineTypeNameAndDebug(fixedTemperatureConstraint, 0)
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
virtual void updateMesh(const mapPolyMesh &)
Update for mesh changes.
virtual wordList addSupFields() const
Return the list of fields for which the fvModel adds source term.
virtual bool read(const dictionary &dict)
Read source dictionary.
Definition: fvModel.C:158
virtual void addSup(fvMatrix< vector > &eqn, const word &fieldName) const
Source term to momentum equation.
Finite volume model abstract base class.
Definition: fvModel.H:55
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
Definition: mapPolyMesh.H:158
Macros for easy insertion into run-time selection tables.
accelerationSource(const word &name, const word &modelType, const dictionary &dict, const fvMesh &mesh)
Construct from components.
A class representing the concept of a GeometricField of 1 used to avoid unnecessary manipulations for...
A class for handling words, derived from string.
Definition: word.H:59
labelList fv(nPoints)
static const word null
An empty word.
Definition: word.H:77
A special matrix type and solver, designed for finite volume solutions of scalar equations. Face addressing is used to make all matrix assembly and solution loops vectorise.
Definition: fvPatchField.H:72
void add(FieldField< Field1, typename typeOfSum< Type1, Type2 >::type > &f, const FieldField< Field1, Type1 > &f1, const FieldField< Field2, Type2 > &f2)
addToRunTimeSelectionTable(fvConstraint, fixedTemperatureConstraint, dictionary)
List< word > wordList
A List of words.
Definition: fileName.H:54
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
virtual bool read(const dictionary &dict)
Read dictionary.
Namespace for OpenFOAM.
static autoPtr< Function1< Type > > New(const word &name, const dictionary &dict)
Selector.
Definition: Function1New.C:32