atmPlantCanopyUSource.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) 2020 ENERCON GmbH
9  Copyright (C) 2020-2022 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
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 \*---------------------------------------------------------------------------*/
28 
29 #include "atmPlantCanopyUSource.H"
31 #include "fvmSup.H"
32 
33 // * * * * * * * * * * * * * Static Member Functions * * * * * * * * * * * * //
34 
35 namespace Foam
36 {
37 namespace fv
38 {
39  defineTypeNameAndDebug(atmPlantCanopyUSource, 0);
40  addToRunTimeSelectionTable(option, atmPlantCanopyUSource, dictionary);
41 }
42 }
43 
44 
45 // * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * * //
46 
47 Foam::volScalarField& Foam::fv::atmPlantCanopyUSource::getOrReadField
48 (
49  const word& fieldName
50 ) const
51 {
52  auto* ptr = mesh_.getObjectPtr<volScalarField>(fieldName);
53 
54  if (!ptr)
55  {
56  ptr = new volScalarField
57  (
58  IOobject
59  (
60  fieldName,
61  mesh_.time().timeName(),
62  mesh_,
66  ),
67  mesh_
68  );
69  mesh_.objectRegistry::store(ptr);
70  }
71 
72  return *ptr;
73 }
74 
75 
76 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
77 
79 (
80  const word& sourceName,
81  const word& modelType,
82  const dictionary& dict,
83  const fvMesh& mesh
84 )
85 :
86  fv::cellSetOption(sourceName, modelType, dict, mesh),
87  CdName_(),
88  LADname_()
89 {
90  read(dict);
91 
92  fieldNames_.resize(1, "U");
93 
95 
96  Log << " Applying atmPlantCanopyUSource to: " << fieldNames_[0] << endl;
97 }
98 
99 
100 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
101 
103 (
104  fvMatrix<vector>& eqn,
105  const label fieldi
106 )
107 {
108  if (V_ < VSMALL)
109  {
110  return;
111  }
112 
113  const volVectorField& U = eqn.psi();
114  const volScalarField& Cd = getOrReadField(CdName_);
115  const volScalarField& LAD = getOrReadField(LADname_);
117  // (SP:Eq. 42), (BSG:Eq. 7)
118  eqn -= fvm::Sp(Cd*LAD*mag(U), U);
119 }
120 
121 
123 (
124  const volScalarField& rho,
125  fvMatrix<vector>& eqn,
126  const label fieldi
127 )
128 {
129  if (V_ < VSMALL)
130  {
131  return;
132  }
133 
134  const volVectorField& U = eqn.psi();
135  const volScalarField& Cd = getOrReadField(CdName_);
136  const volScalarField& LAD = getOrReadField(LADname_);
137 
138  eqn -= fvm::Sp(rho*Cd*LAD*mag(U), U);
139 }
140 
141 
143 (
144  const volScalarField& alpha,
145  const volScalarField& rho,
146  fvMatrix<vector>& eqn,
147  const label fieldi
148 )
149 {
150  if (V_ < VSMALL)
151  {
152  return;
153  }
154 
155  const volVectorField& U = eqn.psi();
156  const volScalarField& Cd = getOrReadField(CdName_);
157  const volScalarField& LAD = getOrReadField(LADname_);
158 
159  eqn -= fvm::Sp(alpha*rho*Cd*LAD*mag(U), U);
160 }
161 
162 
164 {
166  {
167  return false;
168  }
169 
170  CdName_ = dict.getOrDefault<word>("Cd", "Cd");
171  LADname_ = dict.getOrDefault<word>("LAD", "LAD");
172 
173  (void) getOrReadField(CdName_);
174  (void) getOrReadField(LADname_);
175 
176  return true;
177 }
178 
179 
180 // ************************************************************************* //
dictionary dict
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
wordList fieldNames_
Field names to apply source to - populated by derived models.
Definition: fvOption.H:157
void resize(const label len)
Adjust allocated size of list.
Definition: ListI.H:160
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
defineTypeNameAndDebug(atmAmbientTurbSource, 0)
const fvMesh & mesh_
Reference to the mesh database.
Definition: fvOption.H:142
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
atmPlantCanopyUSource(const word &sourceName, const word &modelType, const dictionary &dict, const fvMesh &mesh)
Construct from explicit source name and mesh.
virtual bool read(const dictionary &dict)
Read source dictionary.
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:360
Macros for easy insertion into run-time selection tables.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:81
const GeometricField< Type, fvPatchField, volMesh > & psi(const label i=0) const
Return psi.
Definition: fvMatrix.H:485
dynamicFvMesh & mesh
A class for handling words, derived from Foam::string.
Definition: word.H:63
Type * getObjectPtr(const word &name, const bool recursive=false) const
Return non-const pointer to the object of the given Type, using a const-cast to have it behave like a...
labelList fv(nPoints)
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:64
zeroField Sp(const Foam::zero, const GeometricField< Type, fvPatchField, volMesh > &)
A no-op source.
virtual void addSup(fvMatrix< vector > &eqn, const label fieldi)
Add explicit contribution to momentum equation.
addToRunTimeSelectionTable(option, atmAmbientTurbSource, dictionary)
static word timeName(const scalar t, const int precision=precision_)
Return a time name for the given scalar time value formatted with the given precision.
Definition: Time.C:714
virtual bool read(const dictionary &dict)
Read source dictionary.
U
Definition: pEqn.H:72
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
void resetApplied()
Resize/reset applied flag list for all fieldNames_ entries.
Definition: fvOption.C:41
#define Log
Definition: PDRblock.C:28
Automatically write from objectRegistry::writeObject()
Intermediate abstract class for handling cell-set options for the derived fvOptions.
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
Request registration (bool: true)
Calculate the finiteVolume matrix for implicit and explicit sources.
Namespace for OpenFOAM.