Finite Element Domain Decomposition Library
FEDDLib
Loading...
Searching...
No Matches
LaplaceBlocks_def.hpp
1#ifndef LaplaceBlocks_def_hpp
2#define LaplaceBlocks_def_hpp
3#include "LaplaceBlocks_decl.hpp"
12
13namespace FEDD {
14
15template<class SC,class LO,class GO,class NO>
16LaplaceBlocks<SC,LO,GO,NO>::LaplaceBlocks( const DomainConstPtr_Type &domain1, const DomainConstPtr_Type &domain2, std::string FEType1, std::string FEType2, ParameterListPtr_Type parameterList ):
17Problem<SC,LO,GO,NO>(parameterList, domain1->getComm())
18{
19
20 this->addVariable( domain1 , FEType1 , "u1" , 1);
21 this->addVariable( domain2 , FEType2 , "u2" , 1);
22 this->dim_ = this->getDomain(0)->getDimension();
23}
24
25template<class SC,class LO,class GO,class NO>
26LaplaceBlocks<SC,LO,GO,NO>::~LaplaceBlocks(){
27
28}
29
30template<class SC,class LO,class GO,class NO>
31void LaplaceBlocks<SC,LO,GO,NO>::info(){
32 this->infoProblem();
33}
34
35template<class SC,class LO,class GO,class NO>
36void LaplaceBlocks<SC,LO,GO,NO>::assemble( std::string type ) const{
37
38 if (this->verbose_)
39 std::cout << "-- Assembly LaplaceBlocks ... " << std::flush;
40
41 MatrixPtr_Type A1;
42 MatrixPtr_Type A2;
43
44
45 A1 = Teuchos::rcp(new Matrix_Type( this->domainPtr_vec_.at(0)->getMapUnique(), this->getDomain(0)->getApproxEntriesPerRow() ) );
46 this->feFactory_->assemblyLaplace(this->dim_, this->domain_FEType_vec_.at(0), 2, A1, true, 0);
47
48 A2 = Teuchos::rcp(new Matrix_Type( this->domainPtr_vec_.at(1)->getMapUnique(), this->getDomain(0)->getApproxEntriesPerRow() ) );
49 this->feFactory_->assemblyLaplace(this->dim_, this->domain_FEType_vec_.at(1), 2, A2, true , 1);
50
51 this->system_->addBlock( A1, 0, 0 );
52 this->system_->addBlock( A2, 1, 1 );
53
54 this->assembleSourceTerm( 0. );
55
56 this->addToRhs( this->sourceTerm_ );
57
58 if (this->verbose_)
59 std::cout << "done -- " << std::endl;
60}
61
62//template<class SC,class LO,class GO,class NO>
63//void LaplaceBlocks<SC,LO,GO,NO>::assembleSourceTerm(double time){
64//
65// MultiVectorPtr_Type FERhs1;
66// MultiVectorPtr_Type FERhs2;
67//
68// vec_dbl_Type funcParameter(1,0.);
69// FERhs1 = Teuchos::rcp(new MultiVector_Type( this->domainPtr_vec_.at(0)->getMapRepeated() ) );
70// this->feFactory_->assemblyRHS(this->dim_,
71// this->domain_FEType_vec_.at(0),
72// FERhs1,
73// "Scalar",
74// this->rhsFuncVec_[0],
75// funcParameter,
76// 0);
77//
78// FERhs2 = Teuchos::rcp(new MultiVector_Type( this->domainPtr_vec_.at(1)->getMapRepeated() ) );
79// this->feFactory_->assemblyRHS(this->dim_,
80// this->domain_FEType_vec_.at(1),
81// FERhs2,
82// "Scalar",
83// this->rhsFuncVec_[0],
84// funcParameter,
85// 1);
86//
87//
88// this->sourceTerm_->getBlockNonConst(0)->putScalar(0.);
89//
90// this->sourceTerm_->getBlockNonConst(0)->exportFromVector( FERhs1, false, "Add" );
91//
92// this->sourceTerm_->getBlockNonConst(1)->putScalar(0.);
93//
94// this->sourceTerm_->getBlockNonConst(1)->exportFromVector( FERhs2, false, "Add" );
95//
96//}
97
98//template<class SC,class LO,class GO,class NO>
99//void LaplaceBlocks<SC,LO,GO,NO>::assembleRhs(double time){
100//
101// MultiVectorPtr_Type FERhs1;
102// MultiVectorPtr_Type FERhs2;
103//
104// vec_dbl_Type funcParameter(1,0.);
105// FERhs1 = Teuchos::rcp(new MultiVector_Type( this->domainPtr_vec_.at(0)->getMapRepeated() ) );
106// this->feFactory_->assemblyRHS(this->dim_,
107// this->domain_FEType_vec_.at(0),
108// FERhs1,
109// "Scalar",
110// this->rhsFuncVec_[0],
111// funcParameter,
112// 0);
113//
114// FERhs2 = Teuchos::rcp(new MultiVector_Type( this->domainPtr_vec_.at(1)->getMapRepeated() ) );
115// this->feFactory_->assemblyRHS(this->dim_,
116// this->domain_FEType_vec_.at(1),
117// FERhs2,
118// "Scalar",
119// this->rhsFuncVec_[0],
120// funcParameter,
121// 1);
122//
123//
124// this->sourceTerm_->getBlockNonConst(0)->putScalar(0.);
125//
126// this->sourceTerm_->getBlockNonConst(0)->exportFromVector( FERhs1, false, "Add" );
127//
128// this->sourceTerm_->getBlockNonConst(1)->putScalar(0.);
129//
130// this->sourceTerm_->getBlockNonConst(1)->exportFromVector( FERhs2, false, "Add" );
131//
132//}
133
134//template<class SC,class LO,class GO,class NO>
135//int LaplaceBlocks<SC,LO,GO,NO>::SetupPreconditioner(BMat_ptr_Type systemPrec){
136//
137//#ifdef ASSERTS_WARNINGS
138// MYASSERT(this->boolProblemAssembled_,"Problem not assembled!");
139//#endif
140//
141// Teuchos::RCP<PrecondManagerFROSch<SC,LO,GO,NO> > precondManager(new PrecondManagerFROSch<SC,LO,GO,NO>());
142//
143// precondManager->BuildPreconditioner( this->XpetraPrec_, this->XpetraMatrix_, systemPrec, this->domainPtr_vec_.at(0), this->bcFactory_ , this->parameterList_ , this->comm_);
144//
145//
146// return 0;
147//};
148
149//template<class SC,class LO,class GO,class NO>
150//int LaplaceBlocks<SC,LO,GO,NO>::SetupPreconditioner( BMat_ptr_Type systemPrec, ThyraConstLinOpPtr_Type thyraMatrix, ThyraPrecPtr_Type thyraPreconditioner, LinSolverBuilderPtr_Type linearSolverBuilder ) const{
151//
152//#ifdef ASSERTS_WARNINGS
153// MYASSERT(this->boolProblemAssembled_,"Problem not assembled!");
154//#endif
155//
156// Teuchos::RCP<PrecondManagerFROSch<SC,LO,GO,NO> > precondManager(new PrecondManagerFROSch<SC,LO,GO,NO>());
157//
158// if (!thyraPreconditioner.is_null() && !thyraMatrix.is_null() && !linearSolverBuilder.is_null()) {
159// precondManager->BuildPreconditioner( thyraPreconditioner, thyraMatrix, linearSolverBuilder ,systemPrec, this->domainPtr_vec_, this->bcFactory_ , this->parameterList_ , this->PListThyraSolver_, this->comm_, this->PrecondtionerIsBuilt_);
160// }
161// else{
162// precondManager->BuildPreconditioner( this->ThyraPreconditioner_, this->ThyraMatrix_, this->LinearSolverBuilder_ ,systemPrec, this->domainPtr_vec_, this->bcFactory_ , this->parameterList_ , this->PListThyraSolver_, this->comm_, this->PrecondtionerIsBuilt_);
163// }
164//
165//
166//
167// return 0;
168//};
169//
170//template<class SC,class LO,class GO,class NO>
171//void LaplaceBlocks<SC,LO,GO,NO>::initSolutionWithBC(){
172//
174//
175//}
176}
177#endif
Definition Problem_decl.hpp:38
Adaptive Mesh Refinement.
Definition AdaptiveMeshRefinement.cpp:5