oscillatingDisplacementPointPatchVectorField.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
41oscillatingDisplacementPointPatchVectorField::
42oscillatingDisplacementPointPatchVectorField
43(
44 const pointPatch& p,
46)
47:
49 amplitude_(Zero),
50 omega_(0.0)
51{}
52
53
54oscillatingDisplacementPointPatchVectorField::
55oscillatingDisplacementPointPatchVectorField
56(
57 const pointPatch& p,
59 const dictionary& dict
60)
61:
63 amplitude_(dict.lookup("amplitude")),
64 omega_(dict.get<scalar>("omega"))
65{
66 if (!dict.found("value"))
67 {
69 }
70}
71
72
73oscillatingDisplacementPointPatchVectorField::
74oscillatingDisplacementPointPatchVectorField
75(
77 const pointPatch& p,
79 const pointPatchFieldMapper& mapper
80)
81:
82 fixedValuePointPatchField<vector>(ptf, p, iF, mapper),
83 amplitude_(ptf.amplitude_),
84 omega_(ptf.omega_)
85{}
86
87
88oscillatingDisplacementPointPatchVectorField::
89oscillatingDisplacementPointPatchVectorField
90(
93)
94:
96 amplitude_(ptf.amplitude_),
97 omega_(ptf.omega_)
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 Field<vector>::operator=(amplitude_*sin(omega_*t.value()));
114
116}
117
118
120{
122 os.writeEntry("amplitude", amplitude_);
123 os.writeEntry("omega", omega_);
124 writeEntry("value", os);
125}
126
127
128// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
129
131(
134);
135
136// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
137
138} // End namespace Foam
139
140// ************************************************************************* //
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.
void operator=(const Field< Type > &)
Copy assignment.
Definition: Field.C:641
void writeEntry(const word &keyword, Ostream &os) const
Write the field as a dictionary entry.
Definition: Field.C:614
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
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:80
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
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 updateCoeffs()
Update the coefficients associated with the patch field.
Foam::pointPatchFieldMapper.
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)
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
#define makePointPatchTypeField(PatchTypeField, typePatchTypeField)
dictionary dict