3#include "Stokes_decl.hpp"
14void ZeroDirichlet(
double* x,
double* res,
double t,
double* parameters){
21double OneFunction(
double* x,
int* parameter)
26template<
class SC,
class LO,
class GO,
class NO>
27Stokes<SC,LO,GO,NO>::Stokes(
const DomainConstPtr_Type &domainVelocity, std::string FETypeVelocity,
const DomainConstPtr_Type &domainPressure, std::string FETypePressure, ParameterListPtr_Type parameterList):
28Problem<SC,LO,GO,NO>(parameterList, domainVelocity->getComm())
31 this->addVariable( domainVelocity , FETypeVelocity ,
"u" , domainVelocity->getDimension());
32 this->addVariable( domainPressure , FETypePressure ,
"p" , 1);
34 this->dim_ = this->getDomain(0)->getDimension();
37template<
class SC,
class LO,
class GO,
class NO>
38Stokes<SC,LO,GO,NO>::~Stokes(){
42template<
class SC,
class LO,
class GO,
class NO>
43void Stokes<SC,LO,GO,NO>::info(){
47template<
class SC,
class LO,
class GO,
class NO>
48void Stokes<SC,LO,GO,NO>::assemble( std::string type )
const{
51 std::cout <<
"-- Assembly ... " << std::flush;
53 double viscosity = this->parameterList_->sublist(
"Parameter").get(
"Viscosity",1.);
55 MatrixPtr_Type A(
new Matrix_Type( this->getDomain(0)->getMapVecFieldUnique(), this->getDomain(0)->getApproxEntriesPerRow() ) );
56 MatrixPtr_Type BT(
new Matrix_Type( this->getDomain(0)->getMapVecFieldUnique(), this->getDomain(1)->getDimension() * this->getDomain(1)->getApproxEntriesPerRow() ) );
58 MapConstPtr_Type pressureMap;
59 if ( this->getDomain(1)->getFEType() ==
"P0" )
60 pressureMap = this->getDomain(1)->getElementMap();
62 pressureMap = this->getDomain(1)->getMapUnique();
64 MatrixPtr_Type B(
new Matrix_Type( pressureMap, this->getDomain(0)->getDimension() * this->getDomain(0)->getApproxEntriesPerRow() ) );
68 std::cout <<
" A ... " << std::flush;
70 if ( this->parameterList_->sublist(
"Parameter").get(
"Symmetric gradient",
false) )
71 this->feFactory_->assemblyStress(this->dim_, this->domain_FEType_vec_.at(0), A, OneFunction, dummy,
true);
73 this->feFactory_->assemblyLaplaceVecField(this->dim_, this->domain_FEType_vec_.at(0), 2, A,
true);
76 std::cout <<
"B and B^T ... " << std::flush;
78 this->feFactory_->assemblyDivAndDivT(this->dim_, this->getFEType(0), this->getFEType(1), 2, B, BT, this->getDomain(0)->getMapVecFieldUnique(), pressureMap,
true );
88 A->fillComplete( this->getDomain(0)->getMapVecFieldUnique(), this->getDomain(0)->getMapVecFieldUnique());
89 B->fillComplete( this->getDomain(0)->getMapVecFieldUnique(), pressureMap );
90 BT->fillComplete( pressureMap, this->getDomain(0)->getMapVecFieldUnique() );
92 this->system_.reset(
new BlockMatrix_Type(2));
94 this->system_->addBlock( A, 0, 0 );
95 this->system_->addBlock( BT, 0, 1 );
96 this->system_->addBlock( B, 1, 0 );
100 if ( !this->getFEType(0).compare(
"P1") ) {
101 C.reset(
new Matrix_Type( this->getDomain(1)->getMapUnique(), this->getDomain(1)->getApproxEntriesPerRow() ) );
102 this->feFactory_->assemblyBDStabilization( this->dim_,
"P1", C,
true);
104 C->scale( -1./viscosity );
105 C->fillComplete( pressureMap, pressureMap );
106 this->system_->addBlock( C, 1, 1 );
109 if ( !this->parameterList_->sublist(
"General").get(
"Preconditioner Method",
"Monolithic").compare(
"Teko") ) {
110 if (!this->parameterList_->sublist(
"General").get(
"Assemble Velocity Mass",
false)) {
111 MatrixPtr_Type Mvelocity(
new Matrix_Type( this->getDomain(0)->getMapVecFieldUnique(), this->getDomain(0)->getApproxEntriesPerRow() ) );
113 this->feFactory_->assemblyMass( this->dim_, this->domain_FEType_vec_.at(0),
"Vector", Mvelocity,
true );
115 this->getPreconditionerConst()->setVelocityMassMatrix( Mvelocity );
117 std::cout <<
"\nVelocity mass matrix for LSC block preconditioner is assembled." << std::endl;
120 std::cout <<
"\nVelocity mass matrix for LSC block preconditioner not assembled." << std::endl;
124 std::string precType = this->parameterList_->sublist(
"General").get(
"Preconditioner Method",
"Monolithic");
125 if ( precType ==
"Diagonal" || precType ==
"Triangular" ) {
126 MatrixPtr_Type Mpressure(
new Matrix_Type( this->getDomain(1)->getMapUnique(), this->getDomain(1)->getApproxEntriesPerRow() ) );
128 this->feFactory_->assemblyMass( this->dim_, this->domain_FEType_vec_.at(1),
"Scalar", Mpressure,
true );
130 Mpressure->resumeFill();
131 Mpressure->scale(-1./viscosity);
132 Mpressure->fillComplete( pressureMap, pressureMap );
133 this->getPreconditionerConst()->setPressureMassMatrix( Mpressure );
136 this->assembleSourceTerm( 0. );
137 this->addToRhs( this->sourceTerm_ );
141 std::cout <<
"done -- " << std::endl;
Definition Problem_decl.hpp:38
Adaptive Mesh Refinement.
Definition AdaptiveMeshRefinement.cpp:5