angularOscillatingVelocityPointPatchVectorField.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  axis_(dict.lookup("axis")),
49  origin_(dict.lookup("origin")),
50  angle0_(dict.lookup<scalar>("angle0")),
51  amplitude_(dict.lookup<scalar>("amplitude")),
52  omega_(dict.lookup<scalar>("omega"))
53 {
54  if (!dict.found("value"))
55  {
56  updateCoeffs();
57  }
58 
59  if (dict.found("p0"))
60  {
61  p0_ = vectorField("p0", dict , p.size());
62  }
63  else
64  {
65  p0_ = p.localPoints();
66  }
67 }
68 
69 
72 (
74  const pointPatch& p,
76  const pointPatchFieldMapper& mapper
77 )
78 :
79  fixedValuePointPatchField<vector>(ptf, p, iF, mapper),
80  axis_(ptf.axis_),
81  origin_(ptf.origin_),
82  angle0_(ptf.angle0_),
83  amplitude_(ptf.amplitude_),
84  omega_(ptf.omega_),
85  p0_(ptf.p0_)
86 {}
87 
88 
91 (
94 )
95 :
97  axis_(ptf.axis_),
98  origin_(ptf.origin_),
99  angle0_(ptf.angle0_),
100  amplitude_(ptf.amplitude_),
101  omega_(ptf.omega_),
102  p0_(ptf.p0_)
103 {}
104 
105 
106 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
107 
109 (
110  const pointPatchField<vector>& ptf,
111  const pointPatchFieldMapper& mapper
112 )
113 {
115  refCast<const angularOscillatingVelocityPointPatchVectorField>(ptf);
116 
118 
119  mapper(p0_, aOVptf.p0_);
120 }
121 
122 
124 (
125  const pointPatchField<vector>& ptf
126 )
127 {
129  refCast<const angularOscillatingVelocityPointPatchVectorField>(ptf);
130 
132 
133  p0_.reset(aOVptf.p0_);
134 }
135 
136 
138 {
139  if (this->updated())
140  {
141  return;
142  }
143 
144  const polyMesh& mesh = this->internalField().mesh()();
145  const Time& t = mesh.time();
146  const pointPatch& p = this->patch();
147 
148  scalar angle = angle0_ + amplitude_*sin(omega_*t.value());
149  vector axisHat = axis_/mag(axis_);
150  vectorField p0Rel(p0_ - origin_);
151 
153  (
154  (
155  p0_
156  + p0Rel*(cos(angle) - 1)
157  + (axisHat ^ p0Rel*sin(angle))
158  + (axisHat & p0Rel)*(1 - cos(angle))*axisHat
159  - p.localPoints()
160  )/t.deltaTValue()
161  );
162 
164 }
165 
166 
168 (
169  Ostream& os
170 ) const
171 {
173  writeEntry(os, "axis", axis_);
174  writeEntry(os, "origin", origin_);
175  writeEntry(os, "angle0", angle0_);
176  writeEntry(os, "amplitude", amplitude_);
177  writeEntry(os, "omega", omega_);
178  writeEntry(os, "p0", p0_);
179  writeEntry(os, "value", *this);
180 }
181 
182 
183 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
184 
186 (
189 );
190 
191 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
192 
193 } // End namespace Foam
194 
195 // ************************************************************************* //
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 reset(const Field< Type > &)
Reset the field values to the given field.
Definition: Field.C:432
friend Ostream & operator(Ostream &, const Field< vector > &)
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
angularOscillatingVelocityPointPatchVectorField(const pointPatch &, const DimensionedField< vector, pointMesh > &, const dictionary &)
Construct from patch, internal field and dictionary.
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.
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.
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
dimensioned< scalar > mag(const dimensioned< Type > &)
Field< vector > vectorField
Specialisation of Field<T> for vector.
dimensionedScalar cos(const dimensionedScalar &ds)
dictionary dict
volScalarField & p