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-2023 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 {
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 fvMesh& mesh,
61  const dictionary& dict
62 )
63 :
64  fvModel(name, modelType, mesh, dict),
65  set_(mesh, coeffs()),
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_.movePoints();
117  return true;
118 }
119 
120 
122 {
123  set_.topoChange(map);
124 }
125 
126 
128 {
129  set_.mapMesh(map);
130 }
131 
132 
134 (
135  const polyDistributionMap& map
136 )
137 {
138  set_.distribute(map);
139 }
140 
141 
143 {
144  if (fvModel::read(dict))
145  {
146  set_.read(coeffs());
147  readCoeffs();
148  return true;
149  }
150  else
151  {
152  return false;
153  }
154 }
155 
156 
157 // ************************************************************************* //
Macros for easy insertion into run-time selection tables.
static autoPtr< Function1< Type > > New(const word &name, const dictionary &dict)
Selector.
Definition: Function1New.C:32
Generic GeometricField class.
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:160
T lookupOrDefault(const word &, const T &, bool recursive=false, bool patternMatch=true) const
Find and return a T,.
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition: fvMatrix.H:118
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:101
Finite volume model abstract base class.
Definition: fvModel.H:59
const dictionary & coeffs() const
Return dictionary.
Definition: fvModelI.H:40
virtual bool read(const dictionary &dict)
Read source dictionary.
Definition: fvModel.C:187
This fvModel applies an explicit acceleration force to components of the velocity field.
virtual bool movePoints()
Update for mesh motion.
virtual wordList addSupFields() const
Return the list of fields for which the fvModel adds source term.
virtual void topoChange(const polyTopoChangeMap &)
Update topology using the given map.
virtual void distribute(const polyDistributionMap &)
Redistribute or update using the given distribution map.
virtual bool read(const dictionary &dict)
Read dictionary.
virtual void addSup(fvMatrix< vector > &eqn, const word &fieldName) const
Source term to momentum equation.
virtual void mapMesh(const polyMeshMap &)
Update from another mesh using the given map.
accelerationSource(const word &name, const word &modelType, const fvMesh &mesh, const dictionary &dict)
Construct from components.
A class representing the concept of a GeometricField of 1 used to avoid unnecessary manipulations for...
Class containing mesh-to-mesh mapping information after a mesh distribution where we send parts of me...
Class containing mesh-to-mesh mapping information.
Definition: polyMeshMap.H:51
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
A class for handling words, derived from string.
Definition: word.H:62
volScalarField alpha(IOobject("alpha", runTime.name(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
addToRunTimeSelectionTable(fvConstraint, bound, dictionary)
defineTypeNameAndDebug(bound, 0)
Namespace for OpenFOAM.
List< word > wordList
A List of words.
Definition: fileName.H:54
void add(FieldField< Field1, typename typeOfSum< Type1, Type2 >::type > &f, const FieldField< Field1, Type1 > &f1, const FieldField< Field2, Type2 > &f2)
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
labelList fv(nPoints)
dictionary dict