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