1#ifndef NonLinLaplace_def_hpp
2#define NonLinLaplace_def_hpp
3#include "NonLinLaplace_decl.hpp"
14template <
class SC,
class LO,
class GO,
class NO>
15NonLinLaplace<SC, LO, GO, NO>::NonLinLaplace(
16 const DomainConstPtr_Type &domain, std::string FEType,
17 ParameterListPtr_Type parameterList)
20 this->nonLinearTolerance_ =
21 this->parameterList_->sublist(
"Parameter").get(
"relNonLinTol", 1.0e-6);
22 this->initNOXParameters();
23 this->addVariable(domain, FEType,
"u", 1);
24 this->dim_ = this->getDomain(0)->getDimension();
25 u_rep_ = Teuchos::rcp(
26 new MultiVector_Type(this->getDomain(0)->getMapRepeated()));
29template <
class SC,
class LO,
class GO,
class NO>
30NonLinLaplace<SC, LO, GO, NO>::~NonLinLaplace() {}
32template <
class SC,
class LO,
class GO,
class NO>
33void NonLinLaplace<SC, LO, GO, NO>::info() {
35 this->infoNonlinProblem();
42template <
class SC,
class LO,
class GO,
class NO>
43void NonLinLaplace<SC, LO, GO, NO>::assemble(std::string type)
const {
46 std::cout <<
"-- Assembly nonlinear laplace ... " << std::flush;
50 this->reAssemble(type);
53 std::cout <<
"done -- " << std::endl;
60template <
class SC,
class LO,
class GO,
class NO>
61void NonLinLaplace<SC, LO, GO, NO>::initAssemble()
const {
64 std::cout <<
"-- Initial assembly " << std::flush;
67 if (this->system_.is_null())
68 this->system_.reset(
new BlockMatrix_Type(1));
70 if (this->residualVec_.is_null())
71 this->residualVec_.reset(
new BlockMultiVector_Type(1));
73 this->feFactory_->assemblyNonlinearLaplace(
74 this->dim_, this->getDomain(0)->getFEType(), 2, this->u_rep_,
75 this->system_, this->residualVec_, this->parameterList_,
"Jacobian");
78 this->solution_->putScalar(1.);
81 std::cout <<
"done -- " << std::endl;
89template <
class SC,
class LO,
class GO,
class NO>
90void NonLinLaplace<SC, LO, GO, NO>::reAssemble(std::string type)
const {
93 std::cout <<
"-- Reassembly nonlinear laplace"
94 <<
" (" << type <<
") ... " << std::flush;
96 MultiVectorConstPtr_Type u = this->solution_->getBlock(0);
97 this->u_rep_->importFromVector(u,
true);
101 this->feFactory_->assemblyNonlinearLaplace(
102 this->dim_, this->getDomain(0)->getFEType(), 2, this->u_rep_,
103 this->system_, this->residualVec_, this->parameterList_,
"Rhs");
105 }
else if (type ==
"Newton") {
107 this->feFactory_->assemblyNonlinearLaplace(
108 this->dim_, this->getDomain(0)->getFEType(), 2, this->u_rep_,
109 this->system_, this->residualVec_, this->parameterList_,
113 std::cout <<
"done -- " << std::endl;
116template <
class SC,
class LO,
class GO,
class NO>
117void NonLinLaplace<SC, LO, GO, NO>::reAssembleExtrapolation(
118 BlockMultiVectorPtrArray_Type previousSolutions) {
120 TEUCHOS_TEST_FOR_EXCEPTION(
121 true, std::logic_error,
122 "Only Newton/NOX implemented for nonlinear Laplace!");
129template <
class SC,
class LO,
class GO,
class NO>
130void NonLinLaplace<SC, LO, GO, NO>::calculateNonLinResidualVec(
131 std::string type,
double time)
const {
133 this->reAssemble(
"Rhs");
136 if (!type.compare(
"standard")) {
138 this->residualVec_->update(-1., *this->rhs_, 1.);
139 this->bcFactory_->setVectorMinusBC(this->residualVec_, this->solution_,
141 }
else if (!type.compare(
"reverse")) {
144 this->residualVec_->update(1., *this->rhs_, -1.);
149 this->bcFactory_->setBCMinusVector(this->residualVec_, this->solution_,
152 TEUCHOS_TEST_FOR_EXCEPTION(
true, std::runtime_error,
153 "Unknown type for residual computation.");
Definition NonLinearProblem_decl.hpp:24
Adaptive Mesh Refinement.
Definition AdaptiveMeshRefinement_decl.hpp:33