1#ifndef ASSEMBLEFENONLINLAPLACE_DEF_hpp
2#define ASSEMBLEFENONLINLAPLACE_DEF_hpp
4#include "feddlib/core/FEDDCore.hpp"
5#include <Teuchos_Array.hpp>
6#include <Teuchos_ScalarTraitsDecl.hpp>
11void rhsFunc(
double *x,
double *res,
double *parameters) {
23template <
class SC,
class LO,
class GO,
class NO>
25 int flag, vec2D_dbl_Type nodesRefConfig, ParameterListPtr_Type params,
26 tuple_disk_vec_ptr_Type tuple)
27 :
AssembleFE<SC, LO, GO, NO>(flag, nodesRefConfig, params, tuple) {
29 this->FEType_ = std::get<1>(this->diskTuple_->at(0));
30 this->dofs_ = std::get<2>(this->diskTuple_->at(0));
32 this->numNodes_ = std::get<3>(this->diskTuple_->at(0));
33 this->dofsElement_ = this->numNodes_ * this->dofs_;
36 this->rhsFunc_ = rhsFunc;
44template <
class SC,
class LO,
class GO,
class NO>
47 SmallMatrixPtr_Type elementMatrix =
48 Teuchos::rcp(
new SmallMatrix_Type(this->dofsElement_));
49 assemblyNonLinLaplacian(elementMatrix);
51 this->jacobian_ = elementMatrix;
60template <
class SC,
class LO,
class GO,
class NO>
61void AssembleFENonLinLaplace<SC, LO, GO, NO>::assemblyNonLinLaplacian(
62 SmallMatrixPtr_Type &elementMatrix) {
64 int dim = this->getDim();
67 vec3D_dbl_ptr_Type dPhi;
68 vec2D_dbl_ptr_Type phi;
69 vec_dbl_ptr_Type weights = Teuchos::rcp(
new vec_dbl_Type(0));
81 buildTransformation(B);
82 detB = B.computeInverse(Binv);
83 absDetB = std::fabs(detB);
84 vec3D_dbl_Type dPhiTrans(
86 vec2D_dbl_Type(dPhi->at(0).size(), vec_dbl_Type(dim, 0.)));
88 applyBTinv(dPhi, dPhiTrans, Binv);
90 vec_dbl_Type uLoc(weights->size(), 0.);
93 vec2D_dbl_Type duLoc(weights->size(), vec_dbl_Type(dim, 0));
97 for (
int w = 0; w < phi->size(); w++) {
99 for (
int i = 0; i < phi->at(0).size(); i++) {
100 uLoc[w] += (*this->solution_)[i] * phi->at(w).at(i);
101 for (
int d = 0; d < dim; d++) {
102 duLoc[w][d] += (*this->solution_)[i] * dPhiTrans[w][i][d];
108 auto value = Teuchos::ScalarTraits<SC>::zero();
109 for (UN i = 0; i < this->numNodes_; i++) {
110 for (UN j = 0; j < this->numNodes_; j++) {
112 for (UN w = 0; w < weights->size(); w++) {
113 for (UN d = 0; d < dim; d++) {
114 value += weights->at(w) *
115 ((1 + uLoc[w] * uLoc[w]) * dPhiTrans[w][i][d] +
116 2 * uLoc[w] * phi->at(w).at(j) * duLoc[w][d]) *
121 (*elementMatrix)[i][j] = value;
132template <
class SC,
class LO,
class GO,
class NO>
138 vec3D_dbl_ptr_Type dPhi;
139 vec2D_dbl_ptr_Type phi;
140 vec_dbl_ptr_Type weights = Teuchos::rcp(
new vec_dbl_Type(0));
149 buildTransformation(B);
150 detB = B.computeInverse(Binv);
151 absDetB = std::fabs(detB);
152 vec3D_dbl_Type dPhiTrans(
154 vec2D_dbl_Type(dPhi->at(0).size(), vec_dbl_Type(dim, 0.)));
155 applyBTinv(dPhi, dPhiTrans, Binv);
156 vec_dbl_Type uLoc(weights->size(), 0.);
157 vec_dbl_ptr_Type funcValues;
158 Helper::getFuncAtQuadNodes(funcValues, this->rhsFunc_, dim, this->FEType_,
162 for (
int w = 0; w < phi->size(); w++) {
163 for (
int i = 0; i < phi->at(0).size();
165 uLoc[w] += (*this->solution_)[i] * phi->at(w).at(i);
169 this->rhsVec_.reset(
new vec_dbl_Type(this->dofsElement_, 0.));
171 auto valueOne = Teuchos::ScalarTraits<SC>::zero();
173 auto valueTwo = Teuchos::ScalarTraits<SC>::zero();
175 auto reductionValue = Teuchos::ScalarTraits<SC>::zero();
177 for (UN i = 0; i < this->numNodes_; i++) {
180 for (UN j = 0; j < this->numNodes_; j++) {
183 for (UN w = 0; w < dPhiTrans.size(); w++) {
184 for (UN d = 0; d < dim; d++) {
185 valueTwo += weights->at(w) * (1 + uLoc[w] * uLoc[w]) *
186 dPhiTrans[w][i][d] * dPhiTrans[w][j][d];
190 reductionValue += valueTwo * (*this->solution_)[j];
192 reductionValue *= absDetB;
194 for (UN w = 0; w < dPhiTrans.size(); w++) {
196 valueOne += weights->at(w) * funcValues->at(w) * phi->at(w).at(i);
199 (*this->rhsVec_)[i] = -(valueOne - reductionValue);
208template <
class SC,
class LO,
class GO,
class NO>
209void AssembleFENonLinLaplace<SC, LO, GO, NO>::buildTransformation(
212 TEUCHOS_TEST_FOR_EXCEPTION((B.size() < 2 || B.size() > 3), std::logic_error,
213 "Initialize SmallMatrix for transformation.");
216 for (UN j = 0; j < B.size(); j++) {
218 for (UN i = 0; i < B.size(); i++) {
221 B[i][j] = this->nodesRefConfig_.at(index).at(i) -
222 this->nodesRefConfig_.at(index0).at(i);
233template <
class SC,
class LO,
class GO,
class NO>
234void AssembleFENonLinLaplace<SC, LO, GO, NO>::applyBTinv(
235 vec3D_dbl_ptr_Type &dPhiIn, vec3D_dbl_Type &dPhiOut,
237 UN dim = Binv.size();
238 for (UN w = 0; w < dPhiIn->size(); w++) {
239 for (UN i = 0; i < dPhiIn->at(w).size(); i++) {
240 for (UN d1 = 0; d1 < dim; d1++) {
241 for (UN d2 = 0; d2 < dim; d2++) {
243 dPhiIn->at(w).at(i).at(d2) * Binv[d2][d1];
void assembleRHS() override
Assemble the element right hand side vector.
Definition AssembleFENonLinLaplace_def.hpp:133
void assembleJacobian() override
Assemble the element Jacobian matrix.
Definition AssembleFENonLinLaplace_def.hpp:45
AssembleFENonLinLaplace(int flag, vec2D_dbl_Type nodesRefConfig, ParameterListPtr_Type parameters, tuple_disk_vec_ptr_Type tuple)
Constructor for AssembleFENonLinLaplace.
Definition AssembleFENonLinLaplace_def.hpp:24
int getDim()
Definition AssembleFE_def.hpp:135
AssembleFE(int flag, vec2D_dbl_Type nodesRefConfig, ParameterListPtr_Type parameters, tuple_disk_vec_ptr_Type tuple)
Definition AssembleFE_def.hpp:8
static UN determineDegree(UN dim, std::string FEType, VarType orderOfDerivative)
Determine polynomial degree of a finite element basis function or its gradient that is required to se...
Definition Helper.cpp:68
@ Deriv1
order 1, gradient(f(x))
Definition Helper.hpp:28
static int getPhi(vec2D_dbl_ptr_Type &Phi, vec_dbl_ptr_Type &weightsPhi, int dim, std::string FEType, int Degree, std::string FETypeQuadPoints="")
Get basisfunction phi per quadrature point.
Definition Helper.cpp:921
static int getDPhi(vec3D_dbl_ptr_Type &DPhi, vec_dbl_ptr_Type &weightsDPhi, int Dimension, std::string FEType, int Degree)
Full matrix representation of gradient of a basis function for each quadrature point.
Definition Helper.cpp:215
This class represents a templated small Matrix of type T. Primarily created for 2x2 and 3x3 matrices....
Definition SmallMatrix.hpp:20
Adaptive Mesh Refinement.
Definition AdaptiveMeshRefinement_decl.hpp:36