13#include "feddlib/core/FE/Domain.hpp"
14#include "feddlib/core/FE/FE.hpp"
18template<
class SC,
class LO,
class GO,
class NO>
19Laplace<SC,LO,GO,NO>::Laplace(
const DomainConstPtr_Type &domain, std::string FEType, ParameterListPtr_Type parameterList,
bool vectorLaplace):
20Problem<SC,LO,GO,NO>(parameterList, domain->getComm()),
21vectorLaplace_(vectorLaplace)
24 this->addVariable( domain , FEType ,
"u" , 1);
25 this->dim_ = this->getDomain(0)->getDimension();
28template<
class SC,
class LO,
class GO,
class NO>
29Laplace<SC,LO,GO,NO>::~Laplace(){
33template<
class SC,
class LO,
class GO,
class NO>
34void Laplace<SC,LO,GO,NO>::info(){
38template<
class SC,
class LO,
class GO,
class NO>
39void Laplace<SC,LO,GO,NO>::assemble( std::string type )
const{
42 std::cout <<
"-- Assembly Laplace ... " << std::flush;
45 vec_dbl_Type funcParameter(1,0.);
47 A = Teuchos::rcp(
new Matrix_Type( this->domainPtr_vec_.at(0)->getMapVecFieldUnique(), this->getDomain(0)->getApproxEntriesPerRow() ) );
48 this->feFactory_->assemblyLaplaceVecField(this->dim_, this->domain_FEType_vec_.at(0), 2, A );
51 A = Teuchos::rcp(
new Matrix_Type( this->domainPtr_vec_.at(0)->getMapUnique(), this->getDomain(0)->getApproxEntriesPerRow() ) );
52 this->feFactory_->assemblyLaplace(this->dim_, this->domain_FEType_vec_.at(0), 2, A );
55 this->system_->addBlock(A,0,0);
57 this->assembleSourceTerm( 0. );
59 this->addToRhs( this->sourceTerm_ );
62 std::cout <<
"done -- " << std::endl;
65template<
class SC,
class LO,
class GO,
class NO>
66typename Laplace<SC,LO,GO,NO>::MatrixPtr_Type Laplace<SC,LO,GO,NO>::getMassMatrix()
const{
69 A = Teuchos::rcp(
new Matrix_Type( this->domainPtr_vec_.at(0)->getMapUnique(), this->getDomain(0)->getApproxEntriesPerRow() ) );
70 this->feFactory_->assemblyMass(this->dim_,this->domain_FEType_vec_.at(0),
"Scalar", A);
Definition Problem_decl.hpp:42
Adaptive Mesh Refinement.
Definition AdaptiveMeshRefinement_decl.hpp:36