SRFVelocityFvPatchVectorField.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-2016 OpenFOAM Foundation
6  \\/ M anipulation |
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 
28 #include "volFields.H"
29 
30 #include "SRFModel.H"
31 
32 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
33 
35 (
36  const fvPatch& p,
38 )
39 :
40  fixedValueFvPatchVectorField(p, iF),
41  relative_(0),
42  inletValue_(p.size(), Zero)
43 {}
44 
45 
47 (
49  const fvPatch& p,
51  const fvPatchFieldMapper& mapper
52 )
53 :
54  fixedValueFvPatchVectorField(ptf, p, iF, mapper),
55  relative_(ptf.relative_),
56  inletValue_(ptf.inletValue_, mapper)
57 {}
58 
59 
61 (
62  const fvPatch& p,
64  const dictionary& dict
65 )
66 :
67  fixedValueFvPatchVectorField(p, iF),
68  relative_(dict.lookup("relative")),
69  inletValue_("inletValue", dict, p.size())
70 {
71  fvPatchVectorField::operator=(vectorField("value", dict, p.size()));
72 }
73 
74 
76 (
77  const SRFVelocityFvPatchVectorField& srfvpvf
78 )
79 :
80  fixedValueFvPatchVectorField(srfvpvf),
81  relative_(srfvpvf.relative_),
82  inletValue_(srfvpvf.inletValue_)
83 {}
84 
85 
87 (
88  const SRFVelocityFvPatchVectorField& srfvpvf,
90 )
91 :
92  fixedValueFvPatchVectorField(srfvpvf, iF),
93  relative_(srfvpvf.relative_),
94  inletValue_(srfvpvf.inletValue_)
95 {}
96 
97 
98 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
99 
101 (
102  const fvPatchFieldMapper& m
103 )
104 {
105  vectorField::autoMap(m);
106  inletValue_.autoMap(m);
107 }
108 
109 
111 (
112  const fvPatchVectorField& ptf,
113  const labelList& addr
114 )
115 {
116  fixedValueFvPatchVectorField::rmap(ptf, addr);
117 
118  const SRFVelocityFvPatchVectorField& tiptf =
119  refCast<const SRFVelocityFvPatchVectorField>(ptf);
120 
121  inletValue_.rmap(tiptf.inletValue_, addr);
122 }
123 
124 
126 {
127  if (updated())
128  {
129  return;
130  }
131 
132  // If not relative to the SRF include the effect of the SRF
133  if (!relative_)
134  {
135  // Get reference to the SRF model
136  const SRF::SRFModel& srf =
137  db().lookupObject<SRF::SRFModel>("SRFProperties");
138 
139  // Determine patch velocity due to SRF
140  const vectorField SRFVelocity(srf.velocity(patch().Cf()));
141 
142  operator==(-SRFVelocity + inletValue_);
143  }
144  // If already relative to the SRF simply supply the inlet value as a fixed
145  // value
146  else
147  {
148  operator==(inletValue_);
149  }
150 
151  fixedValueFvPatchVectorField::updateCoeffs();
152 }
153 
154 
156 {
158  os.writeKeyword("relative") << relative_ << token::END_STATEMENT << nl;
159  inletValue_.writeEntry("inletValue", os);
160  writeEntry("value", os);
161 }
162 
163 
164 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
165 
166 namespace Foam
167 {
169  (
172  );
173 }
174 
175 // ************************************************************************* //
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
SRFVelocityFvPatchVectorField(const fvPatch &, const DimensionedField< vector, volMesh > &)
Construct from patch and internal field.
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:61
void writeEntry(const word &keyword, Ostream &os) const
Write the field as a dictionary entry.
Definition: Field.C:719
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:65
volVectorField vectorField(fieldObject, mesh)
Macros for easy insertion into run-time selection tables.
virtual void autoMap(const fvPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
makePatchTypeField(fvPatchVectorField, SRFFreestreamVelocityFvPatchVectorField)
vectorField velocity(const vectorField &positions) const
Return velocity vector from positions.
Definition: SRFModel.C:171
virtual void write(Ostream &) const
Write.
tmp< fvMatrix< Type > > operator==(const fvMatrix< Type > &, const fvMatrix< Type > &)
Foam::fvPatchFieldMapper.
virtual label size() const
Return size.
Definition: fvPatch.H:161
static const zero Zero
Definition: zero.H:91
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:53
static const char nl
Definition: Ostream.H:262
Ostream & writeKeyword(const keyType &)
Write the keyword followed by an appropriate indentation.
Definition: Ostream.C:54
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
virtual void rmap(const fvPatchVectorField &, const labelList &)
Reverse map the given fvPatchField onto this fvPatchField.
Top level model for single rotating frame.
Definition: SRFModel.H:62
Namespace for OpenFOAM.
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:363
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:451
Freestream velocity condition to be used in conjunction with the single rotating frame (SRF) model (s...