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 | www.openfoam.com
6 \\/ M anipulation |
7-------------------------------------------------------------------------------
8 Copyright (C) 2011-2016 OpenFOAM Foundation
9-------------------------------------------------------------------------------
10License
11 This file is part of OpenFOAM.
12
13 OpenFOAM is free software: you can redistribute it and/or modify it
14 under the terms of the GNU General Public License as published by
15 the Free Software Foundation, either version 3 of the License, or
16 (at your option) any later version.
17
18 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21 for more details.
22
23 You should have received a copy of the GNU General Public License
24 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25
26\*---------------------------------------------------------------------------*/
27
30#include "volFields.H"
31
32#include "SRFModel.H"
33
34// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
35
37(
38 const fvPatch& p,
40)
41:
42 fixedValueFvPatchVectorField(p, iF),
43 relative_(false),
44 inletValue_(p.size(), Zero)
45{}
46
47
49(
51 const fvPatch& p,
53 const fvPatchFieldMapper& mapper
54)
55:
56 fixedValueFvPatchVectorField(ptf, p, iF, mapper),
57 relative_(ptf.relative_),
58 inletValue_(ptf.inletValue_, mapper)
59{}
60
61
63(
64 const fvPatch& p,
66 const dictionary& dict
67)
68:
69 fixedValueFvPatchVectorField(p, iF, dict),
70 relative_(dict.get<Switch>("relative")),
71 inletValue_("inletValue", dict, p.size())
72{}
73
74
76(
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{
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.writeEntry("relative", relative_);
159 inletValue_.writeEntry("inletValue", os);
160 writeEntry("value", os);
161}
162
163
164// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
165
166namespace Foam
167{
169 (
172 );
173}
174
175// ************************************************************************* //
Macros for easy insertion into run-time selection tables.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
void autoMap(const FieldMapper &map, const bool applyFlip=true)
Map from self.
Definition: Field.C:403
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:62
Ostream & writeEntry(const keyType &key, const T &value)
Write a keyword/value entry.
Definition: Ostream.H:239
Velocity condition to be used in conjunction with the single rotating frame (SRF) model (see: SRFMode...
virtual void autoMap(const fvPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
virtual void rmap(const fvPatchVectorField &, const labelList &)
Reverse map the given fvPatchField onto this fvPatchField.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Top level model for single rotating frame.
Definition: SRFModel.H:67
vectorField velocity(const vectorField &positions) const
Return velocity vector from positions.
Definition: SRFModel.C:168
A simple wrapper around bool so that it can be read as a word: true/false, on/off,...
Definition: Switch.H:78
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
virtual bool write()
Write the output fields.
A FieldMapper for finite-volume patch fields.
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:71
volScalarField & p
OBJstream os(runTime.globalPath()/outputName)
#define makePatchTypeField(PatchTypeField, typePatchTypeField)
Definition: fvPatchField.H:676
Namespace for OpenFOAM.
tmp< faMatrix< Type > > operator==(const faMatrix< Type > &, const faMatrix< Type > &)
dictionary dict