windDrivenRainFoam.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) 2004-2011 OpenCFD Ltd.
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 Application
25  windDrivenRainFoam
26 
27 Description
28  Solves for wind-driven rain with an Eulerian multiphase model
29  Written by Aytac Kubilay, March 2012, ETH Zurich/Empa
30 
31  Latest Update: 17.02.2015
32 
33 \*---------------------------------------------------------------------------*/
34 
35 #include "fvCFD.H"
36 #include "simpleControl.H"
37 
38 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
39 
40 #include "updateValues.H"
41 
42 int main(int argc, char *argv[])
43 {
44  #include "setRootCase.H"
45  #include "createTime.H"
46  #include "createMesh.H"
47  #include "createFields.H"
48  #include "readGravitationalAcceleration.H"
49 
51 
52  // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
53 
54  #include "createRainFields.H"
55  #include "createTDFields.H"
56 
57  if (solveTD) { Info <<nl<< "Solving for the turbulent dispersion of raindrops" << endl; }
58  else { Info <<nl<< "Turbulent dispersion of raindrops is neglected" << endl; }
59 
60  GET_parameters(temp,rhoa,mua,rhop);
61 
62  Info <<nl<< "Temperature: " << temp.value() << " K" << endl;
63  Info << "Air density: " << rhoa.value() << " kg/m3" << endl;
64  Info << "Air dynamic viscosity: " << mua.value() << " kg/m-s" << endl;
65  Info << "Water density: " << rhop.value() << " kg/m3" << endl;
66 
67  while (simple.loop())
68  {
69  Info<< nl << "Time = " << runTime.timeName() << nl;
70 
71  for (int nonOrth=0; nonOrth<=simple.nNonOrthCorr(); nonOrth++)
72  {
73  for (int phase_no = 0; phase_no < phases.size(); phase_no++)
74  {
75  // phi is used by the inletOutlet boundary condition and courant number calculation
77  (
78  IOobject
79  (
80  "phi",
81  runTime.timeName(),
82  mesh,
83  IOobject::NO_READ,
84  IOobject::NO_WRITE
85  ),
86  phirain[phase_no]
87  );
89 
90  #include "CourantNo.H"
91 
92  #include "alphaEqns.H"
93 
94  dimensionedScalar dp ("dp", dimensionSet(0,1,0,0,0,0,0), phases[phase_no][0]);
95 
96  volScalarField magUr = mag(U - Urain[phase_no]);
97 
98  Re = (magUr*dp*rhoa)/mua;
99  CdRe = GET_CdRe(Re);
100  CdRe.correctBoundaryConditions();
101 
102  if (solveTD)
103  {
104  tfl = 0.2*(k/epsilon);
105  tp = (4*rhop*dp*dp)/(3*mua*CdRe);
106  Ctrain[phase_no] = sqrt( tfl/(tfl+tp) );
107  }
108 
109  nutrain = nut*sqr(Ctrain[phase_no]);
110 
111  #include "UEqns.H"
112 
113  }
114  }
115 
116  if (runTime.outputTime())
117  {
118  for (int phase_no = 0; phase_no < phases.size(); phase_no++)
119  {
120  Urain[phase_no].write();
121  alpharain[phase_no].write();
122  //Ctrain[phase_no].write();
123  }
124  }
125  runTime.write();
126 
127  Info<< "ExecutionTime = " << runTime.elapsedCpuTime() << " s"
128  << " ClockTime = " << runTime.elapsedClockTime() << " s"
129  << nl << endl;
130  }
131 
132  Info<< "Writing final output\n" << endl;
133  for (int phase_no = 0; phase_no < phases.size(); phase_no++)
134  {
135  Urain[phase_no].write();
136  alpharain[phase_no].write();
137  //Ctrain[phase_no].write();
138  }
139 
140  #include "calculateCatchRatio.H"
141 
142  Info<< "ExecutionTime = " << runTime.elapsedCpuTime() << " s"
143  << " ClockTime = " << runTime.elapsedClockTime() << " s"
144  << nl << endl;
145 
146  Info<< "End\n" << endl;
147 
148  return 0;
149 }
150 
151 // ************************************************************************* //
createRainFields.H
Foam::IOobject
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:91
updateValues.H
Foam::simpleControl
SIMPLE control class to supply convergence information/checks for the SIMPLE loop.
Definition: simpleControl.H:46
phi
surfaceScalarField & phi
Definition: setRegionFluidFields.H:8
Ctrain
PtrList< volScalarField > Ctrain
Definition: createTDFields.H:1
Foam::Re
scalarField Re(const UList< complex > &cf)
Definition: complexFields.C:97
phirain
PtrList< surfaceScalarField > phirain
Definition: createRainFields.H:3
tp
volScalarField tp(IOobject("tp", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::NO_WRITE), mesh, dimensionedScalar("tp", dimensionSet(0, 0, 1, 0, 0, 0, 0), 0))
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
Foam::dimensionSet
Dimension set for the base types.
Definition: dimensionSet.H:116
Foam::mag
dimensioned< scalar > mag(const dimensioned< Type > &)
UEqns.H
simple
Simple relative velocity model.
U
U
Definition: pEqn.H:46
nut
nut
Definition: createTDFields.H:71
Foam::nl
static const char nl
Definition: Ostream.H:260
Foam::Info
messageStream Info
simpleControl.H
GET_parameters
void GET_parameters(dimensionedScalar temp_, dimensionedScalar &rhoa_, dimensionedScalar &mua_, dimensionedScalar &rhop_)
Definition: updateValues.H:1
GET_CdRe
volScalarField GET_CdRe(volScalarField Reynolds)
Definition: updateValues.H:56
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:18
Foam::dimensioned
Generic dimensioned Type class.
Definition: dimensionedScalarFwd.H:41
calculateCatchRatio.H
alpharain
PtrList< volScalarField > alpharain
Definition: createRainFields.H:5
main
int main(int argc, char *argv[])
Definition: windDrivenRainFoam.C:39
createFields.H
setRootCase.H
epsilon
epsilon
Definition: createTDFields.H:56
Foam::sqr
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Definition: dimensionedSymmTensor.C:49
simple
const dictionary & simple
Definition: readFluidMultiRegionSIMPLEControls.H:1
Urain
PtrList< volVectorField > Urain
Definition: createRainFields.H:1
tfl
volScalarField tfl(IOobject("tfl", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::NO_WRITE), mesh, dimensionedScalar("tfl", dimensionSet(0, 0, 1, 0, 0, 0, 0), 0))
alphaEqns.H
Foam::sqrt
dimensionedScalar sqrt(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:142
nutrain
volScalarField nutrain(IOobject("nutrain", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::NO_WRITE), mesh, dimensionedScalar("nutrain", dimensionSet(0, 2,-1, 0, 0, 0, 0), 1))
phases
multiphaseSystem::phaseModelList & phases
Definition: createFields.H:11
createMesh.H
k
label k
Boltzmann constant.
Definition: LISASMDCalcMethod2.H:41
createTDFields.H
createTime.H
fvCFD.H
Foam::GeometricField
Generic GeometricField class.
Definition: surfaceFieldsFwd.H:52