oscillatingDisplacementPointPatchVectorField.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2023 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
27 #include "pointPatchFields.H"
29 #include "Time.H"
30 #include "polyMesh.H"
31 
32 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
33 
34 namespace Foam
35 {
36 
37 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
38 
41 (
42  const pointPatch& p,
44  const dictionary& dict
45 )
46 :
48  amplitude_(dict.lookup("amplitude")),
49  omega_(dict.lookup<scalar>("omega"))
50 {
51  if (!dict.found("value"))
52  {
53  updateCoeffs();
54  }
55 }
56 
57 
60 (
62  const pointPatch& p,
64  const pointPatchFieldMapper& mapper
65 )
66 :
67  fixedValuePointPatchField<vector>(ptf, p, iF, mapper),
68  amplitude_(ptf.amplitude_),
69  omega_(ptf.omega_)
70 {}
71 
72 
75 (
78 )
79 :
81  amplitude_(ptf.amplitude_),
82  omega_(ptf.omega_)
83 {}
84 
85 
86 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
87 
89 {
90  if (this->updated())
91  {
92  return;
93  }
94 
95  const polyMesh& mesh = this->internalField().mesh()();
96  const Time& t = mesh.time();
97 
98  Field<vector>::operator=(amplitude_*sin(omega_*t.value()));
99 
101 }
102 
103 
105 {
107  writeEntry(os, "amplitude", amplitude_);
108  writeEntry(os, "omega", omega_);
109  writeEntry(os, "value", *this);
110 }
111 
112 
113 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
114 
116 (
119 );
120 
121 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
122 
123 } // End namespace Foam
124 
125 // ************************************************************************* //
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 operator=(const Field< Type > &)
Definition: Field.C:526
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:57
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:76
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:160
const Type & value() const
Return const reference to value.
A FixedValue boundary condition for pointField.
const Time & time() const
Return time.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
oscillatingDisplacementPointPatchVectorField(const pointPatch &, const DimensionedField< vector, pointMesh > &, const dictionary &)
Construct from patch, internal field and dictionary.
Foam::pointPatchFieldMapper.
Abstract base class for point-mesh patch fields.
virtual void write(Ostream &) const
Write.
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:57
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:80
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Namespace for OpenFOAM.
makePointPatchTypeField(pointPatchVectorField, solidBodyMotionDisplacementPointPatchVectorField)
dimensionedScalar sin(const dimensionedScalar &ds)
void writeEntry(Ostream &os, const HashTable< T, Key, Hash > &ht)
Definition: HashTableIO.C:96
dictionary dict
volScalarField & p