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 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 \*---------------------------------------------------------------------------*/
27 
31 #include "fvPatchFieldMapper.H"
32 #include "volFields.H"
33 #include "surfaceFields.H"
34 #include "psiReactionThermo.H"
35 #include "rhoReactionThermo.H"
36 
37 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
38 
40 Foam::semiPermeableBaffleVelocityFvPatchVectorField::composition() const
41 {
42  const word& name = basicThermo::dictName;
43 
44  if (db().foundObject<psiReactionThermo>(name))
45  {
46  return db().lookupObject<psiReactionThermo>(name).composition();
47  }
48  else if (db().foundObject<rhoReactionThermo>(name))
49  {
50  return db().lookupObject<rhoReactionThermo>(name).composition();
51  }
52  else
53  {
55  << "Could not find a multi-component thermodynamic model."
56  << exit(FatalError);
57 
58  return NullObjectRef<basicSpecieMixture>();
59  }
60 }
61 
62 
63 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
64 
67 (
68  const fvPatch& p,
70 )
71 :
72  fixedValueFvPatchVectorField(p, iF),
73  rhoName_("rho")
74 {}
75 
76 
79 (
80  const fvPatch& p,
82  const dictionary& dict
83 )
84 :
85  fixedValueFvPatchVectorField(p, iF),
86  rhoName_(dict.lookupOrDefault<word>("rho", "rho"))
87 {
89 }
90 
91 
94 (
96  const fvPatch& p,
98  const fvPatchFieldMapper& mapper
99 )
100 :
101  fixedValueFvPatchVectorField(ptf, p, iF, mapper),
102  rhoName_(ptf.rhoName_)
103 {}
104 
105 
108 (
110 )
111 :
112  fixedValueFvPatchVectorField(ptf),
113  rhoName_(ptf.rhoName_)
114 {}
115 
116 
119 (
122 )
123 :
124  fixedValueFvPatchVectorField(ptf, iF),
125  rhoName_(ptf.rhoName_)
126 {}
127 
128 
129 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
130 
132 {
133  if (updated())
134  {
135  return;
136  }
137 
139 
140  const scalarField& rhop =
141  patch().lookupPatchField<volScalarField, scalar>(rhoName_);
142 
143  const PtrList<volScalarField>& Y = composition().Y();
144 
145  scalarField phip(patch().size(), Zero);
146  forAll(Y, i)
147  {
148  const fvPatchScalarField& Yp = Y[i].boundaryField()[patch().index()];
149 
150  if (!isA<YBCType>(Yp))
151  {
153  << "The mass-fraction condition on patch " << patch().name()
154  << " is not of type " << YBCType::typeName << "."
155  << exit(FatalError);
156  }
157 
158  phip += refCast<const YBCType>(Yp).phiY();
159  }
160 
161  this->operator==(patch().nf()*phip/(rhop*patch().magSf()));
162 
163  fixedValueFvPatchVectorField::updateCoeffs();
164 }
165 
166 
168 (
169  Ostream& os
170 ) const
171 {
173  os.writeEntryIfDifferent<word>("rho", "rho", rhoName_);
174  writeEntry("value", os);
175 }
176 
177 
178 // * * * * * * * * * * * * * * Build Macro Function * * * * * * * * * * * * //
179 
180 namespace Foam
181 {
183  (
186  );
187 }
188 
189 // ************************************************************************* //
Foam::fvPatchField< scalar >
volFields.H
Foam::fvPatchField< vector >::write
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:364
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:231
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:62
Foam::basicSpecieMixture
Specialization of basicMultiComponentMixture for a mixture consisting of a number for molecular speci...
Definition: basicSpecieMixture.H:54
Foam::Zero
static constexpr const zero Zero
Global zero.
Definition: zero.H:128
Foam::fvPatchField< vector >::operator==
virtual void operator==(const fvPatchField< vector > &)
Definition: fvPatchField.C:545
surfaceFields.H
Foam::surfaceFields.
fvPatchFieldMapper.H
Foam::basicMultiComponentMixture::Y
PtrList< volScalarField > & Y()
Return the mass-fraction fields.
Definition: basicMultiComponentMixtureI.H:69
Foam::dictionary::lookupOrDefault
T lookupOrDefault(const word &keyword, const T &deflt, enum keyType::option matchOpt=keyType::REGEX) const
Definition: dictionary.H:1241
composition
basicSpecieMixture & composition
Definition: createFieldRefs.H:6
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:290
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:67
rhoReactionThermo.H
Foam::Field< scalar >
Foam::operator==
tmp< faMatrix< Type > > operator==(const faMatrix< Type > &, const faMatrix< Type > &)
semiPermeableBaffleVelocityFvPatchVectorField.H
Foam::name
word name(const complex &c)
Return string representation of complex.
Definition: complex.C:76
semiPermeableBaffleMassFractionFvPatchScalarField.H
Foam::fvPatch
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:63
Foam::PtrList
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition: List.H:65
Foam::volScalarField
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:57
Foam::dictionary::dictName
word dictName() const
The local dictionary name (final part of scoped name)
Definition: dictionary.H:458
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:121
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:355
Foam::semiPermeableBaffleVelocityFvPatchVectorField::write
virtual void write(Ostream &) const
Write.
Definition: semiPermeableBaffleVelocityFvPatchVectorField.C:168
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:131
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::DimensionedField
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: DimensionedField.H:54