varyingVelocityInletOutletFvPatchVectorField.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 | Copyright (C) 2011-2017 OpenFOAM Foundation
6  \\/ M anipulation | Copyright (C) 2017 OpenCFD Ltd.
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 
28 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
29 
32  const fvPatch& p,
33  const DimensionedField<vector, volMesh>& iF)
34  : mixedFvPatchField<vector>(p, iF),
35  phiName_("phi"),
36  U0_(0.0),
37  URate_(0.0),
38  flowComponent_(0),
39  normalComponent_(1),
40  alpha0_(0),
41  alphaRate_(0)
42 {
43  forAll(this->refValue(), idxI)
44  {
45  this->refValue()[idxI] = pTraits<vector>::zero;
46  }
47  this->refGrad() = pTraits<vector>::zero;
48  this->valueFraction() = 0.0;
49 }
50 
54  const fvPatch& p,
55  const DimensionedField<vector, volMesh>& iF,
56  const fvPatchFieldMapper& mapper)
57  : mixedFvPatchField<vector>(ptf, p, iF, mapper),
58  phiName_(ptf.phiName_),
59  U0_(ptf.U0_),
60  URate_(ptf.URate_),
61  flowComponent_(ptf.flowComponent_),
62  normalComponent_(ptf.normalComponent_),
63  alpha0_(ptf.alpha0_),
64  alphaRate_(ptf.alphaRate_)
65 {
66 }
67 
70  const fvPatch& p,
71  const DimensionedField<vector, volMesh>& iF,
72  const dictionary& dict)
73  : mixedFvPatchField<vector>(p, iF),
74  phiName_(dict.lookupOrDefault<word>("phi", "phi")),
75  U0_(dict.lookupOrDefault<scalar>("U0", 0.0)),
76  URate_(dict.lookupOrDefault<scalar>("URate", 0.0)),
77  flowComponent_(dict.lookupOrDefault<label>("flowComponent", 0)),
78  normalComponent_(dict.lookupOrDefault<label>("normalComponent", 1)),
79  alpha0_(dict.lookupOrDefault<scalar>("alpha0", 0.0)),
80  alphaRate_(dict.lookupOrDefault<scalar>("alphaRate", 0.0))
81 {
82  this->patchType() = dict.lookupOrDefault<word>("patchType", word::null);
83 
84  vector UInit = vector::zero;
85  UInit[flowComponent_] = U0_ * cos(alpha0_);
86  UInit[normalComponent_] = U0_ * sin(alpha0_);
87  forAll(this->refValue(), idxI)
88  {
89  this->refValue()[idxI] = UInit;
90  }
91 
92  if (dict.found("value"))
93  {
94  fvPatchField<vector>::operator=(
95  vectorField("value", dict, p.size()));
96  }
97  else
98  {
99  fvPatchField<vector>::operator=(this->refValue());
100  }
101 
102  this->refGrad() = pTraits<vector>::zero;
103  this->valueFraction() = 0.0;
104 }
105 
109  : mixedFvPatchField<vector>(ptf),
110  phiName_(ptf.phiName_),
111  U0_(ptf.U0_),
112  URate_(ptf.URate_),
113  flowComponent_(ptf.flowComponent_),
114  normalComponent_(ptf.normalComponent_),
115  alpha0_(ptf.alpha0_),
116  alphaRate_(ptf.alphaRate_)
117 {
118 }
119 
123  const DimensionedField<vector, volMesh>& iF)
124  : mixedFvPatchField<vector>(ptf, iF),
125  phiName_(ptf.phiName_),
126  U0_(ptf.U0_),
127  URate_(ptf.URate_),
128  flowComponent_(ptf.flowComponent_),
129  normalComponent_(ptf.normalComponent_),
130  alpha0_(ptf.alpha0_),
131  alphaRate_(ptf.alphaRate_)
132 {
133 }
134 
135 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
136 
138 {
139  if (updated())
140  {
141  return;
142  }
143 
144  const scalar t = this->db().time().timeOutputValue();
145  scalar alpha = alpha0_ + t * alphaRate_;
146  scalar U = U0_ + t * URate_;
147 
148  const Field<scalar>& phip =
149  this->patch().template lookupPatchField<surfaceScalarField, scalar>(phiName_);
150 
151  forAll(this->refValue(), faceI)
152  {
153  this->refValue()[faceI][flowComponent_] = U * cos(alpha) * (1.0 - pos0(phip[faceI]));
154  this->refValue()[faceI][normalComponent_] = U * sin(alpha) * (1.0 - pos0(phip[faceI]));
155  }
156 
157  this->valueFraction() = 1.0 - pos0(phip);
158 
159  mixedFvPatchField<vector>::updateCoeffs();
160 }
161 
163  const
164 {
165  fvPatchVectorField::write(os);
166  os.writeEntry("phiName", phiName_);
167  os.writeEntry("U0", U0_);
168  os.writeEntry("URate", URate_);
169  os.writeEntry("flowComponent", flowComponent_);
170  os.writeEntry("normalComponent", normalComponent_);
171  os.writeEntry("alpha0", alpha0_);
172  os.writeEntry("alphaRate", alphaRate_);
173  //writeEntry("value", os);
174 }
175 
176 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
177 
179  const fvPatchVectorField& ptf)
180 {
181  fvPatchVectorField::operator=(
182  this->valueFraction() * this->refValue()
183  + (1 - this->valueFraction()) * ptf);
184 }
185 
186 namespace Foam
187 {
189  fvPatchVectorField,
191 }
192 
193 // ************************************************************************* //
U
volVectorField & U
Definition: createRefsPimpleDyM.H:7
Foam::varyingVelocityInletOutletFvPatchVectorField::operator=
virtual void operator=(const fvPatchVectorField &pvf)
Definition: varyingVelocityInletOutletFvPatchVectorField.C:178
Foam::varyingVelocityInletOutletFvPatchVectorField::write
virtual void write(Ostream &) const
Definition: varyingVelocityInletOutletFvPatchVectorField.C:162
varyingVelocityInletOutletFvPatchVectorField.H
Foam::varyingVelocityInletOutletFvPatchVectorField::varyingVelocityInletOutletFvPatchVectorField
varyingVelocityInletOutletFvPatchVectorField(const fvPatch &, const DimensionedField< vector, volMesh > &)
Definition: varyingVelocityInletOutletFvPatchVectorField.C:31
p
volScalarField & p
Definition: createRefsPimpleDyM.H:6
Foam::varyingVelocityInletOutletFvPatchVectorField
Definition: varyingVelocityInletOutletFvPatchVectorField.H:82
Foam
Definition: checkGeometry.C:32
forAll
forAll(nuTilda1, cellI)
Definition: nuTilda1EqnIrkPimple.H:19
Foam::varyingVelocityInletOutletFvPatchVectorField::updateCoeffs
virtual void updateCoeffs()
Definition: varyingVelocityInletOutletFvPatchVectorField.C:137
Foam::makePatchTypeField
makePatchTypeField(fvPatchScalarField, fixedWallHeatFluxFvPatchScalarField)