BlobsSheetAtomization.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) 2011-2013 OpenFOAM Foundation
9-------------------------------------------------------------------------------
10License
11 This file is part of OpenFOAM.
12
13 OpenFOAM is free software: you can redistribute it and/or modify it
14 under the terms of the GNU General Public License as published by
15 the Free Software Foundation, either version 3 of the License, or
16 (at your option) any later version.
17
18 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21 for more details.
22
23 You should have received a copy of the GNU General Public License
24 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25
26\*---------------------------------------------------------------------------*/
27
29
30// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
31
32template<class CloudType>
34(
35 const dictionary& dict,
36 CloudType& owner
37)
38:
39 AtomizationModel<CloudType>(dict, owner, typeName),
40 B_(this->coeffDict().getScalar("B")),
41 angle_(this->coeffDict().getScalar("angle"))
42{}
43
44
45template<class CloudType>
47(
49)
50:
52 B_(am.B_),
53 angle_(am.angle_)
54{}
55
56
57// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
58
59template<class CloudType>
61{}
62
63
64// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
65
66template<class CloudType>
68{
69 return 1.0;
70}
71
72
73template<class CloudType>
75{
76 return false;
77}
78
79
80template<class CloudType>
82(
83 const scalar dt,
84 scalar& d,
85 scalar& liquidCore,
86 scalar& tc,
87 const scalar rho,
88 const scalar mu,
89 const scalar sigma,
90 const scalar volFlowRate,
91 const scalar rhoAv,
92 const scalar Urel,
93 const vector& pos,
94 const vector& injectionPos,
95 const scalar pAmbient,
96 const scalar chi,
98) const
99{
100 scalar lBU =
101 B_
102 * sqrt
103 (
104 rho*sigma*d*cos(angle_*constant::mathematical::pi/360.0)
105 / sqr(rhoAv*Urel)
106 );
107
108 scalar pWalk = mag(pos - injectionPos);
109
110 if (pWalk > lBU)
111 {
112 liquidCore = 0.0;
113 }
114}
115
116
117// ************************************************************************* //
Templated atomization model class.
Primary Breakup Model for pressure swirl atomizers.
virtual bool calcChi() const
Flag to indicate if chi needs to be calculated.
virtual scalar initLiquidCore() const
Initial value of liquidCore.
virtual ~BlobsSheetAtomization()
Destructor.
Templated base class for dsmc cloud.
Definition: DSMCCloud.H:75
Random number generator.
Definition: Random.H:60
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
virtual bool update()
Update the mesh for both mesh motion and topology change.
const volScalarField & mu
Urel
Definition: pEqn.H:56
constexpr scalar pi(M_PI)
dimensionedScalar pos(const dimensionedScalar &ds)
dimensionedSymmTensor sqr(const dimensionedVector &dv)
dimensionedScalar sqrt(const dimensionedScalar &ds)
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
dimensionedScalar cos(const dimensionedScalar &ds)
dictionary dict
Random rndGen
Definition: createFields.H:23