From b9de9145100a0d1c093940685ad68f3b3ea2d8a3 Mon Sep 17 00:00:00 2001 From: Elias Pipping <elias.pipping@fu-berlin.de> Date: Fri, 23 Dec 2011 16:18:54 +0100 Subject: [PATCH] Clean up construction of TNNMG solver Mostly heap -> stack --- src/one-body-sample.cc | 88 +++++++++++++++++------------------------- 1 file changed, 35 insertions(+), 53 deletions(-) diff --git a/src/one-body-sample.cc b/src/one-body-sample.cc index 8fdb314f..c5347ba0 100644 --- a/src/one-body-sample.cc +++ b/src/one-body-sample.cc @@ -268,51 +268,42 @@ int main(int argc, char *argv[]) { myGlobalNonlinearity, nodalIntegrals); // {{{ Set up TNNMG solver - auto linearBaseSolverStep = - new TruncatedBlockGSStep<OperatorType, VectorType>; - auto const baseEnergyNorm = - new EnergyNorm<OperatorType, VectorType>(*linearBaseSolverStep); - auto linearBaseSolver = new LoopSolver<VectorType>( - linearBaseSolverStep, solver_maxIterations, solver_tolerance, - baseEnergyNorm, Solver::QUIET); - - // {{{ Smoothers - auto linearPresmoother = new TruncatedBlockGSStep<OperatorType, VectorType>; - auto linearPostsmoother = - new TruncatedBlockGSStep<OperatorType, VectorType>; - // }}} - - auto linearIterationStep = new MultigridStep<OperatorType, VectorType>; - linearIterationStep->setNumberOfLevels(levels); - linearIterationStep->setMGType(parset.get<int>("solver.tnnmg.linear.mu"), - parset.get<int>("solver.tnnmg.linear.nu1"), - parset.get<int>("solver.tnnmg.linear.nu2")); - linearIterationStep->basesolver_ = linearBaseSolver; - linearIterationStep->setSmoother(linearPresmoother, linearPostsmoother); - - // {{{ Transfer operators - std::vector<CompressedMultigridTransfer<VectorType> *> transferOperators; - transferOperators.resize(refinements); + // linear iteration step components + TruncatedBlockGSStep<OperatorType, VectorType> linearBaseSolverStep; + EnergyNorm<OperatorType, VectorType> baseEnergyNorm(linearBaseSolverStep); + LoopSolver<VectorType> linearBaseSolver( + &linearBaseSolverStep, solver_maxIterations, solver_tolerance, + &baseEnergyNorm, Solver::QUIET); + TruncatedBlockGSStep<OperatorType, VectorType> linearPresmoother; + TruncatedBlockGSStep<OperatorType, VectorType> linearPostsmoother; + + // linear iteration step + MultigridStep<OperatorType, VectorType> linearIterationStep; + linearIterationStep.setNumberOfLevels(levels); + linearIterationStep.setMGType(parset.get<int>("solver.tnnmg.linear.mu"), + parset.get<int>("solver.tnnmg.linear.nu1"), + parset.get<int>("solver.tnnmg.linear.nu2")); + linearIterationStep.basesolver_ = &linearBaseSolver; + linearIterationStep.setSmoother(&linearPresmoother, &linearPostsmoother); + + // transfer operators + std::vector<CompressedMultigridTransfer<VectorType> *> transferOperators( + refinements); for (auto &x : transferOperators) x = new CompressedMultigridTransfer<VectorType>; + TransferOperatorAssembler<GridType>(grid) + .assembleOperatorPointerHierarchy(transferOperators); + linearIterationStep.setTransferOperators(transferOperators); - TransferOperatorAssembler<GridType> transferOperatorAssembler(grid); - transferOperatorAssembler.assembleOperatorPointerHierarchy( - transferOperators); - linearIterationStep->setTransferOperators(transferOperators); - // }}} - + // tnnmg iteration step typedef GenericNonlinearGS<MyBlockProblemType> NonlinearSmootherType; - typedef TruncatedNonsmoothNewtonMultigrid< - MyBlockProblemType, NonlinearSmootherType> TNNMGStepType; - - auto nonlinearSmoother = new NonlinearSmootherType; - auto multigridStep = - new TNNMGStepType(*linearIterationStep, *nonlinearSmoother); - multigridStep->setSmoothingSteps(parset.get<int>("solver.tnnmg.main.nu1"), - parset.get<int>("solver.tnnmg.main.mu"), - parset.get<int>("solver.tnnmg.main.nu2")); - multigridStep->ignoreNodes_ = &ignoreNodes; + NonlinearSmootherType nonlinearSmoother; + TruncatedNonsmoothNewtonMultigrid<MyBlockProblemType, NonlinearSmootherType> + multigridStep(linearIterationStep, nonlinearSmoother); + multigridStep.setSmoothingSteps(parset.get<int>("solver.tnnmg.main.nu1"), + parset.get<int>("solver.tnnmg.main.mu"), + parset.get<int>("solver.tnnmg.main.nu2")); + multigridStep.ignoreNodes_ = &ignoreNodes; // }}} for (size_t run = 1; run <= parset.get<size_t>("timesteps"); ++run) { @@ -376,10 +367,10 @@ int main(int argc, char *argv[]) { } MyBlockProblemType myBlockProblem(parset, myConvexProblem); - multigridStep->setProblem(u4_diff, myBlockProblem); + multigridStep.setProblem(u4_diff, myBlockProblem); LoopSolver<VectorType> overallSolver( - multigridStep, solver_maxIterations, solver_tolerance, &energyNorm, + &multigridStep, solver_maxIterations, solver_tolerance, &energyNorm, verbosity); overallSolver.solve(); } @@ -437,18 +428,9 @@ int main(int argc, char *argv[]) { } std::cout << std::endl; - // {{{ Clean up after the TNNMG solver - delete linearBaseSolverStep; - delete baseEnergyNorm; - delete linearBaseSolver; - delete linearPresmoother; - delete linearPostsmoother; + // Clean up after the TNNMG solver for (auto &x : transferOperators) delete x; - delete linearIterationStep; - delete nonlinearSmoother; - delete multigridStep; - // }}} if (parset.get<bool>("printDifference")) { VectorType diff2 = u2; -- GitLab