Loading...
Searching...
No Matches
waveDisplacementPointPatchVectorField.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 Copyright (C) 2020 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
12 This file is part of OpenFOAM.
13
14 OpenFOAM is free software: you can redistribute it and/or modify it
15 under the terms of the GNU General Public License as published by
16 the Free Software Foundation, either version 3 of the License, or
17 (at your option) any later version.
18
19 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22 for more details.
23
24 You should have received a copy of the GNU General Public License
25 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26
27\*---------------------------------------------------------------------------*/
28
32#include "Time.H"
33#include "polyMesh.H"
34
35// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
36
39(
40 const pointPatch& p,
42)
43:
45 amplitude_(Zero),
46 omega_(0.0),
47 waveNumber_(Zero)
48{}
49
50
53(
54 const pointPatch& p,
56 const dictionary& dict
57)
58:
60 amplitude_(dict.lookup("amplitude")),
61 omega_(dict.get<scalar>("omega")),
62 waveNumber_(dict.getOrDefault<vector>("waveNumber", Zero))
63{
64 if (!dict.found("value"))
65 {
66 updateCoeffs();
67 }
68}
69
70
73(
75 const pointPatch& p,
77 const pointPatchFieldMapper& mapper
78)
79:
81 amplitude_(ptf.amplitude_),
82 omega_(ptf.omega_),
83 waveNumber_(ptf.waveNumber_)
84{}
85
86
89(
92)
93:
95 amplitude_(ptf.amplitude_),
96 omega_(ptf.omega_),
97 waveNumber_(ptf.waveNumber_)
98{}
99
100
101// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
102
104{
105 if (this->updated())
106 {
107 return;
108 }
109
110 const polyMesh& mesh = this->internalField().mesh()();
111 const Time& t = mesh.time();
112
113 const scalarField points( waveNumber_ & patch().localPoints());
114
115 Field<vector>::operator=
116 (
117 amplitude_*cos(omega_*t.value() - points)
118 );
119
121}
122
123
125{
127 os.writeEntry("amplitude", amplitude_);
128 os.writeEntry("omega", omega_);
129 os.writeEntry("waveNumber", waveNumber_);
130 this->writeValueEntry(os);
132
133
134// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
135
136namespace Foam
137{
139 (
141 waveDisplacementPointPatchVectorField
142 );
143}
144
145// ************************************************************************* //
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...
Generic templated field type that is much like a Foam::List except that it is expected to hold numeri...
Definition Field.H:172
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition Ostream.H:59
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 > &)
const pointPatch & patch() const noexcept
Return the patch.
bool updated() const noexcept
True if the boundary condition has already been updated.
Foam::pointPatchFieldMapper.
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.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
waveDisplacementPointPatchVectorField(const pointPatch &, const DimensionedField< vector, pointMesh > &)
Construct from patch and internal field.
volScalarField & p
dynamicFvMesh & mesh
OBJstream os(runTime.globalPath()/outputName)
const pointField & points
Namespace for OpenFOAM.
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
pointPatchField< vector > pointPatchVectorField
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
Vector< scalar > vector
Definition vector.H:57
dimensionedScalar cos(const dimensionedScalar &ds)
#define makePointPatchTypeField(PatchTypeField, typePatchTypeField)
Define a concrete pointPatchField type and add to run-time tables Example, (pointPatchScalarField,...
dictionary dict