varyingVelocityFvPatchVectorField.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2 
3  DAFoam : Discrete Adjoint with OpenFOAM
4  Version : v4
5 
6  This file is modified from OpenFOAM's source code
7  src/finiteVolume/fields/fvPatchFields/basic/fixedValue/fixedValueFvPatchField.C
8 
9  OpenFOAM: The Open Source CFD Toolbox
10 
11  Copyright (C): 2011-2016 OpenFOAM Foundation
12 
13  OpenFOAM License:
14 
15  OpenFOAM is free software: you can redistribute it and/or modify it
16  under the terms of the GNU General Public License as published by
17  the Free Software Foundation, either version 3 of the License, or
18  (at your option) any later version.
19 
20  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
21  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
22  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
23  for more details.
24 
25  You should have received a copy of the GNU General Public License
26  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
27 
28 \*---------------------------------------------------------------------------*/
29 
31 #include "volFields.H"
32 #include "addToRunTimeSelectionTable.H"
33 
34 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
35 
37  const fvPatch& p,
38  const DimensionedField<vector, volMesh>& iF)
39  : fixedValueFvPatchVectorField(p, iF),
40  U0_(0.0),
41  URate_(0.0),
42  flowComponent_(0),
43  normalComponent_(1),
44  alpha0_(0),
45  alphaRate_(0)
46 {
47 }
48 
50  const fvPatch& p,
51  const DimensionedField<vector, volMesh>& iF,
52  const dictionary& dict)
53  : fixedValueFvPatchVectorField(p, iF),
54  U0_(dict.lookupOrDefault<scalar>("U0", 0.0)),
55  URate_(dict.lookupOrDefault<scalar>("URate", 0.0)),
56  flowComponent_(dict.lookupOrDefault<label>("flowComponent", 0)),
57  normalComponent_(dict.lookupOrDefault<label>("normalComponent", 1)),
58  alpha0_(dict.lookupOrDefault<scalar>("alpha0", 0.0)),
59  alphaRate_(dict.lookupOrDefault<scalar>("alphaRate", 0.0))
60 {
61 
62  if (dict.found("value"))
63  {
64  fvPatchVectorField::operator=(
65  vectorField("value", dict, p.size()));
66  }
67  else
68  {
69  this->evaluate();
70  }
71 }
72 
75  const fvPatch& p,
76  const DimensionedField<vector, volMesh>& iF,
77  const fvPatchFieldMapper& mapper)
78  : fixedValueFvPatchVectorField(p, iF),
79  U0_(ptf.U0_),
80  URate_(ptf.URate_),
81  flowComponent_(ptf.flowComponent_),
82  normalComponent_(ptf.normalComponent_),
83  alpha0_(ptf.alpha0_),
84  alphaRate_(ptf.alphaRate_)
85 {
86  this->evaluate();
87 }
88 
91  : fixedValueFvPatchVectorField(wbppsf),
92  U0_(wbppsf.U0_),
93  URate_(wbppsf.URate_),
94  flowComponent_(wbppsf.flowComponent_),
95  normalComponent_(wbppsf.normalComponent_),
96  alpha0_(wbppsf.alpha0_),
97  alphaRate_(wbppsf.alphaRate_)
98 {
99  this->evaluate();
100 }
101 
103  const varyingVelocityFvPatchVectorField& wbppsf,
104  const DimensionedField<vector, volMesh>& iF)
105  : fixedValueFvPatchVectorField(wbppsf, iF),
106  U0_(wbppsf.U0_),
107  URate_(wbppsf.URate_),
108  flowComponent_(wbppsf.flowComponent_),
109  normalComponent_(wbppsf.normalComponent_),
110  alpha0_(wbppsf.alpha0_),
111  alphaRate_(wbppsf.alphaRate_)
112 {
113  this->evaluate();
114 }
115 
116 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
118 {
119  // calculate patch values
120  const scalar t = this->db().time().timeOutputValue();
121  scalar alpha = alpha0_ + t * alphaRate_;
122  scalar U = U0_ + t * URate_;
123 
124  vectorField& thisPatchRef = *this;
125  vectorField thisPatch = thisPatchRef;
126  forAll(thisPatch, idxI)
127  {
128  thisPatch[idxI][flowComponent_] = U * cos(alpha);
129  thisPatch[idxI][normalComponent_] = U * sin(alpha);
130  }
131 
132  fvPatchVectorField::operator==(thisPatch);
133 
134  fixedValueFvPatchVectorField::updateCoeffs();
135 }
136 
138 {
139  fixedValueFvPatchVectorField::write(os);
140  os.writeEntry("U0", U0_);
141  os.writeEntry("URate", URate_);
142  os.writeEntry("flowComponent", flowComponent_);
143  os.writeEntry("normalComponent", normalComponent_);
144  os.writeEntry("alpha0", alpha0_);
145  os.writeEntry("alphaRate", alphaRate_);
146  //writeEntry("value", os);
147 }
148 
149 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
150 
151 namespace Foam
152 {
154  fvPatchVectorField,
156 }
157 
158 // ************************************************************************* //
U
volVectorField & U
Definition: createRefsPimpleDyM.H:7
Foam::varyingVelocityFvPatchVectorField::updateCoeffs
virtual void updateCoeffs()
Definition: varyingVelocityFvPatchVectorField.C:117
varyingVelocityFvPatchVectorField.H
Foam::varyingVelocityFvPatchVectorField::write
virtual void write(Ostream &) const
Definition: varyingVelocityFvPatchVectorField.C:137
p
volScalarField & p
Definition: createRefsPimpleDyM.H:6
Foam
Definition: checkGeometry.C:32
Foam::varyingVelocityFvPatchVectorField
Definition: varyingVelocityFvPatchVectorField.H:78
forAll
forAll(nuTilda1, cellI)
Definition: nuTilda1EqnIrkPimple.H:19
Foam::makePatchTypeField
makePatchTypeField(fvPatchScalarField, fixedWallHeatFluxFvPatchScalarField)
Foam::varyingVelocityFvPatchVectorField::varyingVelocityFvPatchVectorField
varyingVelocityFvPatchVectorField(const fvPatch &, const DimensionedField< vector, volMesh > &)
Definition: varyingVelocityFvPatchVectorField.C:36