Loading...
Searching...
No Matches
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
43(
44 const pointPatch& p,
46)
47:
49 amplitude_(Zero),
50 omega_(0.0),
51 p0_(p.localPoints())
52{}
53
54
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 {
69 updateCoeffs();
70 }
71
72 if (dict.found("p0"))
73 {
74 p0_ = vectorField("p0", dict , p.size());
75 }
90)
91:
93 amplitude_(ptf.amplitude_),
94 omega_(ptf.omega_),
95 p0_(ptf.p0_, mapper)
96{}
97
98
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)
121
122 p0_.autoMap(m);
123}
124
125
127(
128 const pointPatchField<vector>& ptf,
129 const labelList& addr
130)
131{
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
152 Field<vector>::operator=
153 (
154 (p0_ + amplitude_*sin(omega_*t.value()) - p.localPoints())
163{
165 os.writeEntry("amplitude", amplitude_);
166 os.writeEntry("omega", omega_);
167 p0_.writeEntry("p0", os);
168 this->writeValueEntry(os);
169}
170
171
172// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
173
175(
177 oscillatingVelocityPointPatchVectorField
178);
179
180// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
181
182} // End namespace Foam
183
184// ************************************************************************* //
Macros for easy insertion into run-time selection tables.
if(patchID !=-1)
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
const Mesh & mesh() const noexcept
Return const reference to mesh.
Generic templated field type that is much like a Foam::List except that it is expected to hold numeri...
Definition Field.H:172
void autoMap(const FieldMapper &map, const bool applyFlip=true)
Map from self.
Definition Field.C:465
void writeEntry(const word &keyword, Ostream &os) const
Write the field as a dictionary entry.
Definition Field.C:754
void rmap(const UList< Type > &mapF, const labelUList &mapAddressing)
1 to 1 reverse-map from the given field
Definition Field.C:528
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition Ostream.H:59
scalar deltaTValue() const noexcept
Return time step value.
Definition TimeStateI.H:49
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition Time.H:75
bool get(const label i) const
Definition UList.H:868
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
const Type & value() const noexcept
Return const reference to value.
A FixedValue boundary condition for pointField.
fixedValuePointPatchField(const pointPatch &, const DimensionedField< vector, pointMesh > &)
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.
oscillatingVelocityPointPatchVectorField(const pointPatch &, const DimensionedField< vector, pointMesh > &)
Construct from patch and internal field.
const pointPatch & patch() const noexcept
Return the patch.
bool updated() const noexcept
True if the boundary condition has already been updated.
Foam::pointPatchFieldMapper.
Abstract base class for point-mesh patch fields.
virtual void write(Ostream &os) const
Write.
const DimensionedField< vector, pointMesh > & internalField() const noexcept
Basic pointPatch represents a set of points from the mesh.
Definition pointPatch.H:67
Mesh consisting of general polyhedral cells.
Definition polyMesh.H:79
Lookup type of boundary radiation properties.
Definition lookup.H:60
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
volScalarField & p
dynamicFvMesh & mesh
OBJstream os(runTime.globalPath()/outputName)
Namespace for OpenFOAM.
Type & refCast(U &obj)
A dynamic_cast (for references) to Type reference.
Definition typeInfo.H:172
List< label > labelList
A List of labels.
Definition List.H:62
dimensionedScalar sin(const dimensionedScalar &ds)
pointPatchField< vector > pointPatchVectorField
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
Vector< scalar > vector
Definition vector.H:57
#define makePointPatchTypeField(PatchTypeField, typePatchTypeField)
Define a concrete pointPatchField type and add to run-time tables Example, (pointPatchScalarField,...
dictionary dict