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 -------------------------------------------------------------------------------
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 
29 #include "volMesh.H"
31 #include "fvPatchFieldMapper.H"
32 #include "surfaceFields.H"
33 
34 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
35 
37 (
38  const fvPatch& p,
40 )
41 :
42  fixedValueFvPatchScalarField(p, iF),
43  phiName_("phi"),
44  curTimeIndex_(-1)
45 {}
46 
47 
49 (
50  const fvPatch& p,
52  const dictionary& dict
53 )
54 :
55  fixedValueFvPatchScalarField(p, iF, dict, false),
56  Ap_(dict.get<scalar>("Ap")),
57  Sp_(dict.get<scalar>("Sp")),
58  VsI_(dict.get<scalar>("VsI")),
59  tas_(dict.get<scalar>("tas")),
60  tae_(dict.get<scalar>("tae")),
61  tds_(dict.get<scalar>("tds")),
62  tde_(dict.get<scalar>("tde")),
63  psI_(dict.get<scalar>("psI")),
64  psi_(dict.get<scalar>("psi")),
65  ams_(dict.get<scalar>("ams")),
66  ams0_(ams_),
67  phiName_(dict.lookupOrDefault<word>("phi", "phi")),
68  curTimeIndex_(-1)
69 {
70  scalar ps = (psI_*VsI_ + ams_/psi_)/Vs(db().time().value());
72 }
73 
74 
76 (
78  const fvPatch& p,
80  const fvPatchFieldMapper& mapper
81 )
82 :
83  fixedValueFvPatchScalarField(sppsf, p, iF, mapper),
84  Ap_(sppsf.Ap_),
85  Sp_(sppsf.Sp_),
86  VsI_(sppsf.VsI_),
87  tas_(sppsf.tas_),
88  tae_(sppsf.tae_),
89  tds_(sppsf.tds_),
90  tde_(sppsf.tde_),
91  psI_(sppsf.psI_),
92  psi_(sppsf.psi_),
93  ams_(sppsf.ams_),
94  ams0_(sppsf.ams0_),
95  phiName_(sppsf.phiName_),
96  curTimeIndex_(-1)
97 {}
98 
99 
101 (
104 )
105 :
106  fixedValueFvPatchScalarField(sppsf, iF),
107  Ap_(sppsf.Ap_),
108  Sp_(sppsf.Sp_),
109  VsI_(sppsf.VsI_),
110  tas_(sppsf.tas_),
111  tae_(sppsf.tae_),
112  tds_(sppsf.tds_),
113  tde_(sppsf.tde_),
114  psI_(sppsf.psI_),
115  psi_(sppsf.psi_),
116  ams_(sppsf.ams_),
117  ams0_(sppsf.ams0_),
118  phiName_(sppsf.phiName_),
119  curTimeIndex_(-1)
120 {}
121 
122 
124 (
126 )
127 :
128  fixedValueFvPatchScalarField(sppsf),
129  Ap_(sppsf.Ap_),
130  Sp_(sppsf.Sp_),
131  VsI_(sppsf.VsI_),
132  tas_(sppsf.tas_),
133  tae_(sppsf.tae_),
134  tds_(sppsf.tds_),
135  tde_(sppsf.tde_),
136  psI_(sppsf.psI_),
137  psi_(sppsf.psi_),
138  ams_(sppsf.ams_),
139  ams0_(sppsf.ams0_),
140  phiName_(sppsf.phiName_),
141  curTimeIndex_(-1)
142 {}
143 
144 
145 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
146 
147 Foam::scalar Foam::syringePressureFvPatchScalarField::Vs(const scalar t) const
148 {
149  if (t < tas_)
150  {
151  return VsI_;
152  }
153  else if (t < tae_)
154  {
155  return
156  VsI_
157  + 0.5*Ap_*Sp_*sqr(t - tas_)/(tae_ - tas_);
158  }
159  else if (t < tds_)
160  {
161  return
162  VsI_
163  + 0.5*Ap_*Sp_*(tae_ - tas_)
164  + Ap_*Sp_*(t - tae_);
165  }
166  else if (t < tde_)
167  {
168  return
169  VsI_
170  + 0.5*Ap_*Sp_*(tae_ - tas_)
171  + Ap_*Sp_*(tds_ - tae_)
172  + Ap_*Sp_*(t - tds_)
173  - 0.5*Ap_*Sp_*sqr(t - tds_)/(tde_ - tds_);
174  }
175  else
176  {
177  return
178  VsI_
179  + 0.5*Ap_*Sp_*(tae_ - tas_)
180  + Ap_*Sp_*(tds_ - tae_)
181  + 0.5*Ap_*Sp_*(tde_ - tds_);
182  }
183 }
184 
185 
187 {
188  if (updated())
189  {
190  return;
191  }
192 
193  if (curTimeIndex_ != db().time().timeIndex())
194  {
195  ams0_ = ams_;
196  curTimeIndex_ = db().time().timeIndex();
197  }
198 
199  scalar t = db().time().value();
200  scalar deltaT = db().time().deltaTValue();
201 
202  const surfaceScalarField& phi =
203  db().lookupObject<surfaceScalarField>(phiName_);
204 
205  const fvsPatchField<scalar>& phip =
206  patch().patchField<surfaceScalarField, scalar>(phi);
207 
208  if (phi.dimensions() == dimVelocity*dimArea)
209  {
210  ams_ = ams0_ + deltaT*sum((*this*psi_)*phip);
211  }
212  else if (phi.dimensions() == dimDensity*dimVelocity*dimArea)
213  {
214  ams_ = ams0_ + deltaT*sum(phip);
215  }
216  else
217  {
219  << "dimensions of phi are not correct"
220  << "\n on patch " << this->patch().name()
221  << " of field " << this->internalField().name()
222  << " in file " << this->internalField().objectPath()
223  << exit(FatalError);
224  }
225 
226  scalar ps = (psI_*VsI_ + ams_/psi_)/Vs(t);
227 
228  operator==(ps);
229 
230  fixedValueFvPatchScalarField::updateCoeffs();
231 }
232 
233 
235 {
237 
238  os.writeEntry("Ap", Ap_);
239  os.writeEntry("Sp", Sp_);
240  os.writeEntry("VsI", VsI_);
241  os.writeEntry("tas", tas_);
242  os.writeEntry("tae", tae_);
243  os.writeEntry("tds", tds_);
244  os.writeEntry("tde", tde_);
245  os.writeEntry("psI", psI_);
246  os.writeEntry("psi", psi_);
247  os.writeEntry("ams", ams_);
248  os.writeEntryIfDifferent<word>("phi", "phi", phiName_);
249 
250  writeEntry("value", os);
251 }
252 
253 
254 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
255 
256 namespace Foam
257 {
259  (
262  );
263 }
264 
265 // ************************************************************************* //
Foam::fvPatchField
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: volSurfaceMapping.H:50
Foam::fvPatchField< scalar >::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::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:62
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:234
surfaceFields.H
Foam::surfaceFields.
fvPatchFieldMapper.H
Foam::dimArea
const dimensionSet dimArea(sqr(dimLength))
Definition: dimensionSets.H:60
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:63
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:37
timeIndex
label timeIndex
Definition: getTimeIndex.H:4
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::syringePressureFvPatchScalarField::updateCoeffs
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Definition: syringePressureFvPatchScalarField.C:186
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:355
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:219
syringePressureFvPatchScalarField.H
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::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