semiPermeableBaffleVelocityFvPatchVectorField.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) 2017 OpenFOAM Foundation
9  Copyright (C) 2020 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
32 #include "fvPatchFieldMapper.H"
33 #include "volFields.H"
34 #include "surfaceFields.H"
35 #include "psiReactionThermo.H"
36 #include "rhoReactionThermo.H"
37 
38 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
39 
41 Foam::semiPermeableBaffleVelocityFvPatchVectorField::composition() const
42 {
43  const word& name = basicThermo::dictName;
44 
45  if (db().foundObject<psiReactionThermo>(name))
46  {
47  return db().lookupObject<psiReactionThermo>(name).composition();
48  }
49  else if (db().foundObject<rhoReactionThermo>(name))
50  {
51  return db().lookupObject<rhoReactionThermo>(name).composition();
52  }
53  else
54  {
56  << "Could not find a multi-component thermodynamic model."
57  << exit(FatalError);
58 
59  return NullObjectRef<basicSpecieMixture>();
60  }
61 }
62 
63 
64 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
65 
68 (
69  const fvPatch& p,
71 )
72 :
73  fixedValueFvPatchVectorField(p, iF),
74  rhoName_("rho")
75 {}
76 
77 
80 (
81  const fvPatch& p,
83  const dictionary& dict
84 )
85 :
86  fixedValueFvPatchVectorField(p, iF),
87  rhoName_(dict.getOrDefault<word>("rho", "rho"))
88 {
90 }
91 
92 
95 (
97  const fvPatch& p,
99  const fvPatchFieldMapper& mapper
100 )
101 :
102  fixedValueFvPatchVectorField(ptf, p, iF, mapper),
103  rhoName_(ptf.rhoName_)
104 {}
105 
106 
109 (
111 )
112 :
113  fixedValueFvPatchVectorField(ptf),
114  rhoName_(ptf.rhoName_)
115 {}
116 
117 
120 (
123 )
124 :
125  fixedValueFvPatchVectorField(ptf, iF),
126  rhoName_(ptf.rhoName_)
127 {}
128 
129 
130 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
131 
133 {
134  if (updated())
135  {
136  return;
137  }
138 
140 
141  const scalarField& rhop =
142  patch().lookupPatchField<volScalarField, scalar>(rhoName_);
143 
144  const PtrList<volScalarField>& Y = composition().Y();
145 
146  scalarField phip(patch().size(), Zero);
147  forAll(Y, i)
148  {
149  const fvPatchScalarField& Yp = Y[i].boundaryField()[patch().index()];
150 
151  if (!isA<YBCType>(Yp))
152  {
154  << "The mass-fraction condition on patch " << patch().name()
155  << " is not of type " << YBCType::typeName << "."
156  << exit(FatalError);
157  }
158 
159  phip += refCast<const YBCType>(Yp).phiY();
160  }
161 
162  this->operator==(patch().nf()*phip/(rhop*patch().magSf()));
163 
164  fixedValueFvPatchVectorField::updateCoeffs();
165 }
166 
167 
169 (
170  Ostream& os
171 ) const
172 {
174  os.writeEntryIfDifferent<word>("rho", "rho", rhoName_);
175  writeEntry("value", os);
176 }
177 
178 
179 // * * * * * * * * * * * * * * Build Macro Function * * * * * * * * * * * * //
180 
181 namespace Foam
182 {
184  (
187  );
188 }
189 
190 // ************************************************************************* //
Foam::fvPatchField< scalar >
volFields.H
Foam::fvPatchField< vector >::write
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:384
Foam::Ostream::writeEntryIfDifferent
Ostream & writeEntryIfDifferent(const word &key, const T &value1, const T &value2)
Write a keyword/value entry only when the two values differ.
Definition: Ostream.H:248
p
volScalarField & p
Definition: createFieldRefs.H:8
Foam::semiPermeableBaffleMassFractionFvPatchScalarField
This is a mass-fraction boundary condition for a semi-permeable baffle.
Definition: semiPermeableBaffleMassFractionFvPatchScalarField.H:127
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:65
Foam::basicSpecieMixture
Specialization of basicMultiComponentMixture for a mixture consisting of a number for molecular speci...
Definition: basicSpecieMixture.H:58
Foam::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
Foam::fvPatchField< vector >::operator==
virtual void operator==(const fvPatchField< vector > &)
Definition: fvPatchField.C:570
surfaceFields.H
Foam::surfaceFields.
fvPatchFieldMapper.H
Foam::basicMultiComponentMixture::Y
PtrList< volScalarField > & Y()
Return the mass-fraction fields.
Definition: basicMultiComponentMixtureI.H:69
composition
basicSpecieMixture & composition
Definition: createFieldRefs.H:6
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
Foam::vectorField
Field< vector > vectorField
Specialisation of Field<T> for vector.
Definition: primitiveFieldsFwd.H:54
Foam::semiPermeableBaffleVelocityFvPatchVectorField
This is a velocity boundary condition for a semi-permeable baffle.
Definition: semiPermeableBaffleVelocityFvPatchVectorField.H:83
Foam::semiPermeableBaffleVelocityFvPatchVectorField::semiPermeableBaffleVelocityFvPatchVectorField
semiPermeableBaffleVelocityFvPatchVectorField(const fvPatch &, const DimensionedField< vector, volMesh > &)
Construct from patch and internal field.
Definition: semiPermeableBaffleVelocityFvPatchVectorField.C:68
rhoReactionThermo.H
Foam::Field< scalar >
Foam::operator==
tmp< faMatrix< Type > > operator==(const faMatrix< Type > &, const faMatrix< Type > &)
semiPermeableBaffleVelocityFvPatchVectorField.H
semiPermeableBaffleMassFractionFvPatchScalarField.H
Foam::fvPatch
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:65
Foam::PtrList
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition: List.H:59
Foam::volScalarField
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:57
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::FatalError
error FatalError
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:123
os
OBJstream os(runTime.globalPath()/outputName)
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Y
PtrList< volScalarField > & Y
Definition: createFieldRefs.H:7
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
Foam::semiPermeableBaffleVelocityFvPatchVectorField::write
virtual void write(Ostream &) const
Write.
Definition: semiPermeableBaffleVelocityFvPatchVectorField.C:169
psiReactionThermo.H
Foam::foamVersion::patch
const std::string patch
OpenFOAM patch number as a std::string.
Foam::semiPermeableBaffleVelocityFvPatchVectorField::updateCoeffs
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Definition: semiPermeableBaffleVelocityFvPatchVectorField.C:132
Foam::name
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
Foam::fvPatchFieldMapper
Foam::fvPatchFieldMapper.
Definition: fvPatchFieldMapper.H:47
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:56
Foam::makePatchTypeField
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
Foam::dictionary::getOrDefault
T getOrDefault(const word &keyword, const T &deflt, enum keyType::option matchOpt=keyType::REGEX) const
Definition: dictionaryTemplates.C:148
Foam::DimensionedField
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: DimensionedField.H:54
Foam::dictionary::dictName
word dictName() const
The local dictionary name (final part of scoped name)
Definition: dictionaryI.H:60