Newer
Older

Oliver Sander
committed
#include <config.h>

Oliver Sander
committed
#include <dune/grid/onedgrid.hh>
#include <dune/istl/io.hh>
#include <dune/ag-common/boundarypatch.hh>
#include <dune/ag-common/projectedblockgsstep.hh>
#include <dune/ag-common/iterativesolver.hh>
#include <dune/ag-common/geomestimator.hh>
#include <dune/ag-common/contactobsrestrict.hh>

Oliver Sander
committed
#include "src/rodwriter.hh"

Oliver Sander
committed
// Number of degrees of freedom:
// 3 (x, y, theta) for a planar rod
const int blocksize = 3;
using namespace Dune;
using std::string;
void setTrustRegionObstacles(double trustRegionRadius,
std::vector<BoxConstraint<double,blocksize> >& trustRegionObstacles,
const std::vector<BoxConstraint<double,blocksize> >& trueObstacles,

Oliver Sander
committed
const BitField& dirichletNodes)
{
//std::cout << "True obstacles\n" << trueObstacles << std::endl;
for (int j=0; j<trustRegionObstacles.size(); j++) {
for (int k=0; k<blocksize; k++) {
if (dirichletNodes[j*blocksize+k])
continue;
trustRegionObstacles[j].lower(k) =
(trueObstacles[j].lower(k) < -1e10)
? std::min(-trustRegionRadius, trueObstacles[j].upper(k) - trustRegionRadius)
: trueObstacles[j].lower(k);

Oliver Sander
committed
trustRegionObstacles[j].upper(k) =
(trueObstacles[j].upper(k) > 1e10)
? std::max(trustRegionRadius,trueObstacles[j].lower(k) + trustRegionRadius)
: trueObstacles[j].upper(k);

Oliver Sander
committed
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
}
}
//std::cout << "TrustRegion obstacles\n" << trustRegionObstacles << std::endl;
}
bool refineCondition(const FieldVector<double,1>& pos) {
return pos[2] > -2 && pos[2] < -0.5;
}
bool refineAll(const FieldVector<double,1>& pos) {
return true;
}
int main (int argc, char *argv[]) try
{
// Some types that I need
typedef BCRSMatrix<FieldMatrix<double, blocksize, blocksize> > MatrixType;
typedef BlockVector<FieldVector<double, blocksize> > VectorType;
// parse data file
ConfigParser parameterSet;
parameterSet.parseFile("staticrod2.parset");
// read solver settings
const int minLevel = parameterSet.get("minLevel", int(0));
const int maxLevel = parameterSet.get("maxLevel", int(0));
const int maxNewtonSteps = parameterSet.get("maxNewtonSteps", int(0));
const int numIt = parameterSet.get("numIt", int(0));
const int nu1 = parameterSet.get("nu1", int(0));
const int nu2 = parameterSet.get("nu2", int(0));
const int mu = parameterSet.get("mu", int(0));
const int baseIt = parameterSet.get("baseIt", int(0));
const double tolerance = parameterSet.get("tolerance", double(0));
const double baseTolerance = parameterSet.get("baseTolerance", double(0));
// Problem settings
const int numRodBaseElements = parameterSet.get("numRodBaseElements", int(0));
// ///////////////////////////////////////
// Create the two grids
// ///////////////////////////////////////

Oliver Sander
committed
GridType grid(numRodBaseElements, 0, 1);
std::vector<std::vector<BoxConstraint<double,3> > > trustRegionObstacles(1);
std::vector<BitField> hasObstacle(1);

Oliver Sander
committed
// ////////////////////////////////
// Create a multigrid solver
// ////////////////////////////////
// First create a gauss-seidel base solver
ProjectedBlockGSStep<MatrixType, VectorType> baseSolverStep;
EnergyNorm<MatrixType, VectorType> baseEnergyNorm(baseSolverStep);
IterativeSolver<VectorType> baseSolver(&baseSolverStep,
baseIt,
baseTolerance,
&baseEnergyNorm,
Solver::QUIET);

Oliver Sander
committed
// Make pre and postsmoothers
ProjectedBlockGSStep<MatrixType, VectorType> presmoother;
ProjectedBlockGSStep<MatrixType, VectorType> postsmoother;
MonotoneMGStep<MatrixType, VectorType> multigridStep(1);

Oliver Sander
committed
multigridStep.setMGType(mu, nu1, nu2);
multigridStep.dirichletNodes_ = &dirichletNodes[0];
multigridStep.basesolver_ = &baseSolver;
multigridStep.presmoother_ = &presmoother;
multigridStep.postsmoother_ = &postsmoother;
multigridStep.hasObstacle_ = &hasObstacle;
multigridStep.obstacles_ = &trustRegionObstacles;
multigridStep.verbosity_ = Solver::QUIET;
multigridStep.obstacleRestrictor_ = new ContactObsRestriction<VectorType>;

Oliver Sander
committed
EnergyNorm<MatrixType, VectorType> energyNorm(multigridStep);

Oliver Sander
committed
IterativeSolver<VectorType> solver(&multigridStep,

Oliver Sander
committed
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
184
double trustRegionRadius = 0.1;
VectorType rhs;
VectorType x(grid.size(0,1));
VectorType corr;
// //////////////////////////
// Initial solution
// //////////////////////////
x = 0;
x[x.size()-1][2] = 2*M_PI;
// /////////////////////////////////////////////////////////////////////
// Refinement Loop
// /////////////////////////////////////////////////////////////////////
for (int toplevel=0; toplevel<=maxLevel; toplevel++) {
std::cout << "####################################################" << std::endl;
std::cout << " Solving on level: " << toplevel << std::endl;
std::cout << "####################################################" << std::endl;
dirichletNodes.resize(toplevel+1);
for (int i=0; i<=toplevel; i++) {
dirichletNodes[i].resize( blocksize * grid.size(i,1), false );
for (int j=0; j<blocksize; j++) {
dirichletNodes[i][j] = true;
dirichletNodes[i][dirichletNodes[i].size()-1-j] = true;
}
}
// ////////////////////////////////////////////////////////////
// Create solution and rhs vectors
// ////////////////////////////////////////////////////////////
MatrixType hessianMatrix;

Oliver Sander
committed

Oliver Sander
committed
MatrixIndexSet indices(grid.size(toplevel,1), grid.size(toplevel,1));
rodAssembler.getNeighborsPerVertex(indices);
indices.exportIdx(hessianMatrix);
rhs.resize(grid.size(toplevel,1));
corr.resize(grid.size(toplevel,1));
// //////////////////////////////////////////////////////////
// Create obstacles
// //////////////////////////////////////////////////////////
hasObstacle.resize(toplevel+1);
for (int i=0; i<hasObstacle.size(); i++) {
hasObstacle[i].resize(grid.size(i, 1));
hasObstacle[i].setAll();
}
std::vector<std::vector<BoxConstraint<double,3> > > trueObstacles(toplevel+1);

Oliver Sander
committed
trustRegionObstacles.resize(toplevel+1);
for (int i=0; i<toplevel+1; i++) {
trueObstacles[i].resize(grid.size(i,1));
trustRegionObstacles[i].resize(grid.size(i,1));
}
for (int i=0; i<trueObstacles[toplevel].size(); i++) {
trueObstacles[toplevel][i].clear();
//trueObstacles[toplevel][i].val[0] = - x[i][0];
trueObstacles[toplevel][i].upper(0) = 0.1 - x[i][0];

Oliver Sander
committed
}
trustRegionObstacles.resize(toplevel+1);
for (int i=0; i<=toplevel; i++)
trustRegionObstacles[i].resize(grid.size(i, 1));
// ////////////////////////////////////
// Create the transfer operators
// ////////////////////////////////////
for (int k=0; k<multigridStep.mgTransfer_.size(); k++)
delete(multigridStep.mgTransfer_[k]);

Oliver Sander
committed
multigridStep.mgTransfer_.resize(toplevel);

Oliver Sander
committed
for (int i=0; i<multigridStep.mgTransfer_.size(); i++){

Oliver Sander
committed
TruncatedMGTransfer<VectorType>* newTransferOp = new TruncatedMGTransfer<VectorType>;
newTransferOp->setup(grid,i,i+1);
multigridStep.mgTransfer_[i] = newTransferOp;

Oliver Sander
committed
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
}
// /////////////////////////////////////////////////////
// Trust-Region Solver
// /////////////////////////////////////////////////////
for (int i=0; i<maxNewtonSteps; i++) {
std::cout << "----------------------------------------------------" << std::endl;
std::cout << " Trust-Region Step Number: " << i << std::endl;
std::cout << "----------------------------------------------------" << std::endl;
rhs = 0;
corr = 0;
rodAssembler.assembleGradient(x, rhs);
rodAssembler.assembleMatrix(x, hessianMatrix);
rhs *= -1;
// Create trust-region obstacle on grid0.maxLevel()
setTrustRegionObstacles(trustRegionRadius,
trustRegionObstacles[toplevel],
trueObstacles[toplevel],
dirichletNodes[toplevel]);
dynamic_cast<MultigridStep<MatrixType,VectorType>*>(solver.iterationStep_)->setProblem(hessianMatrix, corr, rhs, toplevel+1);

Oliver Sander
committed
solver.preprocess();
multigridStep.preprocess();

Oliver Sander
committed
// /////////////////////////////
// Solve !
// /////////////////////////////
solver.solve();
corr = multigridStep.getSol();

Oliver Sander
committed
printf("infinity norm of the correction: %g\n", corr.infinity_norm());

Oliver Sander
committed
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
std::cout << "CORRECTION IS SMALL ENOUGH" << std::endl;
break;
}
// ////////////////////////////////////////////////////
// Check whether trust-region step can be accepted
// ////////////////////////////////////////////////////
/** \todo Faster with expression templates */
VectorType newIterate = x; newIterate += corr;
/** \todo Don't always recompute oldEnergy */
double oldEnergy = rodAssembler.computeEnergy(x);
double energy = rodAssembler.computeEnergy(newIterate);
if (energy >= oldEnergy) {
printf("Richtung ist keine Abstiegsrichtung!\n");
// std::cout << "corr[0]\n" << corr[0] << std::endl;
exit(0);
}
// Add correction to the current solution
x += corr;
// Subtract correction from the current obstacle
for (int k=0; k<corr.size(); k++)
trueObstacles[grid.maxLevel()][k] -= corr[k];
}
// //////////////////////////////
// Output result
// //////////////////////////////
// Write Lagrange multiplyers
std::stringstream levelAsAscii;
levelAsAscii << toplevel;
std::string lagrangeFilename = "pressure/lagrange_" + levelAsAscii.str();
std::ofstream lagrangeFile(lagrangeFilename.c_str());
VectorType lagrangeMultipliers;
rodAssembler.assembleGradient(x, lagrangeMultipliers);
lagrangeFile << lagrangeMultipliers << std::endl;
// Write result grid
std::string solutionFilename = "solutions/rod_" + levelAsAscii.str() + ".result";
writeRod(x, solutionFilename);
// ////////////////////////////////////////////////////////////////////////////
// Refine locally and transfer the current solution to the new leaf level
// ////////////////////////////////////////////////////////////////////////////
GeometricEstimator<GridType> estimator;
estimator.estimate(grid, (toplevel<=minLevel) ? refineAll : refineCondition);
P1FunctionManager<GridType,double> functionManager(grid);
LeafP1Function<GridType,double,blocksize> sol(grid);
*sol = x;
grid.preAdapt();
sol.preAdapt();
grid.adapt();
sol.postAdapt(functionManager);
grid.postAdapt();
x = *sol;

Oliver Sander
committed
//writeRod(x, "solutions/rod_1.result");
}
} catch (Exception e) {
std::cout << e << std::endl;
}