1#ifndef LaplaceBlocks_def_hpp
2#define LaplaceBlocks_def_hpp
13#include "feddlib/core/FE/Domain.hpp"
14#include "feddlib/core/FE/FE.hpp"
19template<
class SC,
class LO,
class GO,
class NO>
20LaplaceBlocks<SC,LO,GO,NO>::LaplaceBlocks(
const DomainConstPtr_Type &domain1,
const DomainConstPtr_Type &domain2, std::string FEType1, std::string FEType2, ParameterListPtr_Type parameterList ):
21Problem<SC,LO,GO,NO>(parameterList, domain1->getComm())
24 this->addVariable( domain1 , FEType1 ,
"u1" , 1);
25 this->addVariable( domain2 , FEType2 ,
"u2" , 1);
26 this->dim_ = this->getDomain(0)->getDimension();
29template<
class SC,
class LO,
class GO,
class NO>
30LaplaceBlocks<SC,LO,GO,NO>::~LaplaceBlocks(){
34template<
class SC,
class LO,
class GO,
class NO>
35void LaplaceBlocks<SC,LO,GO,NO>::info(){
39template<
class SC,
class LO,
class GO,
class NO>
40void LaplaceBlocks<SC,LO,GO,NO>::assemble( std::string type )
const{
43 std::cout <<
"-- Assembly LaplaceBlocks ... " << std::flush;
49 A1 = Teuchos::rcp(
new Matrix_Type( this->domainPtr_vec_.at(0)->getMapUnique(), this->getDomain(0)->getApproxEntriesPerRow() ) );
50 this->feFactory_->assemblyLaplace(this->dim_, this->domain_FEType_vec_.at(0), 2, A1,
true, 0);
52 A2 = Teuchos::rcp(
new Matrix_Type( this->domainPtr_vec_.at(1)->getMapUnique(), this->getDomain(0)->getApproxEntriesPerRow() ) );
53 this->feFactory_->assemblyLaplace(this->dim_, this->domain_FEType_vec_.at(1), 2, A2,
true , 1);
55 this->system_->addBlock( A1, 0, 0 );
56 this->system_->addBlock( A2, 1, 1 );
58 this->assembleSourceTerm( 0. );
60 this->addToRhs( this->sourceTerm_ );
63 std::cout <<
"done -- " << std::endl;
Definition Problem_decl.hpp:42
Adaptive Mesh Refinement.
Definition AdaptiveMeshRefinement_decl.hpp:36