oscillatingVelocityPointPatchVectorField.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  if (dict.found("p0"))
57  {
58  p0_ = vectorField("p0", dict , p.size());
59  }
60  else
61  {
62  p0_ = p.localPoints();
63  }
64 }
65 
66 
69 (
71  const pointPatch& p,
73  const pointPatchFieldMapper& mapper
74 )
75 :
76  fixedValuePointPatchField<vector>(ptf, p, iF, mapper),
77  amplitude_(ptf.amplitude_),
78  omega_(ptf.omega_),
79  p0_(mapper(ptf.p0_))
80 {}
81 
82 
85 (
88 )
89 :
91  amplitude_(ptf.amplitude_),
92  omega_(ptf.omega_),
93  p0_(ptf.p0_)
94 {}
95 
96 
97 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
98 
100 (
101  const pointPatchField<vector>& ptf,
102  const pointPatchFieldMapper& mapper
103 )
104 {
106  refCast<const oscillatingVelocityPointPatchVectorField>(ptf);
107 
109 
110  mapper(p0_, oVptf.p0_);
111 }
112 
113 
115 (
116  const pointPatchField<vector>& ptf
117 )
118 {
120  refCast<const oscillatingVelocityPointPatchVectorField>(ptf);
121 
123 
124  p0_.reset(oVptf.p0_);
125 }
126 
127 
129 {
130  if (this->updated())
131  {
132  return;
133  }
134 
135  const polyMesh& mesh = this->internalField().mesh()();
136  const Time& t = mesh.time();
137  const pointPatch& p = this->patch();
138 
140  (
141  (p0_ + amplitude_*sin(omega_*t.value()) - p.localPoints())
142  /t.deltaTValue()
143  );
144 
146 }
147 
148 
150 {
152  writeEntry(os, "amplitude", amplitude_);
153  writeEntry(os, "omega", omega_);
154  writeEntry(os, "p0", p0_);
155  writeEntry(os, "value", *this);
156 }
157 
158 
159 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
160 
162 (
165 );
166 
167 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
168 
169 } // End namespace Foam
170 
171 // ************************************************************************* //
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...
Pre-declare SubField and related Field type.
Definition: Field.H:82
void reset(const Field< Type > &)
Reset the field values to the given field.
Definition: Field.C:432
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:57
scalar deltaTValue() const
Return time step value.
Definition: TimeStateI.H:34
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 reset(const pointPatchField< vector > &)
Reset the pointPatchField to the given pointPatchField.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
virtual void map(const pointPatchField< vector > &, const pointPatchFieldMapper &)
Map the given pointPatchField onto this pointPatchField.
oscillatingVelocityPointPatchVectorField(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.
const pointPatch & patch() const
Return patch.
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.
virtual void reset(const pointPatchField< Type > &)
Reset the pointPatchField to the given pointPatchField.
virtual void map(const pointPatchField< Type > &, const pointPatchFieldMapper &)
Map the given pointPatchField onto this pointPatchField.
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
Field< vector > vectorField
Specialisation of Field<T> for vector.
dictionary dict
volScalarField & p