atmCoriolisUSource.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 CENER
9 Copyright (C) 2020-2021 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
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 "atmCoriolisUSource.H"
30#include "fvMatrices.H"
31#include "unitConversion.H"
34
35using namespace Foam::constant::mathematical;
36
37// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
38
39namespace Foam
40{
41namespace fv
42{
45}
46}
47
48// * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * * //
49
50Foam::vector Foam::fv::atmCoriolisUSource::planetaryRotationVector() const
51{
52 return vector
53 (
54 Zero,
55 twoPi/(planetaryRotationPeriod_*3600.0)*cos(degToRad(latitude_)),
56 twoPi/(planetaryRotationPeriod_*3600.0)*sin(degToRad(latitude_))
57 );
58}
59
60
61// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
62
64(
65 const word& sourceName,
66 const word& modelType,
67 const dictionary& dict,
68 const fvMesh& mesh
69)
70:
71 fv::cellSetOption(sourceName, modelType, dict, mesh),
72 latitude_
73 (
74 coeffs_.getCheckOrDefault<scalar>
75 (
76 "latitude",
77 0.0,
78 [&](const scalar x){ return (90 >= mag(x)) && (mag(x) >= 0); }
79 )
80 ),
81 planetaryRotationPeriod_
82 (
83 coeffs_.getCheckOrDefault<scalar>
84 (
85 "planetaryRotationPeriod",
86 23.9344694,
87 [&](const scalar x){ return x > SMALL; }
88 )
89 ),
90 Omega_
91 (
93 (
95 coeffs_.getOrDefault<vector>
96 (
97 "Omega",
98 planetaryRotationVector()
99 )
100 )
101 )
102{
103 if (mag(Omega_.value()) < SMALL)
104 {
106 << "The magnitude of the rotation vector in atmCoriolisUSource is "
107 << "effectively zero, mag(Omega) = " << mag(Omega_.value()) << nl
108 << "Please check input values in atmCoriolisUSource settings."
109 << endl;
110 }
111
112 fieldNames_.resize(1, "U");
113
115
116 Log << " Applying atmCoriolisUSource to: " << fieldNames_[0] << endl;
117}
118
119
120// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
121
123(
124 fvMatrix<vector>& eqn,
125 const label fieldi
126)
127{
128 const volVectorField& U = eqn.psi();
129
130 if (V_ > VSMALL)
131 {
132 eqn -= (2.0*Omega_)^U;
133 }
134}
135
136
138(
139 const volScalarField& rho,
140 fvMatrix<vector>& eqn,
141 const label fieldi
142)
143{
144 const volVectorField& U = eqn.psi();
145
146 if (V_ > VSMALL)
147 {
148 eqn -= rho*((2.0*Omega_)^U);
149 }
150}
151
152
154(
155 const volScalarField& alpha,
156 const volScalarField& rho,
157 fvMatrix<vector>& eqn,
158 const label fieldi
159)
160{
161 const volVectorField& U = eqn.psi();
162
163 if (V_ > VSMALL)
164 {
165 eqn -= alpha*rho*((2.0*Omega_)^U);
166 }
167}
168
169
170// ************************************************************************* //
#define Log
Definition: PDRblock.C:35
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition: fvMatrix.H:121
const GeometricField< Type, fvPatchField, volMesh > & psi(const label i=0) const
Return psi.
Definition: fvMatrix.H:412
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:91
Applies corrections to incorporate the horizontal and vertical components of the Coriolis force for w...
virtual void addSup(fvMatrix< vector > &eqn, const label fieldi)
Add explicit contribution to incompressible momentum equation.
Intermediate abstract class for handling cell-set options for the derived fvOptions.
Base abstract class for handling finite volume options (i.e. fvOption).
Definition: fvOption.H:127
void resetApplied()
Resize/reset applied flag list for all fieldNames_ entries.
Definition: fvOption.C:48
A class for handling words, derived from Foam::string.
Definition: word.H:68
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition: className.H:121
U
Definition: pEqn.H:72
dynamicFvMesh & mesh
A special matrix type and solver, designed for finite volume solutions of scalar equations.
#define WarningInFunction
Report a warning using Foam::Warning.
Mathematical constants.
constexpr scalar twoPi(2 *M_PI)
Namespace for OpenFOAM.
const dimensionSet dimless
Dimensionless.
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:53
dimensionedScalar sin(const dimensionedScalar &ds)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
constexpr scalar degToRad() noexcept
Multiplication factor for degrees to radians conversion.
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
Vector< scalar > vector
Definition: vector.H:61
dimensionedScalar cos(const dimensionedScalar &ds)
constexpr char nl
The newline '\n' character (0x0a)
Definition: Ostream.H:53
labelList fv(nPoints)
volScalarField & alpha
dictionary dict
Unit conversion functions.