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