volumetricBSplinesMotionSolver.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 | www.openfoam.com
6 \\/ M anipulation |
7-------------------------------------------------------------------------------
8 Copyright (C) 2007-2020 PCOpt/NTUA
9 Copyright (C) 2013-2020 FOSS GP
10 Copyright (C) 2019 OpenCFD Ltd.
11-------------------------------------------------------------------------------
12License
13 This file is part of OpenFOAM.
14
15 OpenFOAM is free software: you can redistribute it and/or modify it
16 under the terms of the GNU General Public License as published by
17 the Free Software Foundation, either version 3 of the License, or
18 (at your option) any later version.
19
20 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
21 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
22 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
23 for more details.
24
25 You should have received a copy of the GNU General Public License
26 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
27
28\*---------------------------------------------------------------------------*/
29
32
33// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
34
35namespace Foam
36{
38
40 (
44 );
45}
46
47
48// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
49
51(
52 const polyMesh& mesh,
53 const IOdictionary& dict
54)
55:
56 motionSolver(mesh, dict, typeName),
57 volBSplinesBase_
58 (
59 const_cast<volBSplinesBase&>
60 (
62 )
63 ),
64 controlPointsMovement_
65 (
66 volBSplinesBase_.getTotalControlPointsNumber(),
67 Zero
68 )
69{}
70
71
72// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
73
76{
77 tmp<vectorField> tPointMovement(new vectorField(mesh().points()));
78 vectorField& pointMovement = tPointMovement.ref();
79
80 label pastControlPoints(0);
81 PtrList<NURBS3DVolume>& boxes = volBSplinesBase_.boxesRef();
82 forAll(boxes, iNURB)
83 {
84 const label nb = boxes[iNURB].getControlPoints().size();
85 vectorField localControlPointsMovement(nb, Zero);
86
87 forAll(localControlPointsMovement, iCP)
88 {
89 localControlPointsMovement[iCP] =
90 controlPointsMovement_[pastControlPoints + iCP];
91 }
92
94 partialMovement
95 (
96 boxes[iNURB].computeNewPoints
97 (
98 localControlPointsMovement
99 )
100 );
101
102 pointMovement += partialMovement() - mesh().points();
103
104 pastControlPoints += nb;
105 }
106
107 return tPointMovement;
108}
109
110
112{
113 // Do nothing
114}
115
116
118{
119 // Do nothing
120}
121
122
124{
125 // Do nothing
126}
127
128
130(
131 const vectorField& controlPointsMovement
132)
133{
134 if (controlPointsMovement_.size() != controlPointsMovement.size())
135 {
137 << "Attempting to replace controlPointsMovement with a set of "
138 << "different size"
139 << exit(FatalError);
140 }
141 controlPointsMovement_ = controlPointsMovement;
142}
143
144
146(
147 vectorField& controlPointsMovement
148)
149{
150 volBSplinesBase_.boundControlPointMovement(controlPointsMovement);
151}
152
153
154// ************************************************************************* //
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:57
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition: PtrList.H:73
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:114
label size() const noexcept
The number of elements in the list.
Definition: UPtrListI.H:106
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:91
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
Definition: mapPolyMesh.H:162
void movePoints()
Update for new mesh geometry.
void updateMesh()
Update for new mesh topology.
Virtual base class for mesh motion solver.
Definition: motionSolver.H:61
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:81
virtual const pointField & points() const
Return raw points.
Definition: polyMesh.C:1083
A class for managing temporary objects.
Definition: tmp.H:65
T & ref() const
Definition: tmpI.H:227
Class constructing a number of volumetric B-Splines boxes, read from dynamicMeshDict....
A mesh motion solver based on volumetric B-Splines.
virtual tmp< pointField > curPoints() const
Return point location obtained from the current motion field.
void boundControlPointMovement(vectorField &controlPointsMovement)
Bound control points movement.
virtual void solve()
Solve for motion/ Does nothing.
void setControlPointsMovement(const vectorField &controlPointsMovement)
Set control points movement.
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition: className.H:121
dynamicFvMesh & mesh
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
const pointField & points
Namespace for OpenFOAM.
To & refCast(From &r)
Reference type cast template function.
Definition: typeInfo.H:131
Field< vector > vectorField
Specialisation of Field<T> for vector.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh > > &tdf1, const word &name, const dimensionSet &dimensions)
Global function forwards to reuseTmpDimensionedField::New.
error FatalError
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333