solidBodyMotionDisplacementPointPatchVectorField.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) 2013 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 "transformField.H"
29 #include "pointPatchFields.H"
30 
31 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
32 
33 namespace Foam
34 {
35 
36 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * * //
37 
40 (
41  const pointPatch& p,
43 )
44 :
45  fixedValuePointPatchVectorField(p, iF),
46  SBMFPtr_()
47 {}
48 
49 
52 (
53  const pointPatch& p,
55  const dictionary& dict
56 )
57 :
58  fixedValuePointPatchVectorField(p, iF, dict, false),
59  SBMFPtr_(solidBodyMotionFunction::New(dict, this->db().time()))
60 {
61  if (!dict.found("value"))
62  {
63  // Determine current local points and offset
64  fixedValuePointPatchVectorField::operator==
65  (
66  transform(SBMFPtr_().transformation(), localPoints0())
67  -localPoints0()
68  );
69  }
70 }
71 
72 
75 (
77  const pointPatch& p,
79  const pointPatchFieldMapper& mapper
80 )
81 :
82  fixedValuePointPatchVectorField(ptf, p, iF, mapper),
83  SBMFPtr_(ptf.SBMFPtr_().clone().ptr())
84 {
85  // For safety re-evaluate
86 
87  fixedValuePointPatchVectorField::operator==
88  (
89  transform(SBMFPtr_().transformation(), localPoints0())
90  -localPoints0()
91  );
92 }
93 
94 
97 (
99 )
100 :
101  fixedValuePointPatchVectorField(ptf),
102  SBMFPtr_(ptf.SBMFPtr_().clone().ptr())
103 {}
104 
105 
108 (
111 )
112 :
113  fixedValuePointPatchVectorField(ptf, iF),
114  SBMFPtr_(ptf.SBMFPtr_().clone().ptr())
115 {
116  // For safety re-evaluate
117 
118  fixedValuePointPatchVectorField::operator==
119  (
120  transform(SBMFPtr_().transformation(), localPoints0())
121  -localPoints0()
122  );
123 }
124 
125 
126 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
127 
128 const pointField&
130 {
131  if (!localPoints0Ptr_.valid())
132  {
133  pointIOField points0
134  (
135  IOobject
136  (
137  "points",
138  this->db().time().constant(),
140  this->db(),
143  false
144  )
145  );
146 
147  localPoints0Ptr_.reset(new pointField(points0, patch().meshPoints()));
148  }
149  return localPoints0Ptr_();
150 }
151 
152 
154 {
155  if (this->updated())
156  {
157  return;
158  }
159 
160  // Determine current local points and offset
161  fixedValuePointPatchVectorField::operator==
162  (
163  transform(SBMFPtr_().transformation(), localPoints0())
164  -localPoints0()
165  );
166 
167  fixedValuePointPatchVectorField::updateCoeffs();
168 }
169 
170 
172 write(Ostream& os) const
173 {
174  // Note: write value
176 
177  os.writeKeyword(solidBodyMotionFunction::typeName) << SBMFPtr_->type()
178  << token::END_STATEMENT << nl;
179  os << indent << word(SBMFPtr_->type() + "Coeffs");
180  SBMFPtr_->writeData(os);
181 }
182 
183 
184 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
185 
187 (
190 );
191 
192 
193 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
194 
195 } // End namespace Foam
196 
197 // ************************************************************************* //
static word meshSubDir
Return the mesh sub-directory name (usually "polyMesh")
Definition: polyMesh.H:309
bool found(const word &, bool recursive=false, bool patternMatch=true) const
Search dictionary for given keyword.
Definition: dictionary.C:306
Basic pointPatch represents a set of points from the mesh.
Definition: pointPatch.H:56
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
static autoPtr< solidBodyMotionFunction > New(const dictionary &SBMFCoeffs, const Time &runTime)
Select constructed from the SBMFCoeffs dictionary and Time.
Ostream & indent(Ostream &os)
Indent stream.
Definition: Ostream.H:221
A class for handling words, derived from string.
Definition: word.H:59
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:42
Namespace for OpenFOAM.
runTime write()
Enables the specification of a fixed value boundary condition using the solid body motion functions...
makePointPatchTypeField(pointPatchVectorField, solidBodyMotionDisplacementPointPatchVectorField)
static const char nl
Definition: Ostream.H:260
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:91
dimensionSet transform(const dimensionSet &)
Definition: dimensionSet.C:465
Spatial transformation functions for primitive fields.
Macros for easy insertion into run-time selection tables.
Ostream & writeKeyword(const keyType &)
Write the keyword followed by an appropriate indentation.
Definition: Ostream.C:59
Foam::pointPatchFieldMapper.
Constant dispersed-phase particle diameter model.
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:53
solidBodyMotionDisplacementPointPatchVectorField(const pointPatch &, const DimensionedField< vector, pointMesh > &)
Construct from patch and internal field.