syringePressureFvPatchScalarField.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) 2011-2016 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 
30 #include "volMesh.H"
32 #include "fvPatchFieldMapper.H"
33 #include "surfaceFields.H"
34 
35 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
36 
38 (
39  const fvPatch& p,
41 )
42 :
43  fixedValueFvPatchScalarField(p, iF),
44  phiName_("phi"),
45  curTimeIndex_(-1)
46 {}
47 
48 
50 (
51  const fvPatch& p,
53  const dictionary& dict
54 )
55 :
56  fixedValueFvPatchScalarField(p, iF, dict, false),
57  Ap_(dict.get<scalar>("Ap")),
58  Sp_(dict.get<scalar>("Sp")),
59  VsI_(dict.get<scalar>("VsI")),
60  tas_(dict.get<scalar>("tas")),
61  tae_(dict.get<scalar>("tae")),
62  tds_(dict.get<scalar>("tds")),
63  tde_(dict.get<scalar>("tde")),
64  psI_(dict.get<scalar>("psI")),
65  psi_(dict.get<scalar>("psi")),
66  ams_(dict.get<scalar>("ams")),
67  ams0_(ams_),
68  phiName_(dict.getOrDefault<word>("phi", "phi")),
69  curTimeIndex_(-1)
70 {
71  scalar ps = (psI_*VsI_ + ams_/psi_)/Vs(db().time().value());
73 }
74 
75 
77 (
79  const fvPatch& p,
81  const fvPatchFieldMapper& mapper
82 )
83 :
84  fixedValueFvPatchScalarField(sppsf, p, iF, mapper),
85  Ap_(sppsf.Ap_),
86  Sp_(sppsf.Sp_),
87  VsI_(sppsf.VsI_),
88  tas_(sppsf.tas_),
89  tae_(sppsf.tae_),
90  tds_(sppsf.tds_),
91  tde_(sppsf.tde_),
92  psI_(sppsf.psI_),
93  psi_(sppsf.psi_),
94  ams_(sppsf.ams_),
95  ams0_(sppsf.ams0_),
96  phiName_(sppsf.phiName_),
97  curTimeIndex_(-1)
98 {}
99 
100 
102 (
105 )
106 :
107  fixedValueFvPatchScalarField(sppsf, iF),
108  Ap_(sppsf.Ap_),
109  Sp_(sppsf.Sp_),
110  VsI_(sppsf.VsI_),
111  tas_(sppsf.tas_),
112  tae_(sppsf.tae_),
113  tds_(sppsf.tds_),
114  tde_(sppsf.tde_),
115  psI_(sppsf.psI_),
116  psi_(sppsf.psi_),
117  ams_(sppsf.ams_),
118  ams0_(sppsf.ams0_),
119  phiName_(sppsf.phiName_),
120  curTimeIndex_(-1)
121 {}
122 
123 
125 (
127 )
128 :
129  fixedValueFvPatchScalarField(sppsf),
130  Ap_(sppsf.Ap_),
131  Sp_(sppsf.Sp_),
132  VsI_(sppsf.VsI_),
133  tas_(sppsf.tas_),
134  tae_(sppsf.tae_),
135  tds_(sppsf.tds_),
136  tde_(sppsf.tde_),
137  psI_(sppsf.psI_),
138  psi_(sppsf.psi_),
139  ams_(sppsf.ams_),
140  ams0_(sppsf.ams0_),
141  phiName_(sppsf.phiName_),
142  curTimeIndex_(-1)
143 {}
144 
145 
146 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
147 
148 Foam::scalar Foam::syringePressureFvPatchScalarField::Vs(const scalar t) const
149 {
150  if (t < tas_)
151  {
152  return VsI_;
153  }
154  else if (t < tae_)
155  {
156  return
157  VsI_
158  + 0.5*Ap_*Sp_*sqr(t - tas_)/(tae_ - tas_);
159  }
160  else if (t < tds_)
161  {
162  return
163  VsI_
164  + 0.5*Ap_*Sp_*(tae_ - tas_)
165  + Ap_*Sp_*(t - tae_);
166  }
167  else if (t < tde_)
168  {
169  return
170  VsI_
171  + 0.5*Ap_*Sp_*(tae_ - tas_)
172  + Ap_*Sp_*(tds_ - tae_)
173  + Ap_*Sp_*(t - tds_)
174  - 0.5*Ap_*Sp_*sqr(t - tds_)/(tde_ - tds_);
175  }
176  else
177  {
178  return
179  VsI_
180  + 0.5*Ap_*Sp_*(tae_ - tas_)
181  + Ap_*Sp_*(tds_ - tae_)
182  + 0.5*Ap_*Sp_*(tde_ - tds_);
183  }
184 }
185 
186 
188 {
189  if (updated())
190  {
191  return;
192  }
193 
194  if (curTimeIndex_ != db().time().timeIndex())
195  {
196  ams0_ = ams_;
197  curTimeIndex_ = db().time().timeIndex();
198  }
199 
200  scalar t = db().time().value();
201  scalar deltaT = db().time().deltaTValue();
202 
203  const surfaceScalarField& phi =
204  db().lookupObject<surfaceScalarField>(phiName_);
205 
206  const fvsPatchField<scalar>& phip =
207  patch().patchField<surfaceScalarField, scalar>(phi);
208 
209  if (phi.dimensions() == dimVelocity*dimArea)
210  {
211  ams_ = ams0_ + deltaT*sum((*this*psi_)*phip);
212  }
213  else if (phi.dimensions() == dimDensity*dimVelocity*dimArea)
214  {
215  ams_ = ams0_ + deltaT*sum(phip);
216  }
217  else
218  {
220  << "dimensions of phi are not correct"
221  << "\n on patch " << this->patch().name()
222  << " of field " << this->internalField().name()
223  << " in file " << this->internalField().objectPath()
224  << exit(FatalError);
225  }
226 
227  scalar ps = (psI_*VsI_ + ams_/psi_)/Vs(t);
228 
229  operator==(ps);
230 
231  fixedValueFvPatchScalarField::updateCoeffs();
232 }
233 
234 
236 {
238 
239  os.writeEntry("Ap", Ap_);
240  os.writeEntry("Sp", Sp_);
241  os.writeEntry("VsI", VsI_);
242  os.writeEntry("tas", tas_);
243  os.writeEntry("tae", tae_);
244  os.writeEntry("tds", tds_);
245  os.writeEntry("tde", tde_);
246  os.writeEntry("psI", psI_);
247  os.writeEntry("psi", psi_);
248  os.writeEntry("ams", ams_);
249  os.writeEntryIfDifferent<word>("phi", "phi", phiName_);
250 
251  writeEntry("value", os);
252 }
253 
254 
255 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
256 
257 namespace Foam
258 {
260  (
263  );
264 }
265 
266 // ************************************************************************* //
Foam::fvPatchField
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: volSurfaceMapping.H:51
Foam::fvPatchField< scalar >::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::syringePressureFvPatchScalarField
This boundary condition provides a pressure condition, obtained from a zero-D model of the cylinder o...
Definition: syringePressureFvPatchScalarField.H:152
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:65
Foam::dimVelocity
const dimensionSet dimVelocity
Foam::dimDensity
const dimensionSet dimDensity
volMesh.H
Foam::fvsPatchField
An abstract base class with a fat-interface to all derived classes covering all possible ways in whic...
Definition: fvsPatchField.H:68
Foam::syringePressureFvPatchScalarField::write
virtual void write(Ostream &) const
Write.
Definition: syringePressureFvPatchScalarField.C:235
surfaceFields.H
Foam::surfaceFields.
fvPatchFieldMapper.H
Foam::dimArea
const dimensionSet dimArea(sqr(dimLength))
Definition: dimensionSets.H:59
Foam::operator==
tmp< faMatrix< Type > > operator==(const faMatrix< Type > &, const faMatrix< Type > &)
Foam::fvPatch
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:65
phi
surfaceScalarField & phi
Definition: setRegionFluidFields.H:8
Foam::syringePressureFvPatchScalarField::syringePressureFvPatchScalarField
syringePressureFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
Definition: syringePressureFvPatchScalarField.C:38
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::syringePressureFvPatchScalarField::updateCoeffs
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Definition: syringePressureFvPatchScalarField.C:187
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
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::sqr
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Definition: dimensionedSymmTensor.C:51
Foam::surfaceScalarField
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
Definition: surfaceFieldsFwd.H:54
Foam::foamVersion::patch
const std::string patch
OpenFOAM patch number as a std::string.
Foam::sum
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &df)
Definition: DimensionedFieldFunctions.C:327
Foam::Ostream::writeEntry
Ostream & writeEntry(const keyType &key, const T &value)
Write a keyword/value entry.
Definition: Ostream.H:236
syringePressureFvPatchScalarField.H
Foam::fvPatchFieldMapper
Foam::fvPatchFieldMapper.
Definition: fvPatchFieldMapper.H:47
timeIndex
label timeIndex
Definition: getTimeIndex.H:30
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::GeometricField< scalar, fvsPatchField, surfaceMesh >
Foam::DimensionedField
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: DimensionedField.H:54