1#ifndef NonLinLaplace_def_hpp
2#define NonLinLaplace_def_hpp
13#include "feddlib/core/FE/Domain.hpp"
14#include "feddlib//core/FE/FE.hpp"
15#include "feddlib/core/General/BCBuilder.hpp"
19template <
class SC,
class LO,
class GO,
class NO>
20NonLinLaplace<SC, LO, GO, NO>::NonLinLaplace(
21 const DomainConstPtr_Type &domain, std::string FEType,
22 ParameterListPtr_Type parameterList)
25 this->nonLinearTolerance_ =
26 this->parameterList_->sublist(
"Parameter").get(
"relNonLinTol", 1.0e-6);
27 this->initNOXParameters();
28 this->addVariable(domain, FEType,
"u", 1);
29 this->dim_ = this->getDomain(0)->getDimension();
30 u_rep_ = Teuchos::rcp(
31 new MultiVector_Type(this->getDomain(0)->getMapRepeated()));
34template <
class SC,
class LO,
class GO,
class NO>
35NonLinLaplace<SC, LO, GO, NO>::~NonLinLaplace() {}
37template <
class SC,
class LO,
class GO,
class NO>
38void NonLinLaplace<SC, LO, GO, NO>::info() {
40 this->infoNonlinProblem();
47template <
class SC,
class LO,
class GO,
class NO>
51 std::cout <<
"-- Assembly nonlinear laplace ... " << std::flush;
55 this->reAssemble(type);
58 std::cout <<
"done -- " << std::endl;
65template <
class SC,
class LO,
class GO,
class NO>
66void NonLinLaplace<SC, LO, GO, NO>::initAssemble()
const {
69 std::cout <<
"-- Initial assembly " << std::flush;
72 if (this->system_.is_null())
73 this->system_.reset(
new BlockMatrix_Type(1));
75 if (this->residualVec_.is_null())
76 this->residualVec_.reset(
new BlockMultiVector_Type(1));
78 this->feFactory_->assemblyNonlinearLaplace(
79 this->dim_, this->getDomain(0)->getFEType(), 2, this->u_rep_,
80 this->system_, this->residualVec_, this->parameterList_,
"Jacobian");
83 this->solution_->putScalar(1.);
86 std::cout <<
"done -- " << std::endl;
94template <
class SC,
class LO,
class GO,
class NO>
95void NonLinLaplace<SC, LO, GO, NO>::reAssemble(std::string type)
const {
98 std::cout <<
"-- Reassembly nonlinear laplace"
99 <<
" (" << type <<
") ... " << std::flush;
101 MultiVectorConstPtr_Type u = this->solution_->getBlock(0);
102 this->u_rep_->importFromVector(u,
true);
106 this->feFactory_->assemblyNonlinearLaplace(
107 this->dim_, this->getDomain(0)->getFEType(), 2, this->u_rep_,
108 this->system_, this->residualVec_, this->parameterList_,
"Rhs");
110 }
else if (type ==
"Newton") {
112 this->feFactory_->assemblyNonlinearLaplace(
113 this->dim_, this->getDomain(0)->getFEType(), 2, this->u_rep_,
114 this->system_, this->residualVec_, this->parameterList_,
118 std::cout <<
"done -- " << std::endl;
121template <
class SC,
class LO,
class GO,
class NO>
122void NonLinLaplace<SC, LO, GO, NO>::reAssembleExtrapolation(
123 BlockMultiVectorPtrArray_Type previousSolutions) {
125 TEUCHOS_TEST_FOR_EXCEPTION(
126 true, std::logic_error,
127 "Only Newton/NOX implemented for nonlinear Laplace!");
134template <
class SC,
class LO,
class GO,
class NO>
136 std::string type,
double time)
const {
138 this->reAssemble(
"Rhs");
141 if (!type.compare(
"standard")) {
143 this->residualVec_->update(-1., *this->rhs_, 1.);
144 this->bcFactory_->setVectorMinusBC(this->residualVec_, this->solution_,
146 }
else if (!type.compare(
"reverse")) {
149 this->residualVec_->update(1., *this->rhs_, -1.);
154 this->bcFactory_->setBCMinusVector(this->residualVec_, this->solution_,
157 TEUCHOS_TEST_FOR_EXCEPTION(
true, std::runtime_error,
158 "Unknown type for residual computation.");
virtual void calculateNonLinResidualVec(std::string type, double time=0.) const
Virtual function which is implemented in the specific non-linear problem classes to calculate the non...
Definition NonLinLaplace_def.hpp:135
virtual void assemble(std::string type="") const
assemble of type exectuted by the derived specific non-linear problem classes
Definition NonLinLaplace_def.hpp:48
Definition NonLinearProblem_decl.hpp:28
Adaptive Mesh Refinement.
Definition AdaptiveMeshRefinement_decl.hpp:36