3#include "Laplace_decl.hpp"
15template<
class SC,
class LO,
class GO,
class NO>
16Laplace<SC,LO,GO,NO>::Laplace(
const DomainConstPtr_Type &domain, std::string FEType, ParameterListPtr_Type parameterList,
bool vectorLaplace):
17Problem<SC,LO,GO,NO>(parameterList, domain->getComm()),
18vectorLaplace_(vectorLaplace)
21 this->addVariable( domain , FEType ,
"u" , 1);
22 this->dim_ = this->getDomain(0)->getDimension();
25template<
class SC,
class LO,
class GO,
class NO>
26Laplace<SC,LO,GO,NO>::~Laplace(){
30template<
class SC,
class LO,
class GO,
class NO>
31void Laplace<SC,LO,GO,NO>::info(){
35template<
class SC,
class LO,
class GO,
class NO>
36void Laplace<SC,LO,GO,NO>::assemble( std::string type )
const{
39 std::cout <<
"-- Assembly Laplace ... " << std::flush;
42 vec_dbl_Type funcParameter(1,0.);
44 A = Teuchos::rcp(
new Matrix_Type( this->domainPtr_vec_.at(0)->getMapVecFieldUnique(), this->getDomain(0)->getApproxEntriesPerRow() ) );
45 this->feFactory_->assemblyLaplaceVecField(this->dim_, this->domain_FEType_vec_.at(0), 2, A );
48 A = Teuchos::rcp(
new Matrix_Type( this->domainPtr_vec_.at(0)->getMapUnique(), this->getDomain(0)->getApproxEntriesPerRow() ) );
49 this->feFactory_->assemblyLaplace(this->dim_, this->domain_FEType_vec_.at(0), 2, A );
52 this->system_->addBlock(A,0,0);
54 this->assembleSourceTerm( 0. );
56 this->addToRhs( this->sourceTerm_ );
59 std::cout <<
"done -- " << std::endl;
62template<
class SC,
class LO,
class GO,
class NO>
63typename Laplace<SC,LO,GO,NO>::MatrixPtr_Type Laplace<SC,LO,GO,NO>::getMassMatrix()
const{
66 A = Teuchos::rcp(
new Matrix_Type( this->domainPtr_vec_.at(0)->getMapUnique(), this->getDomain(0)->getApproxEntriesPerRow() ) );
67 this->feFactory_->assemblyMass(this->dim_,this->domain_FEType_vec_.at(0),
"Scalar", A);
Definition Problem_decl.hpp:38
Adaptive Mesh Refinement.
Definition AdaptiveMeshRefinement.cpp:5