dalpeMassonCanopySource.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 | Copyright (C) 2015-2016 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 
27 #include "fvMesh.H"
28 #include "fvMatrices.H"
29 #include "fvCFD.H"
31 #include "fvOption.H"
32 
33 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
34 
35 namespace Foam
36 {
37 namespace fv
38 {
39  defineTypeNameAndDebug(dalpeMassonCanopySource, 0);
41  (
42  option,
43  dalpeMassonCanopySource,
45  );
46 }
47 }
48 
49 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
50 
52 (
53  const word& name,
54  const word& modelType,
55  const dictionary& dict,
56  const fvMesh& mesh
57 )
58 :
59  canopySource(name, modelType, dict, mesh),
60  betaP_
61  (
62  dimensionedScalar::lookupOrAddToDict
63  (
64  "betaP",
65  coeffs_,
66  1.0
67  )
68  ),
69  betaD_
70  (
71  dimensionedScalar::lookupOrAddToDict
72  (
73  "betaD",
74  coeffs_,
75  5.03
76  )
77  ),
78  C4_
79  (
80  dimensionedScalar::lookupOrAddToDict
81  (
82  "C4",
83  coeffs_,
84  0.78
85  )
86  ),
87  C5_
88  (
89  dimensionedScalar::lookupOrAddToDict
90  (
91  "C5",
92  coeffs_,
93  0.78
94  )
95  )
96 {
97 }
98 
99 
100 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
101 
102 
104 (
105  fvMatrix<vector>& eqn,
106  const label fieldi
107 )
108 {
109 
110  if (eqn.psi().name() == word("U")) {
111 
112  const volVectorField& U = eqn.psi();
113  const volScalarField& canopy = canopy_;
114 
115  fvMatrix<vector> S_canopy
116  (
117  fvm::Sp(canopy * mag(U), U)
118  );
119 
120  eqn -= S_canopy;
121  }
122 }
123 
124 
126  (
127  const volScalarField& rho,
128  fvMatrix<vector>& eqn,
129  const label fieldi
130  )
131  {
132 
133  if (eqn.psi().name() == word("U")) {
134 
135  const volVectorField& U = eqn.psi();
136  const volScalarField& canopy = canopy_;
137 
138  fvMatrix<vector> S_canopy
139  (
140  fvm::Sp(rho*canopy * mag(U), U)
141  );
142 
143  eqn -= S_canopy;
144  }
145  }
146 
147 
149 (
150  fvMatrix<scalar>& eqn,
151  const label fieldi
152 )
153 {
154  const volScalarField& canopy = canopy_;
155 
156  if (eqn.psi().name() == word("k")) {
157 
158  const volScalarField& k = eqn.psi();
159  const volVectorField& U = mesh_.lookupObject<volVectorField>("U");
160 
162  (
163  betaP_*canopy*pow(mag(U),3) - fvm::Sp(betaD_*canopy*mag(U), k)
164  );
165 
166  eqn += Sk;
167  }
168  else if (eqn.psi().name() == word("epsilon")) {
169 
170  const volScalarField& epsilon = eqn.psi();
171  const volScalarField& k = mesh_.lookupObject<volScalarField>("k");
172  const volVectorField& U = mesh_.lookupObject<volVectorField>("U");
173 
174 
175  fvMatrix<scalar> Sepsilon
176  (
177  fvm::Sp(canopy/k*(C4_*betaP_*pow(mag(U),3) - C5_*betaD_*k*mag(U)), epsilon)
178  );
179 
180  eqn += Sepsilon;
181  // canopy*epsilon/k*(C4_*betaP_*pow(mag(U),3) - C5_*betaD_*k*mag(U));
182  }
183 }
184 
185 
187  (
188  const volScalarField& rho,
189  fvMatrix<scalar>& eqn,
190  const label fieldi
191  )
192  {
193 
194  const volScalarField& canopy = canopy_;
195 
196  if (eqn.psi().name() == word("k")) {
197 
198  const volScalarField& k = eqn.psi();
199  const volVectorField& U = mesh_.lookupObject<volVectorField>("U");
200 
202  (
203  betaP_*rho*canopy*pow(mag(U),3) - fvm::Sp(betaD_*rho*canopy*mag(U), k)
204  );
205 
206  eqn += Sk;
207  }
208  else if (eqn.psi().name() == word("epsilon")) {
209 
210  const volScalarField& epsilon = eqn.psi();
211  const volScalarField& k = mesh_.lookupObject<volScalarField>("k");
212  const volVectorField& U = mesh_.lookupObject<volVectorField>("U");
213 
214 
215  fvMatrix<scalar> Sepsilon
216  (
217  fvm::Sp(rho*canopy/k*(C4_*betaP_*pow(mag(U),3) - C5_*betaD_*k*mag(U)), epsilon)
218  );
219 
220 
221  eqn += Sepsilon;
222  }
223 }
224 
225 
226 
227 // ************************************************************************* //
Foam::word
A class for handling words, derived from string.
Definition: word.H:59
Foam::mag
dimensioned< scalar > mag(const dimensioned< Type > &)
Foam::fvMatrix::psi
const GeometricField< Type, fvPatchField, volMesh > & psi() const
Definition: fvMatrix.H:281
U
U
Definition: pEqn.H:46
fvMatrices.H
A special matrix type and solver, designed for finite volume solutions of scalar equations.
dalpeMassonCanopySource.H
Foam::label
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
Foam::fv::canopySource
Base class for momentum and turbulence source/sink-terms for tree canopy.
Definition: canopySource.H:119
Foam::fv::dalpeMassonCanopySource::addSup
void addSup(const RhoFieldType &rho, fvMatrix< vector > &eqn, const label fieldi)
Source term to momentum equation.
Foam::pow
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
Definition: dimensionedScalar.C:73
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::dictionary
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:137
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:18
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
fvMesh.H
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
rho
rho
Definition: pEqn.H:3
fv
labelList fv(nPoints)
Foam::fv::defineTypeNameAndDebug
defineTypeNameAndDebug(option, 0)
Foam::fv::dalpeMassonCanopySource::dalpeMassonCanopySource
dalpeMassonCanopySource(const dalpeMassonCanopySource &)
Disallow default bitwise copy construct.
epsilon
epsilon
Definition: createTDFields.H:56
Foam::fv::addToRunTimeSelectionTable
addToRunTimeSelectionTable(option, fixedTemperatureConstraint, dictionary)
fvOption.H
k
label k
Boltzmann constant.
Definition: LISASMDCalcMethod2.H:41
Foam::fvMatrix
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition: fvPatchField.H:68
fvCFD.H
Foam::fvc::Sp
tmp< GeometricField< Type, fvPatchField, volMesh > > Sp(const volScalarField &sp, const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: fvcSup.C:67
Foam::GeometricField
Generic GeometricField class.
Definition: surfaceFieldsFwd.H:52
Foam::name
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47