35#ifndef __H__UG__LIB_DISC__OPERATOR__LINEAR_OPERATOR__GPUJACOBI__
36#define __H__UG__LIB_DISC__OPERATOR__LINEAR_OPERATOR__GPUJACOBI__
52CUDA_JacobiApply(
const double *diagInv,
double *corr,
const double *defect,
const int N);
58class GPUJacobi :
public IPreconditioner<GPUAlgebra>
62 typedef GPUAlgebra TAlgebra;
63 typedef TAlgebra algebra_type;
72 typedef IPreconditioner<TAlgebra>::matrix_operator_type matrix_operator_type;
75 typedef IPreconditioner<TAlgebra> base_type;
78 using base_type::set_debug;
79 using base_type::debug_writer;
80 using base_type::write_debug;
81 using base_type::damping;
85 GPUJacobi() {this->set_damp(1.0);};
88 GPUJacobi(
number damp) {this->set_damp(damp);};
91 GPUJacobi(
const GPUJacobi &parent )
99 return make_sp(
new GPUJacobi(*
this));
107 virtual bool supports_parallel()
const {
return false; }
110 virtual const char*
name()
const {
return "GPUJacobi";}
113 virtual bool preprocess(
SmartPtr<MatrixOperator<matrix_type, vector_type> > pOp)
119 size_t size = mat.num_rows();
120 if(size != mat.num_cols())
122 UG_LOG(
"Square Matrix needed for GPUJacobi Iteration.\n");
128 if(damping()->constant_damping())
129 damp = damping()->damping();
136 m_diagInv.resize(mat.num_rows());
138 for(
size_t i = 0; i < mat.num_rows(); ++i)
139 m_diagInv[i] = damp/mat(i,i);
141 CUDAHelper::get_instance();
146 virtual bool step(
SmartPtr<MatrixOperator<matrix_type, vector_type> > pOp, vector_type& c,
const vector_type& d)
150 c.set_storage_type(PST_ADDITIVE);
153 if(!c.change_storage_type(PST_CONSISTENT))
155 UG_LOG(
"ERROR in 'JacobiPreconditioner::apply': "
156 "Cannot change parallel status of correction to consistent.\n");
161 const double *devDiagInv = m_diagInv.get_dev_ptr();
162 double *devC = c.get_dev_ptr();
163 const double *devD = d.get_dev_ptr();
167 CUDA_JacobiApply(devDiagInv, devC, devD, N);
174 virtual bool postprocess() {
return true;}
177 GPUVector<double> m_diagInv;
182class Jacobi<GPUAlgebra> :
public GPUJacobi
184 using GPUJacobi::set_damp;
location name
Definition checkpoint_util.lua:128
Definition smart_pointer.h:108
void set_damp(SmartPtr< IDamping< TAlgebra::vector_type, TAlgebra::vector_type > > spScaling)
sets a scaling for the correction
Definition linear_iterator.h:163
Jacobi()
default constructor
Definition jacobi.h:117
void set_block(bool b)
sets if blocked jacobi is used (inverting block-diagonal), or plain (scalar) diagonal if false
Definition jacobi.h:145
bool CheckDiagonalInvertible(const TSparseMatrix &A)
Definition sparsematrix_util.h:1094
#define UG_LOG(msg)
Definition log.h:367
double number
Definition types.h:124
CPUAlgebra::matrix_type matrix_type
CPUAlgebra::vector_type vector_type
#define PROFILE_BEGIN_GROUP(name, groups)
Definition profiler.h:255
SmartPtr< T, FreePolicy > make_sp(T *inst)
returns a SmartPtr for the passed raw pointer
Definition smart_pointer.h:836