uncoupledSixDoFRigidBodyDisplacementPointPatchVectorField.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 | Copyright (C) 2011-2016 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 "fvMesh.H"
31 #include "volFields.H"
33 
34 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
35 
36 namespace Foam
37 {
38 
39 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
40 
43 (
44  const pointPatch& p,
46 )
47 :
49  motion_(),
50  initialPoints_(p.localPoints()),
51  curTimeIndex_(-1)
52 {}
53 
54 
57 (
58  const pointPatch& p,
60  const dictionary& dict
61 )
62 :
64  motion_(dict, dict),
65  curTimeIndex_(-1)
66 {
67  if (!dict.found("value"))
68  {
69  updateCoeffs();
70  }
71 
72  if (dict.found("initialPoints"))
73  {
74  initialPoints_ = vectorField("initialPoints", dict , p.size());
75  }
76  else
77  {
78  initialPoints_ = p.localPoints();
79  }
80 }
81 
82 
85 (
87  const pointPatch& p,
89  const pointPatchFieldMapper& mapper
90 )
91 :
92  fixedValuePointPatchField<vector>(ptf, p, iF, mapper),
93  motion_(ptf.motion_),
94  initialPoints_(ptf.initialPoints_, mapper),
95  curTimeIndex_(-1)
96 {}
97 
98 
101 (
104 )
105 :
107  motion_(ptf.motion_),
108  initialPoints_(ptf.initialPoints_),
109  curTimeIndex_(-1)
110 {}
111 
112 
113 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
114 
116 (
117  const pointPatchFieldMapper& m
118 )
119 {
121 
122  initialPoints_.autoMap(m);
123 }
124 
125 
127 (
128  const pointPatchField<vector>& ptf,
129  const labelList& addr
130 )
131 {
133  refCast
134  <
136  >(ptf);
137 
139 
140  initialPoints_.rmap(uSDoFptf.initialPoints_, addr);
141 }
142 
143 
145 {
146  if (this->updated())
147  {
148  return;
149  }
150 
151  const polyMesh& mesh = this->internalField().mesh()();
152  const Time& t = mesh.time();
153 
154  // Store the motion state at the beginning of the time-step
155  bool firstIter = false;
156  if (curTimeIndex_ != t.timeIndex())
157  {
158  motion_.newTime();
159  curTimeIndex_ = t.timeIndex();
160  firstIter = true;
161  }
162 
163  vector gravity = Zero;
164 
165  if (db().foundObject<uniformDimensionedVectorField>("g"))
166  {
169 
170  gravity = g.value();
171  }
172 
173  // Do not modify the accelerations, except with gravity, where available
174  motion_.update
175  (
176  firstIter,
177  gravity*motion_.mass(),
178  Zero,
179  t.deltaTValue(),
180  t.deltaT0Value()
181  );
182 
184  (
185  motion_.transform(initialPoints_) - initialPoints_
186  );
187 
189 }
190 
191 
193 (
194  Ostream& os
195 ) const
196 {
198  motion_.write(os);
199  initialPoints_.writeEntry("initialPoints", os);
200  writeEntry("value", os);
201 }
202 
203 
204 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
205 
207 (
210 );
211 
212 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
213 
214 } // End namespace Foam
215 
216 // ************************************************************************* //
label timeIndex() const
Return current time index.
Definition: TimeStateI.H:35
const Time & time() const
Return time.
dictionary dict
bool updated() const
Return true if the boundary condition has already been updated.
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
To & refCast(From &r)
Reference type cast template function.
Definition: typeInfo.H:106
scalar mass() const
Return the mass.
scalar deltaT0Value() const
Return old time step value.
Definition: TimeStateI.H:47
Foam::pointPatchFieldMapper.
void writeEntry(const word &keyword, Ostream &os) const
Write the field as a dictionary entry.
Definition: Field.C:719
const Type & value() const
Return const reference to value.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:68
Macros for easy insertion into run-time selection tables.
point transform(const point &initialPoints) const
Transform the given initial state point by the current motion.
const objectRegistry & db() const
Return local objectRegistry.
const Type & lookupObject(const word &name) const
Lookup and return the object of the given Type.
dynamicFvMesh & mesh
Pre-declare SubField and related Field type.
Definition: Field.H:57
bool found(const word &, bool recursive=false, bool patternMatch=true) const
Search dictionary for given keyword.
Definition: dictionary.C:306
static const zero Zero
Definition: zero.H:91
virtual void autoMap(const pointPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:53
const dimensionedVector & g
void newTime()
Store the motion state at the beginning of the time-step.
scalar deltaTValue() const
Return time step value.
Definition: TimeStateI.H:41
Basic pointPatch represents a set of points from the mesh.
Definition: pointPatch.H:56
void write(Ostream &) const
Write.
virtual void autoMap(const pointPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
makePointPatchTypeField(pointPatchVectorField, solidBodyMotionDisplacementPointPatchVectorField)
virtual void rmap(const pointPatchField< vector > &, const labelList &)
Reverse map the given pointPatchField onto this pointPatchField.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Field< vector > vectorField
Specialisation of Field<T> for vector.
const DimensionedField< vector, pointMesh > & internalField() const
Return dimensioned internal field reference.
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
virtual void write(Ostream &) const
Write.
volScalarField & p
uncoupledSixDoFRigidBodyDisplacementPointPatchVectorField(const pointPatch &, const DimensionedField< vector, pointMesh > &)
Construct from patch and internal field.
virtual void rmap(const pointPatchField< Type > &, const labelList &)
Reverse map the given PointPatchField onto.
void update(bool firstIter, const vector &fGlobal, const vector &tauGlobal, scalar deltaT, scalar deltaT0)
Symplectic integration of velocities, orientation and position.
Namespace for OpenFOAM.
label size() const
Return size.