oscillatingVelocityPointPatchVectorField.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
29#include "pointPatchFields.H"
31#include "Time.H"
32#include "polyMesh.H"
33
34// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
35
36namespace Foam
37{
38
39// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
40
41oscillatingVelocityPointPatchVectorField::
42oscillatingVelocityPointPatchVectorField
43(
44 const pointPatch& p,
46)
47:
49 amplitude_(Zero),
50 omega_(0.0),
51 p0_(p.localPoints())
52{}
53
54
55oscillatingVelocityPointPatchVectorField::
56oscillatingVelocityPointPatchVectorField
57(
58 const pointPatch& p,
60 const dictionary& dict
61)
62:
64 amplitude_(dict.lookup("amplitude")),
65 omega_(dict.get<scalar>("omega"))
66{
67 if (!dict.found("value"))
68 {
70 }
71
72 if (dict.found("p0"))
73 {
74 p0_ = vectorField("p0", dict , p.size());
75 }
76 else
77 {
78 p0_ = p.localPoints();
79 }
80}
81
82
83oscillatingVelocityPointPatchVectorField::
84oscillatingVelocityPointPatchVectorField
85(
87 const pointPatch& p,
89 const pointPatchFieldMapper& mapper
90)
91:
92 fixedValuePointPatchField<vector>(ptf, p, iF, mapper),
93 amplitude_(ptf.amplitude_),
94 omega_(ptf.omega_),
95 p0_(ptf.p0_, mapper)
96{}
97
98
99oscillatingVelocityPointPatchVectorField::
100oscillatingVelocityPointPatchVectorField
101(
104)
105:
107 amplitude_(ptf.amplitude_),
108 omega_(ptf.omega_),
109 p0_(ptf.p0_)
110{}
111
112
113// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
114
116(
117 const pointPatchFieldMapper& m
118)
119{
121
122 p0_.autoMap(m);
123}
124
125
127(
128 const pointPatchField<vector>& ptf,
129 const labelList& addr
130)
131{
133 refCast<const oscillatingVelocityPointPatchVectorField>(ptf);
134
136
137 p0_.rmap(oVptf.p0_, addr);
138}
139
140
142{
143 if (this->updated())
144 {
145 return;
146 }
147
148 const polyMesh& mesh = this->internalField().mesh()();
149 const Time& t = mesh.time();
150 const pointPatch& p = this->patch();
151
153 (
154 (p0_ + amplitude_*sin(omega_*t.value()) - p.localPoints())
155 /t.deltaTValue()
156 );
157
159}
160
161
163{
165 os.writeEntry("amplitude", amplitude_);
166 os.writeEntry("omega", omega_);
167 p0_.writeEntry("p0", os);
168 writeEntry("value", os);
169}
170
171
172// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
173
175(
178);
179
180// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
181
182} // End namespace Foam
183
184// ************************************************************************* //
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...
const Mesh & mesh() const
Return mesh.
Generic templated field type.
Definition: Field.H:82
void autoMap(const FieldMapper &map, const bool applyFlip=true)
Map from self.
Definition: Field.C:403
void writeEntry(const word &keyword, Ostream &os) const
Write the field as a dictionary entry.
Definition: Field.C:614
void rmap(const UList< Type > &mapF, const labelUList &mapAddressing)
1 to 1 reverse-map from the given field
Definition: Field.C:466
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
scalar deltaTValue() const noexcept
Return time step value.
Definition: TimeStateI.H:43
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:80
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:114
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
void autoMap(const fvPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
void rmap(const atmBoundaryLayer &, const labelList &)
Reverse map the given fvPatchField onto this fvPatchField.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
bool found(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Search for an entry (const access) with the given keyword.
Definition: dictionaryI.H:87
const Type & value() const
Return const reference to value.
A FixedValue boundary condition for pointField.
virtual bool write()
Write the output fields.
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:290
virtual void autoMap(const pointPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
virtual void rmap(const pointPatchField< vector > &, const labelList &)
Reverse map the given pointPatchField onto this pointPatchField.
Foam::pointPatchFieldMapper.
Abstract base class for point-mesh patch fields.
const pointPatch & patch() const
Return patch.
const DimensionedField< Type, pointMesh > & internalField() const
Return dimensioned internal field reference.
bool updated() const
Return true if the boundary condition has already been updated.
Basic pointPatch represents a set of points from the mesh.
Definition: pointPatch.H:64
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:81
Lookup type of boundary radiation properties.
Definition: lookup.H:66
volScalarField & p
dynamicFvMesh & mesh
OBJstream os(runTime.globalPath()/outputName)
Namespace for OpenFOAM.
dimensionedScalar sin(const dimensionedScalar &ds)
Field< vector > vectorField
Specialisation of Field<T> for vector.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
#define makePointPatchTypeField(PatchTypeField, typePatchTypeField)
dictionary dict