Newer
Older
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
#include <config.h>
#include <fenv.h>
//#define LAPLACE_DEBUG
//#define HARMONIC_ENERGY_FD_GRADIENT
#define RIGIDBODYMOTION3
#include <dune/common/bitsetvector.hh>
#include <dune/common/parametertree.hh>
#include <dune/common/parametertreeparser.hh>
#include <dune/grid/uggrid.hh>
#include <dune/grid/onedgrid.hh>
#include <dune/grid/geometrygrid.hh>
#include <dune/grid/utility/structuredgridfactory.hh>
#include <dune/grid/io/file/amirameshreader.hh>
#include <dune/grid/io/file/amirameshwriter.hh>
#include <dune/solvers/solvers/iterativesolver.hh>
#include <dune/solvers/norms/energynorm.hh>
#include <dune/gfe/rotation.hh>
#include <dune/gfe/unitvector.hh>
#include <dune/gfe/realtuple.hh>
#include <dune/gfe/harmonicenergystiffness.hh>
#include <dune/gfe/geodesicfeassembler.hh>
#include <dune/gfe/riemanniantrsolver.hh>
// grid dimension
const int dim = 2;
// Image space of the geodesic fe functions
#ifdef RIGIDBODYMOTION3
typedef RigidBodyMotion<3> TargetSpace;
#endif
// Tangent vector of the image space
const int blocksize = TargetSpace::TangentVector::size;
using namespace Dune;
template <class HostGridView>
class DeformationFunction
: public Dune :: DiscreteCoordFunction< double, 3, DeformationFunction<HostGridView> >
{
typedef DeformationFunction<HostGridView> This;
typedef Dune :: DiscreteCoordFunction< double, 3, This > Base;
public:
DeformationFunction(const HostGridView& gridView,
const std::vector<RigidBodyMotion<3> >& deformedPosition)
: gridView_(gridView),
deformedPosition_(deformedPosition)
{}
void evaluate ( const typename HostGridView::template Codim<dim>::Entity& hostEntity, unsigned int corner,
FieldVector<double,3> &y ) const
{
const typename HostGridView::IndexSet& indexSet = gridView_.indexSet();
int idx = indexSet.index(hostEntity);
y = deformedPosition_[idx].r;
}
void evaluate ( const typename HostGridView::template Codim<0>::Entity& hostEntity, unsigned int corner,
FieldVector<double,3> &y ) const
{
const typename HostGridView::IndexSet& indexSet = gridView_.indexSet();
int idx = indexSet.subIndex(hostEntity, corner,dim);
y = deformedPosition_[idx].r;
}
private:
HostGridView gridView_;
const std::vector<RigidBodyMotion<3> > deformedPosition_;
};
int main (int argc, char *argv[]) try
{
//feenableexcept(FE_INVALID);
typedef std::vector<TargetSpace> SolutionType;
// parse data file
ParameterTree parameterSet;
if (argc==2)
ParameterTreeParser::readINITree(argv[1], parameterSet);
else
ParameterTreeParser::readINITree("cosserat-continuum.parset", parameterSet);
// read solver settings
const int numLevels = parameterSet.get<int>("numLevels");
const double tolerance = parameterSet.get<double>("tolerance");
const int maxTrustRegionSteps = parameterSet.get<int>("maxTrustRegionSteps");
const double initialTrustRegionRadius = parameterSet.get<double>("initialTrustRegionRadius");
const int multigridIterations = parameterSet.get<int>("numIt");
const int nu1 = parameterSet.get<int>("nu1");
const int nu2 = parameterSet.get<int>("nu2");
const int mu = parameterSet.get<int>("mu");
const int baseIterations = parameterSet.get<int>("baseIt");
const double mgTolerance = parameterSet.get<double>("mgTolerance");
const double baseTolerance = parameterSet.get<double>("baseTolerance");
const bool instrumented = parameterSet.get<bool>("instrumented");
std::string resultPath = parameterSet.get("resultPath", "");
// read problem settings
std::string path = parameterSet.get<std::string>("path");
std::string gridFile = parameterSet.get<std::string>("gridFile");
// ///////////////////////////////////////
// Create the grid
// ///////////////////////////////////////
typedef std::conditional<dim==1,OneDGrid,UGGrid<dim> >::type GridType;
array<unsigned int,dim> elements;
elements.fill(3);
shared_ptr<GridType> gridPtr = StructuredGridFactory<GridType>::createSimplexGrid(FieldVector<double,dim>(0),
FieldVector<double,dim>(1),
elements);
GridType& grid = *gridPtr.get();
grid.globalRefine(numLevels-1);
SolutionType x(grid.size(dim));
// /////////////////////////////////////////
// Read Dirichlet values
// /////////////////////////////////////////
BitSetVector<1> allNodes(grid.size(dim));
allNodes.setAll();
LeafBoundaryPatch<GridType> dirichletBoundary(grid, allNodes);
BitSetVector<blocksize> dirichletNodes(grid.size(dim));
for (int i=0; i<dirichletNodes.size(); i++)
dirichletNodes[i] = dirichletBoundary.containsVertex(i);
// //////////////////////////
// Initial solution
// //////////////////////////
FieldVector<double,3> yAxis(0);
yAxis[1] = 1;
GridType::LeafGridView::Codim<dim>::Iterator vIt = grid.leafbegin<dim>();
GridType::LeafGridView::Codim<dim>::Iterator vEndIt = grid.leafend<dim>();
for (; vIt!=vEndIt; ++vIt) {
int idx = grid.leafIndexSet().index(*vIt);
x[idx].r = 0;
for (int i=0; i<dim; i++)
x[idx].r[i] = vIt->geometry().corner(0)[i];
// x[idx].q is the identity, set by the default constructor
if (dirichletNodes[idx][0]) {
// Only the positions have Dirichlet values
x[idx].r[2] = vIt->geometry().corner(0)[0];
}
}
// backup for error measurement later
SolutionType initialIterate = x;
// ////////////////////////////////////////////////////////////
// Create an assembler for the Harmonic Energy Functional
// ////////////////////////////////////////////////////////////
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
HarmonicEnergyLocalStiffness<GridType::LeafGridView,TargetSpace> harmonicEnergyLocalStiffness;
GeodesicFEAssembler<GridType::LeafGridView,TargetSpace> assembler(grid.leafView(),
&harmonicEnergyLocalStiffness);
// /////////////////////////////////////////////////
// Create a Riemannian trust-region solver
// /////////////////////////////////////////////////
RiemannianTrustRegionSolver<GridType,TargetSpace> solver;
solver.setup(grid,
&assembler,
x,
dirichletNodes,
tolerance,
maxTrustRegionSteps,
initialTrustRegionRadius,
multigridIterations,
mgTolerance,
mu, nu1, nu2,
baseIterations,
baseTolerance,
instrumented);
// /////////////////////////////////////////////////////
// Solve!
// /////////////////////////////////////////////////////
std::cout << "Energy: " << assembler.computeEnergy(x) << std::endl;
//exit(0);
solver.setInitialSolution(x);
solver.solve();
x = solver.getSol();
// //////////////////////////////
// Output result
// //////////////////////////////
typedef GeometryGrid<GridType,DeformationFunction<GridType::LeafGridView> > DeformedGridType;
DeformationFunction<GridType::LeafGridView> deformationFunction(grid.leafView(),x);
DeformedGridType deformedGrid(grid, deformationFunction);
LeafAmiraMeshWriter<DeformedGridType>::writeSurfaceGrid(deformedGrid.leafView(), "cosseratGrid");
// Make three vector fields containing the directors
// I don't think there is a simpler way to get the data into vanilla Amira
for (int i=0; i<3; i++) {
std::vector<FieldVector<double,3> > director(x.size());
for (size_t j=0; j<x.size(); j++)
director[j] = x[j].q.director(i);
LeafAmiraMeshWriter<DeformedGridType> amiramesh;
amiramesh.addVertexData(director, deformedGrid.leafView());
std::stringstream iAsAscii;
iAsAscii << i;
amiramesh.write("cosseratOrientation"+iAsAscii.str(), true);
}
// //////////////////////////////
} catch (Exception e) {
std::cout << e << std::endl;
}