Chi-Tech
lbsmip_TGDSA_preconditioner.cc
Go to the documentation of this file.
2
6
7//###################################################################
8/**Applies TGDSA to the given input vector.*/
9int lbs::MIP_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& solver = gs_context_ptr->lbs_solver_;
18 LBSGroupset& groupset = gs_context_ptr->groupset_;
19
20 //============================================= Copy PETSc vector to STL
21 auto& phi_delta = gs_context_ptr->lbs_solver_.PhiNewLocal();
22 solver.SetPrimarySTLvectorFromGSPETScVec(groupset, phi_input,
23 PhiSTLOption::PHI_NEW);
24
25 //============================================= Apply TGDSA
26 if (groupset.apply_tgdsa_)
27 {
28 std::vector<double> delta_phi_local;
29 solver.AssembleTGDSADeltaPhiVector(groupset, phi_delta, delta_phi_local);
30 groupset.tgdsa_solver_->Assemble_b(delta_phi_local);
31 groupset.tgdsa_solver_->Solve(delta_phi_local);
32 solver.DisAssembleTGDSADeltaPhiVector(groupset, delta_phi_local, phi_delta);
33 }
34
35 //============================================= Copy STL vector to PETSc Vec
37 pc_output,
38 PhiSTLOption::PHI_NEW);
39
40 return 0;
41}
std::shared_ptr< lbs::acceleration::DiffusionMIPSolver > tgdsa_solver_
Definition: lbs_groupset.h:73
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)
virtual void SetPrimarySTLvectorFromGSPETScVec(LBSGroupset &groupset, Vec x_src, PhiSTLOption which_phi)
int MIP_TGDSA_PreConditionerMult(PC pc, Vec phi_input, Vec pc_output)
struct _p_Vec * Vec