SRFVelocityFvPatchVectorField.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2022 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 (
48  const fvPatch& p,
50  const dictionary& dict
51 )
52 :
53  fixedValueFvPatchVectorField(p, iF, dict),
54  relative_(dict.lookup("relative")),
55  inletValue_("inletValue", dict, p.size())
56 {}
57 
58 
60 (
62  const fvPatch& p,
64  const fvPatchFieldMapper& mapper
65 )
66 :
67  fixedValueFvPatchVectorField(ptf, p, iF, mapper),
68  relative_(ptf.relative_),
69  inletValue_(mapper(ptf.inletValue_))
70 {}
71 
72 
74 (
75  const SRFVelocityFvPatchVectorField& srfvpvf,
77 )
78 :
79  fixedValueFvPatchVectorField(srfvpvf, iF),
80  relative_(srfvpvf.relative_),
81  inletValue_(srfvpvf.inletValue_)
82 {}
83 
84 
85 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
86 
88 (
89  const fvPatchFieldMapper& m
90 )
91 {
92  m(*this, *this);
93  m(inletValue_, inletValue_);
94 }
95 
96 
98 (
99  const fvPatchVectorField& ptf,
100  const labelList& addr
101 )
102 {
103  fixedValueFvPatchVectorField::rmap(ptf, addr);
104 
105  const SRFVelocityFvPatchVectorField& tiptf =
106  refCast<const SRFVelocityFvPatchVectorField>(ptf);
107 
108  inletValue_.rmap(tiptf.inletValue_, addr);
109 }
110 
111 
113 (
114  const fvPatchVectorField& ptf
115 )
116 {
117  fixedValueFvPatchVectorField::reset(ptf);
118 
119  const SRFVelocityFvPatchVectorField& tiptf =
120  refCast<const SRFVelocityFvPatchVectorField>(ptf);
121 
122  inletValue_.reset(tiptf.inletValue_);
123 }
124 
125 
127 {
128  if (updated())
129  {
130  return;
131  }
132 
133  // If not relative to the SRF include the effect of the SRF
134  if (!relative_)
135  {
136  // Get reference to the SRF model
137  const SRF::SRFModel& srf =
138  db().lookupObject<SRF::SRFModel>("SRFProperties");
139 
140  // Determine patch velocity due to SRF
141  const vectorField SRFVelocity(srf.velocity(patch().Cf()));
142 
143  operator==(-SRFVelocity + inletValue_);
144  }
145  // If already relative to the SRF simply supply the inlet value as a fixed
146  // value
147  else
148  {
149  operator==(inletValue_);
150  }
151 
152  fixedValueFvPatchVectorField::updateCoeffs();
153 }
154 
155 
157 {
159  writeEntry(os, "relative", relative_);
160  writeEntry(os, "inletValue", inletValue_);
161  writeEntry(os, "value", *this);
162 }
163 
164 
165 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
166 
167 namespace Foam
168 {
170  (
173  );
174 }
175 
176 // ************************************************************************* //
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
SRFVelocityFvPatchVectorField(const fvPatch &, const DimensionedField< vector, volMesh > &)
Construct from patch and internal field.
virtual void reset(const fvPatchVectorField &)
Reset the fvPatchField to the given fvPatchField.
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:63
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:66
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:243
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.
tmp< fvMatrix< Type > > operator==(const fvMatrix< Type > &, const fvMatrix< Type > &)
Foam::fvPatchFieldMapper.
static const zero Zero
Definition: zero.H:97
virtual label size() const
Return size.
Definition: fvPatch.H:157
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:54
virtual void write(Ostream &) const
Write.
void writeEntry(Ostream &os, const HashTable< T, Key, Hash > &ht)
Definition: HashTableIO.C:96
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
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
vectorField velocity(const vectorField &positions) const
Return velocity vector from positions.
Definition: SRFModel.C:151
Namespace for OpenFOAM.
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:864
Velocity condition to be used in conjunction with the single rotating frame (SRF) model (see: SRFMode...