sixDoFRigidBodyDisplacementPointPatchVectorField.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 #include "forces.H"
34 
35 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
36 
37 namespace Foam
38 {
39 
40 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
41 
44 (
45  const pointPatch& p,
47 )
48 :
50  motion_(),
51  initialPoints_(p.localPoints()),
52  rhoInf_(1.0),
53  rhoName_("rho"),
54  lookupGravity_(-1),
55  g_(Zero),
56  curTimeIndex_(-1)
57 {}
58 
59 
62 (
63  const pointPatch& p,
65  const dictionary& dict
66 )
67 :
69  motion_(dict, dict),
70  rhoInf_(1.0),
71  rhoName_(dict.lookupOrDefault<word>("rho", "rho")),
72  lookupGravity_(-1),
73  g_(Zero),
74  curTimeIndex_(-1)
75 {
76  if (rhoName_ == "rhoInf")
77  {
78  rhoInf_ = readScalar(dict.lookup("rhoInf"));
79  }
80 
81  if (dict.readIfPresent("g", g_))
82  {
83  lookupGravity_ = -2;
84  }
85 
86  if (!dict.found("value"))
87  {
88  updateCoeffs();
89  }
90 
91  if (dict.found("initialPoints"))
92  {
93  initialPoints_ = vectorField("initialPoints", dict , p.size());
94  }
95  else
96  {
97  initialPoints_ = p.localPoints();
98  }
99 }
100 
101 
104 (
106  const pointPatch& p,
108  const pointPatchFieldMapper& mapper
109 )
110 :
111  fixedValuePointPatchField<vector>(ptf, p, iF, mapper),
112  motion_(ptf.motion_),
113  initialPoints_(ptf.initialPoints_, mapper),
114  rhoInf_(ptf.rhoInf_),
115  rhoName_(ptf.rhoName_),
116  lookupGravity_(ptf.lookupGravity_),
117  g_(ptf.g_),
118  curTimeIndex_(-1)
119 {}
120 
121 
124 (
127 )
128 :
130  motion_(ptf.motion_),
131  initialPoints_(ptf.initialPoints_),
132  rhoInf_(ptf.rhoInf_),
133  rhoName_(ptf.rhoName_),
134  lookupGravity_(ptf.lookupGravity_),
135  g_(ptf.g_),
136  curTimeIndex_(-1)
137 {}
138 
139 
140 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
141 
143 (
144  const pointPatchFieldMapper& m
145 )
146 {
148 
149  initialPoints_.autoMap(m);
150 }
151 
152 
154 (
155  const pointPatchField<vector>& ptf,
156  const labelList& addr
157 )
158 {
160  refCast<const sixDoFRigidBodyDisplacementPointPatchVectorField>(ptf);
161 
163 
164  initialPoints_.rmap(sDoFptf.initialPoints_, addr);
165 }
166 
167 
169 {
170  if (this->updated())
171  {
172  return;
173  }
174 
175  if (lookupGravity_ < 0)
176  {
177  if (db().foundObject<uniformDimensionedVectorField>("g"))
178  {
179  if (lookupGravity_ == -2)
180  {
182  << "Specifying the value of g in this boundary condition "
183  << "when g is available from the database is considered "
184  << "a fatal error to avoid the possibility of inconsistency"
185  << exit(FatalError);
186  }
187  else
188  {
189  lookupGravity_ = 1;
190  }
191  }
192  else
193  {
194  lookupGravity_ = 0;
195  }
196  }
197 
198  const polyMesh& mesh = this->internalField().mesh()();
199  const Time& t = mesh.time();
200  const pointPatch& ptPatch = this->patch();
201 
202  // Store the motion state at the beginning of the time-step
203  bool firstIter = false;
204  if (curTimeIndex_ != t.timeIndex())
205  {
206  motion_.newTime();
207  curTimeIndex_ = t.timeIndex();
208  firstIter = true;
209  }
210 
211  dictionary forcesDict;
212 
213  forcesDict.add("type", functionObjects::forces::typeName);
214  forcesDict.add("patches", wordList(1, ptPatch.name()));
215  forcesDict.add("rhoInf", rhoInf_);
216  forcesDict.add("rho", rhoName_);
217  forcesDict.add("CofR", motion_.centreOfRotation());
218 
219  functionObjects::forces f("forces", db(), forcesDict);
220 
221  f.calcForcesMoment();
222 
223  // Get the forces on the patch faces at the current positions
224 
225  if (lookupGravity_ == 1)
226  {
229 
230  g_ = g.value();
231  }
232 
233  // scalar ramp = min(max((t.value() - 5)/10, 0), 1);
234  scalar ramp = 1.0;
235 
236  motion_.update
237  (
238  firstIter,
239  ramp*(f.forceEff() + motion_.mass()*g_),
240  ramp*(f.momentEff() + motion_.mass()*(motion_.momentArm() ^ g_)),
241  t.deltaTValue(),
242  t.deltaT0Value()
243  );
244 
246  (
247  motion_.transform(initialPoints_) - initialPoints_
248  );
249 
251 }
252 
253 
255 {
257 
258  os.writeKeyword("rho") << rhoName_ << token::END_STATEMENT << nl;
259 
260  if (rhoName_ == "rhoInf")
261  {
262  os.writeKeyword("rhoInf") << rhoInf_ << token::END_STATEMENT << nl;
263  }
264 
265  if (lookupGravity_ == 0 || lookupGravity_ == -2)
266  {
267  os.writeKeyword("g") << g_ << token::END_STATEMENT << nl;
268  }
269 
270  motion_.write(os);
271 
272  initialPoints_.writeEntry("initialPoints", os);
273 
274  writeEntry("value", os);
275 }
276 
277 
278 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
279 
281 (
284 );
285 
286 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
287 
288 } // End namespace Foam
289 
290 // ************************************************************************* //
label timeIndex() const
Return current time index.
Definition: TimeStateI.H:35
const Time & time() const
Return time.
dictionary dict
virtual void autoMap(const pointPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
bool updated() const
Return true if the boundary condition has already been updated.
error FatalError
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:319
sixDoFRigidBodyDisplacementPointPatchVectorField(const pointPatch &, const DimensionedField< vector, pointMesh > &)
Construct from patch and internal field.
virtual vector forceEff() const
Return the total force.
Definition: forces.C:887
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.
bool add(entry *, bool mergeEntry=false)
Add a new entry.
Definition: dictionary.C:737
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
A class for handling words, derived from string.
Definition: word.H:59
virtual void rmap(const pointPatchField< vector > &, const labelList &)
Reverse map the given pointPatchField onto this pointPatchField.
bool readIfPresent(const word &, T &, bool recursive=false, bool patternMatch=true) const
Find an entry if present, and assign to T.
bool found(const word &, bool recursive=false, bool patternMatch=true) const
Search dictionary for given keyword.
Definition: dictionary.C:306
virtual void calcForcesMoment()
Calculate the forces and moments.
Definition: forces.C:747
static const zero Zero
Definition: zero.H:91
virtual const word & name() const =0
Return name.
bool readScalar(const char *buf, doubleScalar &s)
Read whole of buf as a scalar. Return true if succesful.
Definition: doubleScalar.H:63
virtual void autoMap(const pointPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
const point & centreOfRotation() const
Return the current centre of rotation.
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:53
static const char nl
Definition: Ostream.H:262
const dimensionedVector & g
Ostream & writeKeyword(const keyType &)
Write the keyword followed by an appropriate indentation.
Definition: Ostream.C:54
labelList f(nPoints)
void newTime()
Store the motion state at the beginning of the time-step.
scalar deltaTValue() const
Return time step value.
Definition: TimeStateI.H:41
virtual vector momentEff() const
Return the total moment.
Definition: forces.C:893
const pointPatch & patch() const
Return patch.
List< word > wordList
A List of words.
Definition: fileName.H:54
Basic pointPatch represents a set of points from the mesh.
Definition: pointPatch.H:56
void write(Ostream &) const
Write.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
makePointPatchTypeField(pointPatchVectorField, solidBodyMotionDisplacementPointPatchVectorField)
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.
This function object calculates the forces and moments by integrating the pressure and skin-friction ...
Definition: forces.H:196
volScalarField & p
T lookupOrDefault(const word &, const T &, bool recursive=false, bool patternMatch=true) const
Find and return a T,.
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.
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:451