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 (
35  surfaceScalarField& phi,
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-8)
105  {
107  (
108  "adjustPhi"
109  "("
110  "surfaceScalarField&, "
111  "const volVectorField&,"
112  "volScalarField&"
113  ")"
114  ) << "Continuity error cannot be removed by adjusting the"
115  " outflow.\nPlease check the velocity boundary conditions"
116  " and/or run potentialFoam to initialise the outflow." << nl
117  << "Total flux : " << totalFlux << nl
118  << "Specified mass inflow : " << massIn << nl
119  << "Specified mass outflow : " << fixedMassOut << nl
120  << "Adjustable mass outflow : " << adjustableMassOut << nl
121  << exit(FatalError);
122  }
123 
124  forAll(bphi, patchi)
125  {
126  const fvPatchVectorField& Up = U.boundaryField()[patchi];
127  fvsPatchScalarField& phip = bphi[patchi];
128 
129  if (!phip.coupled())
130  {
131  if
132  (
133  !Up.fixesValue()
134  || isA<inletOutletFvPatchVectorField>(Up)
135  )
136  {
137  forAll(phip, i)
138  {
139  if (phip[i] > 0.0)
140  {
141  phip[i] *= massCorr;
142  }
143  }
144  }
145  }
146  }
147 
148  return mag(massIn)/totalFlux < SMALL
149  && mag(fixedMassOut)/totalFlux < SMALL
150  && mag(adjustableMassOut)/totalFlux < SMALL;
151  }
152  else
153  {
154  return false;
155  }
156 }
157 
158 
159 // ************************************************************************* //
void reduce(const List< UPstream::commsStruct > &comms, T &Value, const BinaryOp &bop, const int tag, const label comm)
virtual bool coupled() const
Return true if this patch field is coupled.
bool needReference() const
Does the field need a reference level for solution.
Foam::surfaceFields.
dimensioned< scalar > mag(const dimensioned< Type > &)
GeometricBoundaryField & boundaryField()
Return reference to GeometricBoundaryField.
bool adjustPhi(surfaceScalarField &phi, const volVectorField &U, volScalarField &p)
Adjust the balance of fluxes to obey continuity.
Definition: adjustPhi.C:34
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
static const char nl
Definition: Ostream.H:260
const dimensionedScalar e
Elementary charge.
Definition: doubleFloat.H:78
#define forAll(list, i)
Definition: UList.H:421
label patchi
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &df)
virtual bool fixesValue() const
Return true if this patch field fixes a value.
Definition: fvPatchField.H:336
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:65
#define FatalErrorIn(functionName)
Report an error message using Foam::FatalError.
Definition: error.H:314
error FatalError
For cases which do no have a pressure boundary adjust the balance of fluxes to obey continuity...
An abstract base class with a fat-interface to all derived classes covering all possible ways in whic...
Definition: fvsPatchField.H:65