All Classes Namespaces Files Functions Variables Typedefs Friends Macros Groups
cosVelocityFvPatchVectorField.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) 2019 OpenCFD Ltd.
10  Copyright (C) 2016-2019 ISP RAS (www.ispras.ru) UniCFD Group (www.unicfd.ru)
11 -------------------------------------------------------------------------------
12 License
13  This file is part of QGDsolver library, based on OpenFOAM+.
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 
30 #include "addToRunTimeSelectionTable.H"
31 #include "fvPatchFieldMapper.H"
32 #include "volFields.H"
33 #include "surfaceFields.H"
34 
35 
36 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
37 
39 (
40  const fvPatch& p,
41  const DimensionedField<vector, volMesh>& iF
42 )
43 :
44  mixedFvPatchVectorField(p, iF),
45  A_(0.0, 0.0, 0.0),
46  H_(1.0),
47  Hdirection_(0.0, 0.0, 0.0),
48  omega0_(0.0),
49  phi0_(0.0),
50  minZ_(0.0)
51 {
52  this->refValue() = pTraits<vector>::zero;
53  this->refGrad() = pTraits<vector>::zero;
54  this->valueFraction() = 0.0;
55  this->setPatchVelocities(this->refValue());
56 
57  Info << "Executing cosVelocityFvPatchVectorField(fvPatch, DimensionedField) " << endl;
58 }
59 
60 
62 (
64  const fvPatch& p,
65  const DimensionedField<vector, volMesh>& iF,
66  const fvPatchFieldMapper& mapper
67 )
68 :
69  mixedFvPatchField<vector>(ptf, p, iF, mapper),
70  A_(ptf.A_),
71  H_(ptf.H_),
72  Hdirection_(ptf.Hdirection_),
73  omega0_(ptf.omega0_),
74  phi0_(ptf.phi0_),
75  minZ_(ptf.minZ_)
76 {
77  this->setPatchVelocities(this->refValue());
78  mixedFvPatchVectorField::updateCoeffs();
79 
80  Info << "Executing cosVelocityFvPatchVectorField(cosVelocity, fvPatch, DimensionedField, fvPatchFieldMapper) " << endl;
81 }
82 
83 
85 (
86  const fvPatch& p,
87  const DimensionedField<vector, volMesh>& iF,
88  const dictionary& dict
89 )
90 :
91  mixedFvPatchField<vector>(p, iF)
92 {
93 
94  dict.lookup("A") >> A_;
95  dict.lookup("H") >> H_;
96  dict.lookup("Hdirection") >> Hdirection_;
97  dict.lookup("omega0") >> omega0_;
98  dict.lookup("phi0") >> phi0_;
99  dict.lookup("minZ") >> minZ_;
100 
101  if (dict.found("value"))
102  {
103  fvPatchField<vector>::operator=
104  (
105  Field<vector>("value", dict, p.size())
106  );
107  }
108  else
109  {
110  fvPatchField<vector>::operator=(this->patchInternalField());
111  }
112 
113  this->refValue() = pTraits<vector>::zero;
114  this->refGrad() = pTraits<vector>::zero;
115  this->valueFraction() = 1.0;
116  this->setPatchVelocities(this->refValue());
117  mixedFvPatchVectorField::updateCoeffs();
118 
119  Info << "Executing cosVelocityFvPatchVectorField(fvPatch, DimensionedField, dictionary) " << endl;
120 }
121 
122 
124 (
125  const cosVelocityFvPatchVectorField& ptpsf
126 )
127 :
128  mixedFvPatchVectorField(ptpsf),
129  A_(ptpsf.A_),
130  H_(ptpsf.H_),
131  Hdirection_(ptpsf.Hdirection_),
132  omega0_(ptpsf.omega0_),
133  phi0_(ptpsf.phi0_),
134  minZ_(ptpsf.minZ_)
135 {
136  this->setPatchVelocities(this->refValue());
137  Info << "Executing cosVelocityFvPatchVectorField(const cosVelocity&) " << endl;
138 }
139 
140 
142 (
143  const cosVelocityFvPatchVectorField& ptpsf,
144  const DimensionedField<vector, volMesh>& iF
145 )
146 :
147  mixedFvPatchVectorField(ptpsf, iF),
148  A_(ptpsf.A_),
149  H_(ptpsf.H_),
150  Hdirection_(ptpsf.Hdirection_),
151  omega0_(ptpsf.omega0_),
152  phi0_(ptpsf.phi0_),
153  minZ_(ptpsf.minZ_)
154 {
155  this->setPatchVelocities(this->refValue());
157  Info << "Executing cosVelocityFvPatchVectorField(const cosVelocity, const DimensionedField) " << endl;
158 }
159 
160 
161 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
162 
163 void Foam::cosVelocityFvPatchVectorField::autoMap(const fvPatchFieldMapper& m)
164 {
165  mixedFvPatchField<vector>::autoMap(m);
166 
167  Info << "Executing autoMap " << endl;
168 }
170 void Foam::cosVelocityFvPatchVectorField::rmap(const fvPatchVectorField& pf, const labelList& ll)
171 {
172  mixedFvPatchField<vector>::rmap(pf, ll);
173  Info << "Executing rmap " << endl;
174 }
175 
177 {
178  scalarField z = (this->patch().Cf() & Hdirection_) - minZ_;
179 
180  forAll(pV, iFace)
181  {
182  scalar t = this->patch().boundaryMesh().mesh().time().value();
183  this->refGrad()[iFace] = vector (0.0, 0.0, 0.0);
184  this->valueFraction()[iFace] = 1.0;
185  pV[iFace] =
186  A_ * cos(Foam::constant::mathematical::pi * z[iFace] / H_) * (-omega0_) * sin(omega0_ * t + phi0_);
187  if (z[iFace] < 0.0 || z[iFace] > H_)
188  {
189  pV[iFace] = A_*0.0;
190  }
191  }
192 }
193 
195 {
196  if (this->updated())
197  {
198  return;
199  }
200 
201  this->setPatchVelocities(this->refValue());
202 
203  Info << "max/min velocity: " << gMax(this->refValue()) << "/" << gMin(this->refValue()) << endl;
204 
205  mixedFvPatchField<vector>::updateCoeffs();
206 }
207 
208 
209 void Foam::cosVelocityFvPatchVectorField::write(Ostream& os) const
210 {
211 
212  mixedFvPatchVectorField::write(os);
213 
214  os.writeKeyword("A") << A_ << token::END_STATEMENT << nl;
215 
216  os.writeKeyword("H") << H_ << token::END_STATEMENT << nl;
217 
218  os.writeKeyword("minZ") << minZ_ << token::END_STATEMENT << nl;
219 
220  os.writeKeyword("Hdirection") << Hdirection_ << token::END_STATEMENT << nl;
221 
222  os.writeKeyword("omega0") << omega0_ << token::END_STATEMENT << nl;
223 
224  os.writeKeyword("phi0") << phi0_ << token::END_STATEMENT << nl;
225 
226 }
227 
228 namespace Foam
229 {
231  (
232  fvPatchVectorField,
233  cosVelocityFvPatchVectorField
234  );
235 }
236 
237 // ************************************************************************* //
makePatchTypeField(fvPatchVectorField, cosVelocityFvPatchVectorField)
pf
Definition: updateFields.H:21
virtual void rmap(const fvPatchVectorField &pf, const labelList &ll)
This boundary condition sets the pressure gradient to the provided value such that the flux on the bo...
forAll(Y, i)
Definition: QGDYEqn.H:36
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
virtual void write(Ostream &) const
Write.
cosVelocityFvPatchVectorField(const fvPatch &, const DimensionedField< vector, volMesh > &)
Construct from patch and internal field.
volScalarField & p
Definition: createFields.H:52
virtual void autoMap(const fvPatchFieldMapper &)