syringePressureFvPatchScalarField.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) 2011-2021 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 
27 #include "volMesh.H"
29 #include "fvPatchFieldMapper.H"
30 #include "surfaceFields.H"
31 
32 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
33 
35 (
36  const fvPatch& p,
38 )
39 :
40  fixedValueFvPatchScalarField(p, iF),
41  phiName_("phi"),
42  curTimeIndex_(-1)
43 {}
44 
45 
47 (
48  const fvPatch& p,
50  const dictionary& dict
51 )
52 :
53  fixedValueFvPatchScalarField(p, iF, dict, false),
54  Ap_(dict.lookup<scalar>("Ap")),
55  Sp_(dict.lookup<scalar>("Sp")),
56  VsI_(dict.lookup<scalar>("VsI")),
57  tas_(dict.lookup<scalar>("tas")),
58  tae_(dict.lookup<scalar>("tae")),
59  tds_(dict.lookup<scalar>("tds")),
60  tde_(dict.lookup<scalar>("tde")),
61  psI_(dict.lookup<scalar>("psI")),
62  psi_(dict.lookup<scalar>("psi")),
63  ams_(dict.lookup<scalar>("ams")),
64  ams0_(ams_),
65  phiName_(dict.lookupOrDefault<word>("phi", "phi")),
66  curTimeIndex_(-1)
67 {
68  scalar ps = (psI_*VsI_ + ams_/psi_)/Vs(db().time().value());
70 }
71 
72 
74 (
76  const fvPatch& p,
78  const fvPatchFieldMapper& mapper
79 )
80 :
81  fixedValueFvPatchScalarField(sppsf, p, iF, mapper),
82  Ap_(sppsf.Ap_),
83  Sp_(sppsf.Sp_),
84  VsI_(sppsf.VsI_),
85  tas_(sppsf.tas_),
86  tae_(sppsf.tae_),
87  tds_(sppsf.tds_),
88  tde_(sppsf.tde_),
89  psI_(sppsf.psI_),
90  psi_(sppsf.psi_),
91  ams_(sppsf.ams_),
92  ams0_(sppsf.ams0_),
93  phiName_(sppsf.phiName_),
94  curTimeIndex_(-1)
95 {}
96 
97 
99 (
102 )
103 :
104  fixedValueFvPatchScalarField(sppsf, iF),
105  Ap_(sppsf.Ap_),
106  Sp_(sppsf.Sp_),
107  VsI_(sppsf.VsI_),
108  tas_(sppsf.tas_),
109  tae_(sppsf.tae_),
110  tds_(sppsf.tds_),
111  tde_(sppsf.tde_),
112  psI_(sppsf.psI_),
113  psi_(sppsf.psi_),
114  ams_(sppsf.ams_),
115  ams0_(sppsf.ams0_),
116  phiName_(sppsf.phiName_),
117  curTimeIndex_(-1)
118 {}
119 
120 
121 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
122 
123 Foam::scalar Foam::syringePressureFvPatchScalarField::Vs(const scalar t) const
124 {
125  if (t < tas_)
126  {
127  return VsI_;
128  }
129  else if (t < tae_)
130  {
131  return
132  VsI_
133  + 0.5*Ap_*Sp_*sqr(t - tas_)/(tae_ - tas_);
134  }
135  else if (t < tds_)
136  {
137  return
138  VsI_
139  + 0.5*Ap_*Sp_*(tae_ - tas_)
140  + Ap_*Sp_*(t - tae_);
141  }
142  else if (t < tde_)
143  {
144  return
145  VsI_
146  + 0.5*Ap_*Sp_*(tae_ - tas_)
147  + Ap_*Sp_*(tds_ - tae_)
148  + Ap_*Sp_*(t - tds_)
149  - 0.5*Ap_*Sp_*sqr(t - tds_)/(tde_ - tds_);
150  }
151  else
152  {
153  return
154  VsI_
155  + 0.5*Ap_*Sp_*(tae_ - tas_)
156  + Ap_*Sp_*(tds_ - tae_)
157  + 0.5*Ap_*Sp_*(tde_ - tds_);
158  }
159 }
160 
161 
163 {
164  if (updated())
165  {
166  return;
167  }
168 
169  if (curTimeIndex_ != db().time().timeIndex())
170  {
171  ams0_ = ams_;
172  curTimeIndex_ = db().time().timeIndex();
173  }
174 
175  scalar t = db().time().value();
176  scalar deltaT = db().time().deltaTValue();
177 
178  const surfaceScalarField& phi =
179  db().lookupObject<surfaceScalarField>(phiName_);
180 
181  const fvsPatchField<scalar>& phip =
182  patch().patchField<surfaceScalarField, scalar>(phi);
183 
184  if (phi.dimensions() == dimFlux)
185  {
186  ams_ = ams0_ + deltaT*sum((*this*psi_)*phip);
187  }
188  else if (phi.dimensions() == dimMassFlux)
189  {
190  ams_ = ams0_ + deltaT*sum(phip);
191  }
192  else
193  {
195  << "dimensions of phi are not correct"
196  << "\n on patch " << this->patch().name()
197  << " of field " << this->internalField().name()
198  << " in file " << this->internalField().objectPath()
199  << exit(FatalError);
200  }
201 
202  scalar ps = (psI_*VsI_ + ams_/psi_)/Vs(t);
203 
204  operator==(ps);
205 
206  fixedValueFvPatchScalarField::updateCoeffs();
207 }
208 
209 
211 {
213 
214  writeEntry(os, "Ap", Ap_);
215  writeEntry(os, "Sp", Sp_);
216  writeEntry(os, "VsI", VsI_);
217  writeEntry(os, "tas", tas_);
218  writeEntry(os, "tae", tae_);
219  writeEntry(os, "tds", tds_);
220  writeEntry(os, "tde", tde_);
221  writeEntry(os, "psI", psI_);
222  writeEntry(os, "psi", psi_);
223  writeEntry(os, "ams", ams_);
224  writeEntryIfDifferent<word>(os, "phi", "phi", phiName_);
225 
226  writeEntry(os, "value", *this);
227 }
228 
229 
230 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
231 
232 namespace Foam
233 {
235  (
238  );
239 }
240 
241 // ************************************************************************* //
Foam::surfaceFields.
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:156
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:306
dimensionedSymmTensor sqr(const dimensionedVector &dv)
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:63
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:243
Macros for easy insertion into run-time selection tables.
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &df)
const dimensionSet & dimensions() const
Return dimensions.
tmp< fvMatrix< Type > > operator==(const fvMatrix< Type > &, const fvMatrix< Type > &)
A class for handling words, derived from string.
Definition: word.H:59
Foam::fvPatchFieldMapper.
const dimensionSet dimFlux
phi
Definition: correctPhi.H:3
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:54
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
void writeEntry(Ostream &os, const HashTable< T, Key, Hash > &ht)
Definition: HashTableIO.C:96
syringePressureFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
const dimensionSet dimMassFlux
T lookupOrDefault(const word &, const T &, bool recursive=false, bool patternMatch=true) const
Find and return a T,.
label timeIndex
Definition: getTimeIndex.H:4
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
This boundary condition provides a pressure condition, obtained from a zero-D model of the cylinder o...
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
An abstract base class with a fat-interface to all derived classes covering all possible ways in whic...
Definition: fvsPatchField.H:65
Namespace for OpenFOAM.
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:864