semiPermeableBaffleVelocityFvPatchVectorField.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2017-2018 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 
29 #include "fvPatchFieldMapper.H"
30 #include "volFields.H"
31 #include "surfaceFields.H"
32 #include "psiReactionThermo.H"
33 #include "rhoReactionThermo.H"
34 
35 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
36 
38 Foam::semiPermeableBaffleVelocityFvPatchVectorField::composition() const
39 {
40  const word& name = basicThermo::dictName;
41 
42  if (db().foundObject<psiReactionThermo>(name))
43  {
44  return db().lookupObject<psiReactionThermo>(name).composition();
45  }
46  else if (db().foundObject<rhoReactionThermo>(name))
47  {
48  return db().lookupObject<rhoReactionThermo>(name).composition();
49  }
50  else
51  {
53  << "Could not find a multi-component thermodynamic model."
54  << exit(FatalError);
55 
56  return NullObjectRef<basicSpecieMixture>();
57  }
58 }
59 
60 
61 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
62 
65 (
66  const fvPatch& p,
68 )
69 :
70  fixedValueFvPatchVectorField(p, iF),
71  rhoName_("rho")
72 {}
73 
74 
77 (
78  const fvPatch& p,
80  const dictionary& dict
81 )
82 :
83  fixedValueFvPatchVectorField(p, iF),
84  rhoName_(dict.lookupOrDefault<word>("rho", "rho"))
85 {
86  fvPatchVectorField::operator==(vectorField("value", dict, p.size()));
87 }
88 
89 
92 (
94  const fvPatch& p,
96  const fvPatchFieldMapper& mapper
97 )
98 :
99  fixedValueFvPatchVectorField(ptf, p, iF, mapper),
100  rhoName_(ptf.rhoName_)
101 {}
102 
103 
106 (
108 )
109 :
110  fixedValueFvPatchVectorField(ptf),
111  rhoName_(ptf.rhoName_)
112 {}
113 
114 
117 (
120 )
121 :
122  fixedValueFvPatchVectorField(ptf, iF),
123  rhoName_(ptf.rhoName_)
124 {}
125 
126 
127 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
128 
130 {
131  if (updated())
132  {
133  return;
134  }
135 
137 
138  const scalarField& rhop =
139  patch().lookupPatchField<volScalarField, scalar>(rhoName_);
140 
141  const PtrList<volScalarField>& Y = composition().Y();
142 
143  scalarField phip(patch().size(), Zero);
144  forAll(Y, i)
145  {
146  const fvPatchScalarField& Yp = Y[i].boundaryField()[patch().index()];
147 
148  if (!isA<YBCType>(Yp))
149  {
151  << "The mass-fraction condition on patch " << patch().name()
152  << " is not of type " << YBCType::typeName << "."
153  << exit(FatalError);
154  }
155 
156  phip += refCast<const YBCType>(Yp).phiY();
157  }
158 
159  this->operator==(patch().nf()*phip/(rhop*patch().magSf()));
160 
161  fixedValueFvPatchVectorField::updateCoeffs();
162 }
163 
164 
166 (
167  Ostream& os
168 ) const
169 {
171  writeEntryIfDifferent<word>(os, "rho", "rho", rhoName_);
172  writeEntry("value", os);
173 }
174 
175 
176 // * * * * * * * * * * * * * * Build Macro Function * * * * * * * * * * * * //
177 
178 namespace Foam
179 {
181  (
184  );
185 }
186 
187 // ************************************************************************* //
Foam::surfaceFields.
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:428
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
error FatalError
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:319
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:61
This is a mass-fraction boundary condition for a semi-permeable baffle.
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:66
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:362
Specialization of basicMultiComponentMixture for a mixture consisting of a number for molecular speci...
PtrList< volScalarField > & Y()
Return the mass-fraction fields.
Macros for easy insertion into run-time selection tables.
virtual void operator==(const fvPatchField< Type > &)
Definition: fvPatchField.C:561
const word dictName() const
Return the local dictionary name (final part of scoped name)
Definition: dictionary.H:115
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
tmp< fvMatrix< Type > > operator==(const fvMatrix< Type > &, const fvMatrix< Type > &)
A class for handling words, derived from string.
Definition: word.H:59
Foam::fvPatchFieldMapper.
static const zero Zero
Definition: zero.H:97
virtual label size() const
Return size.
Definition: fvPatch.H:161
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:53
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
semiPermeableBaffleVelocityFvPatchVectorField(const fvPatch &, const DimensionedField< vector, volMesh > &)
Construct from patch and internal field.
This is a velocity boundary condition for a semi-permeable baffle.
T lookupOrDefault(const word &, const T &, bool recursive=false, bool patternMatch=true) const
Find and return a T,.
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: List.H:63
PtrList< volScalarField > & Y
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Field< vector > vectorField
Specialisation of Field<T> for vector.
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
Namespace for OpenFOAM.