optMeshMovementVolumetricBSplines.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-2019 PCOpt/NTUA
9  Copyright (C) 2013-2019 FOSS GP
10  Copyright (C) 2019 OpenCFD Ltd.
11 -------------------------------------------------------------------------------
12 License
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 
34 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
35 
36 namespace Foam
37 {
38  defineTypeNameAndDebug(optMeshMovementVolumetricBSplines, 0);
40  (
41  optMeshMovement,
42  optMeshMovementVolumetricBSplines,
43  dictionary
44  );
45 }
46 
47 
48 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
49 
51 (
52  const scalarField& correction
53 )
54 {
55  const label nControlPoints(correction.size()/3);
56  vectorField cpMovement(nControlPoints, Zero);
57 
58  for (label iCP = 0; iCP < nControlPoints; ++iCP)
59  {
60  cpMovement[iCP].x() = correction[3*iCP];
61  cpMovement[iCP].y() = correction[3*iCP + 1];
62  cpMovement[iCP].z() = correction[3*iCP + 2];
63  }
64  displMethodPtr_->boundControlField(cpMovement);
65 
66  return cpMovement;
67 }
68 
69 
70 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
71 
72 Foam::optMeshMovementVolumetricBSplines::optMeshMovementVolumetricBSplines
73 (
74  fvMesh& mesh,
75  const dictionary& dict,
76  const labelList& patchIDs
77 )
78 :
79  optMeshMovement(mesh, dict, patchIDs),
80  volBSplinesBase_
81  (
82  const_cast<volBSplinesBase&>(volBSplinesBase::New(mesh))
83  ),
84  cpsInit_(volBSplinesBase_.getNumberOfBoxes())
85 {
87 
88  forAll(boxes, boxI)
89  {
90  cpsInit_[boxI].setSize
91  (
92  boxes[boxI].getControlPoints().size()
93  );
94  cpsInit_[boxI] = boxes[boxI].getControlPoints();
95  }
96 }
97 
98 
99 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
100 
102 {
103  // Get controlPoint movement from correction
104  vectorField cpMovement = controlPointMovement(correction_);
105 
106  // Set movement of the B-Splines control points
107  displMethodPtr_->setControlField(cpMovement);
108 
109  // Move the mesh and check quality
111 }
112 
113 
115 {
117  const PtrList<NURBS3DVolume>& boxes = volBSplinesBase_.boxes();
118  forAll(boxes, boxI)
119  {
120  cpsInit_[boxI] = boxes[boxI].getControlPoints();
121  }
122 }
123 
124 
126 {
127  // Reset mesh points
129 
130  DebugInfo
131  << "optMeshMovementVolumetricBSplines:: resetting control points"
132  << endl;
133 
134  PtrList<NURBS3DVolume>& boxes = volBSplinesBase_.boxesRef();
135  forAll(boxes, boxI)
136  {
137  boxes[boxI].setControlPoints(cpsInit_[boxI]);
138  }
139 }
140 
141 
143 (
144  const scalarField& correction
145 )
146 {
147  const vectorField cpMovement(controlPointMovement(correction));
148  const scalar maxDisplacement
149  (
150  volBSplinesBase_.computeMaxBoundaryDisplacement
151  (
152  cpMovement,
153  patchIDs_
154  )
155  );
156 
157  Info<< "maxAllowedDisplacement/maxDisplacement of boundary\t"
158  << getMaxAllowedDisplacement() << "/" << maxDisplacement << endl;
159  scalar eta = getMaxAllowedDisplacement() / maxDisplacement;
160  Info<< "Setting eta value to " << eta << endl;
161 
162  return eta;
163 }
164 
165 
168 {
169  return volBSplinesBase_.getActiveDesignVariables();
170 }
171 
172 
173 // ************************************************************************* //
Foam::addToRunTimeSelectionTable
addToRunTimeSelectionTable(decompositionMethod, kahipDecomp, dictionary)
Foam::scalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Definition: primitiveFieldsFwd.H:46
Foam::optMeshMovement::displMethodPtr_
autoPtr< displacementMethod > displMethodPtr_
Definition: optMeshMovement.H:87
Foam::optMeshMovementVolumetricBSplines::computeEta
virtual scalar computeEta(const scalarField &correction)
Definition: optMeshMovementVolumetricBSplines.C:136
Foam::optMeshMovementVolumetricBSplines::resetDesignVariables
virtual void resetDesignVariables()
Definition: optMeshMovementVolumetricBSplines.C:118
Foam::volBSplinesBase
Class constructing a number of volumetric B-Splines boxes, read from dynamicMeshDict....
Definition: volBSplinesBase.H:55
Foam::Zero
static constexpr const zero Zero
Definition: zero.H:131
Foam::optMeshMovement::storeDesignVariables
virtual void storeDesignVariables()
Definition: optMeshMovement.C:175
Foam::optMeshMovementVolumetricBSplines::storeDesignVariables
virtual void storeDesignVariables()
Definition: optMeshMovementVolumetricBSplines.C:107
Foam::endl
Ostream & endl(Ostream &os)
Definition: Ostream.H:381
Foam::correction
tmp< fvMatrix< Type > > correction(const fvMatrix< Type > &)
Foam::optMeshMovementVolumetricBSplines::cpsInit_
List< vectorField > cpsInit_
Definition: optMeshMovementVolumetricBSplines.H:62
forAll
#define forAll(list, i)
Definition: stdFoam.H:349
optMeshMovementVolumetricBSplines.H
Foam::vectorField
Field< vector > vectorField
Specialisation of Field<T> for vector.
Definition: primitiveFieldsFwd.H:48
Foam::Field
Generic templated field type.
Definition: Field.H:59
Foam::Info
messageStream Info
Foam::List::setSize
void setSize(const label n)
Definition: List.H:218
Foam::PtrList
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition: List.H:55
Foam::volBSplinesBase::boxesRef
PtrList< NURBS3DVolume > & boxesRef()
Definition: volBSplinesBase.C:114
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::optMeshMovementVolumetricBSplines::volBSplinesBase_
volBSplinesBase & volBSplinesBase_
Definition: optMeshMovementVolumetricBSplines.H:59
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:119
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:81
Foam
Definition: atmBoundaryLayer.C:26
Foam::optMeshMovementVolumetricBSplines::controlPointMovement
vectorField controlPointMovement(const scalarField &correction)
Definition: optMeshMovementVolumetricBSplines.C:44
Foam::New
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
Definition: DimensionedFieldReuseFunctions.H:100
Foam::optMeshMovement::moveMesh
virtual void moveMesh()
Definition: optMeshMovement.C:121
DebugInfo
#define DebugInfo
Definition: messageStream.H:436
Foam::optMeshMovementVolumetricBSplines::getActiveDesignVariables
virtual labelList getActiveDesignVariables() const
Definition: optMeshMovementVolumetricBSplines.C:160
Foam::List
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: BitOps.H:58
Foam::optMeshMovement::correction_
scalarField correction_
Definition: optMeshMovement.H:77
Foam::optMeshMovement
Abstract base class for translating an update of the design variables into mesh movement.
Definition: optMeshMovement.H:50
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(combustionModel, 0)
Foam::optMeshMovement::resetDesignVariables
virtual void resetDesignVariables()
Definition: optMeshMovement.C:181
Foam::optMeshMovementVolumetricBSplines::moveMesh
void moveMesh()
Definition: optMeshMovementVolumetricBSplines.C:94