svenssonHaggkvistCanopySource.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(svenssonHaggkvistCanopySource, 0);
41  (
42  option,
43  svenssonHaggkvistCanopySource,
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  CpEps1_
61  (
62  dimensionedScalar::lookupOrAddToDict
63  (
64  "CpEps1",
65  coeffs_,
66  1.8
67  )
68  )
69 {
70 }
71 
72 
73 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
74 
75 
76 
78 (
79  fvMatrix<vector>& eqn,
80  const label fieldi
81 )
82 {
83 
84  const volVectorField& U = eqn.psi();
85  const volScalarField& canopy = canopy_;
86 
87  fvMatrix<vector> S_canopy
88  (
89  fvm::Sp(canopy * mag(U), U)
90  );
91 
92  eqn -= S_canopy;
93 
94 }
95 
96 
98  (
99  const volScalarField& rho,
100  fvMatrix<vector>& eqn,
101  const label fieldi
102  )
103  {
104 
105  if (eqn.psi().name() == word("U")) {
106 
107  const volVectorField& U = eqn.psi();
108  const volScalarField& canopy = canopy_;
109 
110  fvMatrix<vector> S_canopy
111  (
112  fvm::Sp(rho*canopy * mag(U), U)
113  );
114 
115  eqn -= S_canopy;
116  }
117 
118 
119  }
120 
121 
123 (
124  fvMatrix<scalar>& eqn,
125  const label fieldi
126 )
127 {
128  const volScalarField& canopy = canopy_;
129 
130  if (eqn.psi().name() == word("k")) {
131 
132  const volVectorField& U = mesh_.lookupObject<volVectorField>("U");
133  eqn += canopy*pow(mag(U),3);
134  }
135  else if (eqn.psi().name() == word("epsilon")) {
136 
137  const volScalarField& epsilon = eqn.psi();
138  const volScalarField& k = mesh_.lookupObject<volScalarField>("k");
139  const volVectorField& U = mesh_.lookupObject<volVectorField>("U");
140 
141 
142  fvMatrix<scalar> Sepsilon
143  (
144  fvm::Sp(canopy*CpEps1_/k*pow(mag(U),3), epsilon)
145  );
146 
147  eqn += Sepsilon;
148  }
149 }
150 
151 
153  (
154  const volScalarField& rho,
155  fvMatrix<scalar>& eqn,
156  const label fieldi
157  )
158  {
159 
160  const volScalarField& canopy = canopy_;
161 
162  if (eqn.psi().name() == word("k")) {
163  const volVectorField& U = mesh_.lookupObject<volVectorField>("U");
164 
165  eqn += rho*canopy*pow(mag(U),3);
166  }
167  else if (eqn.psi().name() == word("epsilon")) {
168 
169  const volScalarField& epsilon = eqn.psi();
170  const volScalarField& k = mesh_.lookupObject<volScalarField>("k");
171  const volVectorField& U = mesh_.lookupObject<volVectorField>("U");
172 
173 
174  fvMatrix<scalar> Sepsilon
175  (
176  fvm::Sp(rho*canopy*CpEps1_/k*pow(mag(U),3), epsilon)
177  );
178 
179 
180  eqn += Sepsilon;
181  }
182 }
183 
184 
185 
186 // ************************************************************************* //
svenssonHaggkvistCanopySource.H
Foam::word
A class for handling words, derived from string.
Definition: word.H:59
Foam::fv::svenssonHaggkvistCanopySource::addSup
void addSup(const RhoFieldType &rho, fvMatrix< vector > &eqn, const label fieldi)
Source term to momentum equation.
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.
Foam::fv::svenssonHaggkvistCanopySource::svenssonHaggkvistCanopySource
svenssonHaggkvistCanopySource(const svenssonHaggkvistCanopySource &)
Disallow default bitwise copy construct.
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::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)
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