Finite Element Domain Decomposition Library
FEDDLib
Loading...
Searching...
No Matches
NonLinLaplace_def.hpp
1#ifndef NonLinLaplace_def_hpp
2#define NonLinLaplace_def_hpp
3#include "NonLinLaplace_decl.hpp"
12
13namespace FEDD {
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)
18 : NonLinearProblem<SC, LO, GO, NO>(parameterList, domain->getComm()),
19 u_rep_() {
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()));
27}
28
29template <class SC, class LO, class GO, class NO>
30NonLinLaplace<SC, LO, GO, NO>::~NonLinLaplace() {}
31
32template <class SC, class LO, class GO, class NO>
33void NonLinLaplace<SC, LO, GO, NO>::info() {
34 this->infoProblem();
35 this->infoNonlinProblem();
36}
37
38/*
39 * General assemble function that is called whenever assembly of any kind is
40 * required type specifies whether the residual or the Jacobian is needed
41 */
42template <class SC, class LO, class GO, class NO>
43void NonLinLaplace<SC, LO, GO, NO>::assemble(std::string type) const {
44 if (type == "") {
45 if (this->verbose_) {
46 std::cout << "-- Assembly nonlinear laplace ... " << std::flush;
47 }
48 this->initAssemble();
49 } else {
50 this->reAssemble(type);
51 }
52 if (this->verbose_) {
53 std::cout << "done -- " << std::endl;
54 }
55}
56
57/*
58 * Initial assembly of the system matrix (Jacobian)
59 */
60template <class SC, class LO, class GO, class NO>
61void NonLinLaplace<SC, LO, GO, NO>::initAssemble() const {
62
63 if (this->verbose_) {
64 std::cout << "-- Initial assembly " << std::flush;
65 }
66
67 if (this->system_.is_null())
68 this->system_.reset(new BlockMatrix_Type(1));
69
70 if (this->residualVec_.is_null())
71 this->residualVec_.reset(new BlockMultiVector_Type(1));
72
73 this->feFactory_->assemblyNonlinearLaplace(
74 this->dim_, this->getDomain(0)->getFEType(), 2, this->u_rep_,
75 this->system_, this->residualVec_, this->parameterList_, "Jacobian");
76
77 // Initialise solution to 1 everywhere
78 this->solution_->putScalar(1.);
79
80 if (this->verbose_) {
81 std::cout << "done -- " << std::endl;
82 }
83}
84
85/*
86 * Reassembly of the residual or Jacobian
87 * Required by the nonlinear solver
88 */
89template <class SC, class LO, class GO, class NO>
90void NonLinLaplace<SC, LO, GO, NO>::reAssemble(std::string type) const {
91
92 if (this->verbose_)
93 std::cout << "-- Reassembly nonlinear laplace"
94 << " (" << type << ") ... " << std::flush;
95 // Update the locally stored solution to the problem
96 MultiVectorConstPtr_Type u = this->solution_->getBlock(0);
97 this->u_rep_->importFromVector(u, true);
98
99 if (type == "Rhs") {
100
101 this->feFactory_->assemblyNonlinearLaplace(
102 this->dim_, this->getDomain(0)->getFEType(), 2, this->u_rep_,
103 this->system_, this->residualVec_, this->parameterList_, "Rhs");
104
105 } else if (type == "Newton") {
106
107 this->feFactory_->assemblyNonlinearLaplace(
108 this->dim_, this->getDomain(0)->getFEType(), 2, this->u_rep_,
109 this->system_, this->residualVec_, this->parameterList_,
110 "Jacobian");
111 }
112 if (this->verbose_)
113 std::cout << "done -- " << std::endl;
114}
115
116template <class SC, class LO, class GO, class NO>
117void NonLinLaplace<SC, LO, GO, NO>::reAssembleExtrapolation(
118 BlockMultiVectorPtrArray_Type previousSolutions) {
119
120 TEUCHOS_TEST_FOR_EXCEPTION(
121 true, std::logic_error,
122 "Only Newton/NOX implemented for nonlinear Laplace!");
123}
124
125
126/*
127 * Updates the residual vector by reassembling and applying boundary conditions
128 */
129template <class SC, class LO, class GO, class NO>
130void NonLinLaplace<SC, LO, GO, NO>::calculateNonLinResidualVec(
131 std::string type, double time) const {
132
133 this->reAssemble("Rhs");
134 // rhs_ contains the dirichlet boundary conditions
135 // Apparently so does bcFactory_. Not sure why both do.
136 if (!type.compare("standard")) {
137 // this = 1*this - 1*rhs
138 this->residualVec_->update(-1., *this->rhs_, 1.);
139 this->bcFactory_->setVectorMinusBC(this->residualVec_, this->solution_,
140 time);
141 } else if (!type.compare("reverse")) {
142 // this = -1*this + 1*rhs
143
144 this->residualVec_->update(1., *this->rhs_, -1.);
145
146 // Sets the residualVec_ at the boundary = boundary condition - solution
147 // Necessary since reAssemble("Rhs") only assembles the residual on
148 // internal nodes
149 this->bcFactory_->setBCMinusVector(this->residualVec_, this->solution_,
150 time);
151 } else {
152 TEUCHOS_TEST_FOR_EXCEPTION(true, std::runtime_error,
153 "Unknown type for residual computation.");
154 }
155}
156} // namespace FEDD
157#endif
Definition NonLinearProblem_decl.hpp:24
Adaptive Mesh Refinement.
Definition AdaptiveMeshRefinement_decl.hpp:33