adjustPhi.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) 2011-2015 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 
26 #include "adjustPhi.H"
27 #include "volFields.H"
28 #include "surfaceFields.H"
30 
31 // * * * * * * * * * * * * * * * Global Functions * * * * * * * * * * * * * //
32 
33 bool Foam::adjustPhi
34 (
36  const volVectorField& U,
38 )
39 {
40  if (p.needReference())
41  {
42  scalar massIn = 0.0;
43  scalar fixedMassOut = 0.0;
44  scalar adjustableMassOut = 0.0;
45 
46  surfaceScalarField::GeometricBoundaryField& bphi = phi.boundaryField();
47 
48  forAll(bphi, patchi)
49  {
50  const fvPatchVectorField& Up = U.boundaryField()[patchi];
51  const fvsPatchScalarField& phip = bphi[patchi];
52 
53  if (!phip.coupled())
54  {
55  if (Up.fixesValue() && !isA<inletOutletFvPatchVectorField>(Up))
56  {
57  forAll(phip, i)
58  {
59  if (phip[i] < 0.0)
60  {
61  massIn -= phip[i];
62  }
63  else
64  {
65  fixedMassOut += phip[i];
66  }
67  }
68  }
69  else
70  {
71  forAll(phip, i)
72  {
73  if (phip[i] < 0.0)
74  {
75  massIn -= phip[i];
76  }
77  else
78  {
79  adjustableMassOut += phip[i];
80  }
81  }
82  }
83  }
84  }
85 
86  // Calculate the total flux in the domain, used for normalisation
87  scalar totalFlux = VSMALL + sum(mag(phi)).value();
88 
89  reduce(massIn, sumOp<scalar>());
90  reduce(fixedMassOut, sumOp<scalar>());
91  reduce(adjustableMassOut, sumOp<scalar>());
92 
93  scalar massCorr = 1.0;
94  scalar magAdjustableMassOut = mag(adjustableMassOut);
95 
96  if
97  (
98  magAdjustableMassOut > VSMALL
99  && magAdjustableMassOut/totalFlux > SMALL
100  )
101  {
102  massCorr = (massIn - fixedMassOut)/adjustableMassOut;
103  }
104  else if (mag(fixedMassOut - massIn)/totalFlux > 1e-6)
105  {
107  << "Continuity error cannot be removed by adjusting the"
108  " outflow.\nPlease check the velocity boundary conditions"
109  " and/or run potentialFoam to initialise the outflow." << nl
110  << "Total flux : " << totalFlux << nl
111  << "Specified mass inflow : " << massIn << nl
112  << "Specified mass outflow : " << fixedMassOut << nl
113  << "Adjustable mass outflow : " << adjustableMassOut << nl
114  << exit(FatalError);
115  }
116 
117  forAll(bphi, patchi)
118  {
119  const fvPatchVectorField& Up = U.boundaryField()[patchi];
120  fvsPatchScalarField& phip = bphi[patchi];
121 
122  if (!phip.coupled())
123  {
124  if
125  (
126  !Up.fixesValue()
127  || isA<inletOutletFvPatchVectorField>(Up)
128  )
129  {
130  forAll(phip, i)
131  {
132  if (phip[i] > 0.0)
133  {
134  phip[i] *= massCorr;
135  }
136  }
137  }
138  }
139  }
140 
141  return mag(massIn)/totalFlux < SMALL
142  && mag(fixedMassOut)/totalFlux < SMALL
143  && mag(adjustableMassOut)/totalFlux < SMALL;
144  }
145  else
146  {
147  return false;
148  }
149 }
150 
151 
152 // ************************************************************************* //
Foam::fvPatchField
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:65
volFields.H
Foam::adjustPhi
bool adjustPhi(surfaceScalarField &phi, const volVectorField &U, volScalarField &p)
Adjust the balance of fluxes to obey continuity.
Definition: adjustPhi.C:34
p
p
Definition: pEqn.H:62
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:406
phi
surfaceScalarField & phi
Definition: setRegionFluidFields.H:8
Foam::fvsPatchField
An abstract base class with a fat-interface to all derived classes covering all possible ways in whic...
Definition: fvsPatchField.H:65
surfaceFields.H
Foam::surfaceFields.
Foam::dimensioned::value
const Type & value() const
Return const reference to value.
Definition: dimensionedType.C:261
Foam::mag
dimensioned< scalar > mag(const dimensioned< Type > &)
Foam::fvPatchField::fixesValue
virtual bool fixesValue() const
Return true if this patch field fixes a value.
Definition: fvPatchField.H:336
adjustPhi.H
For cases which do no have a pressure boundary adjust the balance of fluxes to obey continuity....
U
U
Definition: pEqn.H:46
Foam::reduce
void reduce(const List< UPstream::commsStruct > &comms, T &Value, const BinaryOp &bop, const int tag, const label comm)
Definition: PstreamReduceOps.H:43
Foam::nl
static const char nl
Definition: Ostream.H:260
Foam::FatalError
error FatalError
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:318
Foam::sumOp
Definition: ops.H:162
patchi
label patchi
Definition: getPatchFieldScalar.H:1
Foam::constant::electromagnetic::e
const dimensionedScalar e
Elementary charge.
Definition: doubleFloat.H:94
Foam::sum
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &df)
Definition: DimensionedFieldFunctions.C:333
Foam::GeometricField
Generic GeometricField class.
Definition: surfaceFieldsFwd.H:52
inletOutletFvPatchFields.H
Foam::fvsPatchField::coupled
virtual bool coupled() const
Return true if this patch field is coupled.
Definition: fvsPatchField.H:308
Foam::GeometricField::GeometricBoundaryField
Definition: GeometricField.H:105