Chi-Tech
lbs_DSA_preconditioner.cc
Go to the documentation of this file.
2
6
7//###################################################################
8/**Applies WGDSA or TGDSA to the given input vector.*/
9int lbs::WGDSA_TGDSA_PreConditionerMult(PC pc, Vec phi_input, Vec pc_output)
10{
11 void* context;
12 PCShellGetContext(pc,&context);
13
14 auto gs_context_ptr = (lbs::WGSContext<Mat,Vec,KSP>*)(context);
15
16 //Shorten some names
17 lbs::LBSSolver& lbs_solver = gs_context_ptr->lbs_solver_;
18 LBSGroupset& groupset = gs_context_ptr->groupset_;
19
20 //============================================= Copy PETSc vector to STL
21 auto& phi_new_local = gs_context_ptr->lbs_solver_.PhiNewLocal();
22 lbs_solver.SetPrimarySTLvectorFromGSPETScVec(groupset, phi_input,
23 PhiSTLOption::PHI_NEW);
24
25 //============================================= Apply WGDSA
26 if (groupset.apply_wgdsa_)
27 {
28 std::vector<double> delta_phi_local;
29 lbs_solver.AssembleWGDSADeltaPhiVector(groupset,phi_new_local, //From
30 delta_phi_local); //To
31
32 groupset.wgdsa_solver_->Assemble_b(delta_phi_local);
33 groupset.wgdsa_solver_->Solve(delta_phi_local);
34
35 lbs_solver.DisAssembleWGDSADeltaPhiVector(groupset, delta_phi_local, //From
36 phi_new_local); //To
37 }
38 //============================================= Apply TGDSA
39 if (groupset.apply_tgdsa_)
40 {
41 std::vector<double> delta_phi_local;
42 lbs_solver.AssembleTGDSADeltaPhiVector(groupset, phi_new_local, //From
43 delta_phi_local); //To
44
45 groupset.tgdsa_solver_->Assemble_b(delta_phi_local);
46 groupset.tgdsa_solver_->Solve(delta_phi_local);
47
48 lbs_solver.DisAssembleTGDSADeltaPhiVector(groupset, delta_phi_local, //From
49 phi_new_local); //To
50 }
51
52 //============================================= Copy STL vector to PETSc Vec
53 lbs_solver.SetGSPETScVecFromPrimarySTLvector(groupset, pc_output,
54 PhiSTLOption::PHI_NEW);
55
56 return 0;
57}
58
59//###################################################################
60/**Applies WGDSA or TGDSA to the given input vector.*/
62 lbs::WGSContext<Mat,Vec,KSP>& gs_context_ptr,
63 Vec phi_input, Vec pc_output)
64{
65 //Shorten some names
66 lbs::LBSSolver& lbs_solver = gs_context_ptr.lbs_solver_;
67 LBSGroupset& groupset = gs_context_ptr.groupset_;
68
69 //============================================= Copy PETSc vector to STL
70 auto& phi_new_local = gs_context_ptr.lbs_solver_.PhiNewLocal();
71 lbs_solver.SetPrimarySTLvectorFromGSPETScVec(groupset, phi_input,
72 PhiSTLOption::PHI_NEW);
73
74 //============================================= Apply WGDSA
75 if (groupset.apply_wgdsa_)
76 {
77 std::vector<double> delta_phi_local;
78 lbs_solver.AssembleWGDSADeltaPhiVector(groupset,phi_new_local, //From
79 delta_phi_local); //To
80
81 groupset.wgdsa_solver_->Assemble_b(delta_phi_local);
82 groupset.wgdsa_solver_->Solve(delta_phi_local);
83
84 lbs_solver.DisAssembleWGDSADeltaPhiVector(groupset, delta_phi_local, //From
85 phi_new_local); //To
86 }
87 //============================================= Apply TGDSA
88 if (groupset.apply_tgdsa_)
89 {
90 std::vector<double> delta_phi_local;
91 lbs_solver.AssembleTGDSADeltaPhiVector(groupset, phi_new_local, //From
92 delta_phi_local); //To
93
94 groupset.tgdsa_solver_->Assemble_b(delta_phi_local);
95 groupset.tgdsa_solver_->Solve(delta_phi_local);
96
97 lbs_solver.DisAssembleTGDSADeltaPhiVector(groupset, delta_phi_local, //From
98 phi_new_local); //To
99 }
100
101 //============================================= Copy STL vector to PETSc Vec
102 lbs_solver.SetGSPETScVecFromPrimarySTLvector(groupset, pc_output,
103 PhiSTLOption::PHI_NEW);
104
105 return 0;
106}
std::shared_ptr< lbs::acceleration::DiffusionMIPSolver > wgdsa_solver_
Definition: lbs_groupset.h:72
std::shared_ptr< lbs::acceleration::DiffusionMIPSolver > tgdsa_solver_
Definition: lbs_groupset.h:73
std::vector< double > & PhiNewLocal()
void DisAssembleWGDSADeltaPhiVector(const LBSGroupset &groupset, const std::vector< double > &delta_phi_local, std::vector< double > &ref_phi_new)
void AssembleTGDSADeltaPhiVector(const LBSGroupset &groupset, const std::vector< double > &phi_in, std::vector< double > &delta_phi_local)
void DisAssembleTGDSADeltaPhiVector(const LBSGroupset &groupset, const std::vector< double > &delta_phi_local, std::vector< double > &ref_phi_new)
virtual void SetGSPETScVecFromPrimarySTLvector(LBSGroupset &groupset, Vec x, PhiSTLOption which_phi)
void AssembleWGDSADeltaPhiVector(const LBSGroupset &groupset, const std::vector< double > &phi_in, std::vector< double > &delta_phi_local)
virtual void SetPrimarySTLvectorFromGSPETScVec(LBSGroupset &groupset, Vec x_src, PhiSTLOption which_phi)
int WGDSA_TGDSA_PreConditionerMult(PC pc, Vec phi_input, Vec pc_output)
int WGDSA_TGDSA_PreConditionerMult2(lbs::WGSContext< Mat, Vec, KSP > &gs_context_ptr, Vec phi_input, Vec pc_output)
struct _p_Vec * Vec
LBSGroupset & groupset_
Definition: wgs_context.h:28
LBSSolver & lbs_solver_
Definition: wgs_context.h:27