34 #ifndef __H__UG__PLUGINS__NAVIER_STOKES__INCOMPRESSIBLE__FVCR__PCRILUT__
35 #define __H__UG__PLUGINS__NAVIER_STOKES__INCOMPRESSIBLE__FVCR__PCRILUT__
58 template <
typename TAlgebra>
122 newInst->set_damp(this->
damping());
124 newInst->set_info(this->
m_info);
174 virtual const char*
name()
const {
return "PCRILUTPreconditioner";}
179 STATIC_ASSERT(matrix_type::rows_sorted, Matrix_has_to_have_sorted_rows);
192 std::vector<IndexLayout::Element> vIndex;
199 typedef typename matrix_type::connection connection;
200 m_L.resize_and_clear(
m_A.num_rows(),
m_A.num_cols());
201 m_U.resize_and_clear(
m_A.num_rows(),
m_A.num_cols());
204 std::vector<typename matrix_type::connection> con;
208 static const size_t velocity=0;
209 static const size_t pressure=1;
212 for(
typename matrix_type::row_iterator i_it =
m_A.begin_row(0); i_it !=
m_A.end_row(0); ++i_it)
213 con.push_back(connection(i_it.index(), i_it.value()));
214 m_U.set_matrix_row(0, &con[0], con.size());
216 size_t totalentries=0;
219 for(
size_t i=1; i<
m_A.num_rows(); i++)
222 if (
m_A(i,i)!=0) itype=velocity;
else itype=pressure;
231 for(
typename matrix_type::row_iterator i_it =
m_A.begin_row(i); i_it !=
m_A.end_row(i); ++i_it)
233 con.push_back(connection(i_it.index(), i_it.value()));
238 typename matrix_type::row_iterator i_it =
m_A.begin_row(i);
239 con.push_back(connection(i_it.index(), i_it.value()));
241 for(; i_it !=
m_A.end_row(i); ++i_it)
243 if (i_it.value()==0)
continue;
244 con.push_back(connection(i_it.index(), i_it.value()));
250 for(
size_t i_it = 0; i_it < con.size(); ++i_it)
252 size_t k = con[i_it].iIndex;
259 if(con[i_it].dValue == 0.0)
continue;
260 UG_ASSERT(
m_U.num_connections(k) != 0 &&
m_U.begin_row(k).index() == k,
"");
266 block_type &d = con[i_it].dValue = con[i_it].dValue / ukk;
268 typename matrix_type::row_iterator k_it =
m_U.begin_row(k);
271 while(k_it !=
m_U.end_row(k) && j < con.size())
275 if(k_it.index() == con[j].iIndex)
278 con[j].dValue -= k_it.value() * d;
281 else if(k_it.index() < con[j].iIndex)
286 typename matrix_type::connection
287 c(k_it.index(), k_it.value() * d * -1.0);
292 if (
m_A(k_it.index(),k_it.index())!=0) kitType=velocity;
else kitType=pressure;
294 if ((itype==velocity)&&(kitType==velocity)){
298 con.insert(con.begin()+j, c);
302 if ((itype==velocity)&&(kitType==pressure)){
306 con.insert(con.begin()+j, c);
310 if ((itype==pressure)&&(kitType==velocity)){
314 con.insert(con.begin()+j, c);
318 if ((itype==pressure)&&(kitType==pressure)){
322 con.insert(con.begin()+j, c);
337 if (k_it!=
m_U.end_row(k)){
338 for (;k_it!=
m_U.end_row(k);++k_it){
339 typename matrix_type::connection c(k_it.index(),-k_it.value()*d);
342 if (
m_A(k_it.index(),k_it.index())!=0) kitType=velocity;
else kitType=pressure;
343 if ((itype==velocity)&&(kitType==velocity)){
349 if ((itype==velocity)&&(kitType==pressure)){
355 if ((itype==pressure)&&(kitType==velocity)){
376 totalentries+=con.size();
377 if(maxentries < con.size()) maxentries = con.size();
380 m_L.set_matrix_row(i, &con[0], u_part);
381 m_U.set_matrix_row(i, &con[u_part], con.size()-u_part);
390 UG_LOG(
"\n PCRILUT storage information:\n");
391 UG_LOG(
" A nr of connections: " <<
m_A.total_num_connections() <<
"\n");
392 UG_LOG(
" L+U nr of connections: " <<
m_L.total_num_connections()+
m_U.total_num_connections() <<
"\n");
393 UG_LOG(
" Increase factor: " << (
float)(
m_L.total_num_connections() +
m_U.total_num_connections() )/
m_A.total_num_connections() <<
"\n");
407 for(
size_t i=0; i <
m_L.num_rows(); i++)
411 for(
typename matrix_type::row_iterator it =
m_L.begin_row(i); it !=
m_L.end_row(i); ++it)
412 MatMultAdd(c[i], 1.0, c[i], -1.0, it.value(), c[it.index()] );
420 size_t i=
m_U.num_rows()-1;
421 typename matrix_type::row_iterator it =
m_U.begin_row(i);
430 UG_LOG(
"Warning: zero entry in last row of U with corresponding non-zero rhs entry (" <<
BlockNorm(
s) <<
")\n");
438 for(
size_t i=
m_U.num_rows()-2; ; i--)
440 typename matrix_type::row_iterator it =
m_U.begin_row(i);
446 for(; it !=
m_U.end_row(i); ++it)
448 MatMultAdd(
s, 1.0,
s, -1.0, it.value(), c[it.index()] );
461 for(
size_t i=0; i <
m_L.num_rows(); i++)
465 for(
typename matrix_type::row_iterator it =
m_L.begin_row(i); it !=
m_L.end_row(i); ++it)
466 MatMultAdd(c[i], 1.0, c[i], -1.0, it.value(), c[it.index()] );
474 size_t i=
m_U.num_rows()-1;
475 typename matrix_type::row_iterator it =
m_U.begin_row(i);
484 UG_LOG(
"Warning: zero entry in last row of U with corresponding non-zero rhs entry (" <<
BlockNorm(
s) <<
")\n");
492 for(
size_t i=
m_U.num_rows()-2; ; i--)
494 typename matrix_type::row_iterator it =
m_U.begin_row(i);
500 for(; it !=
m_U.end_row(i); ++it)
502 MatMultAdd(
s, 1.0,
s, -1.0, it.value(), c[it.index()] );
531 template <
typename TAlgebra>
533 template <
typename TAlgebra>
TAlgebra::vector_type vector_type
SmartPtr< IDebugWriter< algebra_type > > debug_writer()
TAlgebra::matrix_type matrix_type
SmartPtr< IDamping< TAlgebra::vector_type, TAlgebra::vector_type > > damping()
PCRILUTPreconditioner class.
Definition: pcr_ilut.h:60
void set_threshold(number threshvv, number threshvp, number threshpv, number threshpp)
sets threshold for incomplete LU factorisation
Definition: pcr_ilut.h:137
TAlgebra::vector_type vector_type
Definition: pcr_ilut.h:66
bool m_info
Definition: pcr_ilut.h:525
number m_eps
Definition: pcr_ilut.h:520
number m_eps_vv
Definition: pcr_ilut.h:521
bool m_remove_zeros
Definition: pcr_ilut.h:526
PCRILUTPreconditioner(number threshvv, number threshvp, number threshpv, number threshpp, bool info)
Definition: pcr_ilut.h:110
virtual bool preprocess(SmartPtr< MatrixOperator< matrix_type, vector_type > > pOp)
Definition: pcr_ilut.h:177
PCRILUTPreconditioner(number threshvv, number thresh_vp_pv_pp)
Definition: pcr_ilut.h:87
virtual bool postprocess()
Definition: pcr_ilut.h:514
PCRILUTPreconditioner(double eps=1e-6)
Definition: pcr_ilut.h:81
number m_eps_pp
Definition: pcr_ilut.h:524
IPreconditioner< TAlgebra >::matrix_operator_type matrix_operator_type
Matrix Operator type.
Definition: pcr_ilut.h:72
void set_threshold(number threshvv, number thresh_vp_pv_pp)
sets threshold for incomplete LU factorisation
Definition: pcr_ilut.h:147
PCRILUTPreconditioner(double eps, bool info)
Definition: pcr_ilut.h:99
TAlgebra::matrix_type matrix_type
Definition: pcr_ilut.h:69
matrix_type m_U
Definition: pcr_ilut.h:518
matrix_type m_A
Definition: pcr_ilut.h:519
void set_info(bool info)
sets storage information output to true or false
Definition: pcr_ilut.h:166
static const number m_small_lower
Definition: pcr_ilut.h:527
matrix_type m_L
Definition: pcr_ilut.h:517
number m_eps_pv
Definition: pcr_ilut.h:523
number m_eps_vp
Definition: pcr_ilut.h:522
void set_threshold(number thresh)
Definition: pcr_ilut.h:156
TAlgebra algebra_type
Definition: pcr_ilut.h:63
virtual const char * name() const
Definition: pcr_ilut.h:174
SmartPtr< ILinearIterator< vector_type > > clone()
Definition: pcr_ilut.h:118
PCRILUTPreconditioner(number threshvv, number threshvp, number threshpv, number threshpp)
Definition: pcr_ilut.h:93
matrix_type::value_type block_type
Definition: pcr_ilut.h:75
virtual bool supports_parallel() const
returns if parallel solving is supported
Definition: pcr_ilut.h:134
virtual bool step(SmartPtr< MatrixOperator< matrix_type, vector_type > > pOp, vector_type &c, const vector_type &d)
Definition: pcr_ilut.h:400
virtual ~PCRILUTPreconditioner()
Definition: pcr_ilut.h:129
static const number m_small_upper
Definition: pcr_ilut.h:528
PCRILUTPreconditioner(number threshvv, number thresh_vp_pv_pp, bool info)
Definition: pcr_ilut.h:104
TVector::value_type value_type
bool change_storage_type(ParallelStorageType type)
ParallelMatrix< SparseMatrix< double > > matrix_type
ParallelVector< Vector< double > > vector_type
void MatAddSlaveRowsToMasterRowOverlap0(TMatrix &mat)
void CollectUniqueElements(std::vector< typename TLayout::Element > &elemsOut, const TLayout &layout)
bool InverseMatMult(DenseVector< FixedArray1< double, 1 > > &dest, double beta1, const DenseMatrix< FixedArray2< double, 1, 1 > > &A1, const DenseVector< FixedArray1< double, 1 > > &w1)
void MatMultAdd(DenseVector< vector_t > &dest, const number &alpha1, const DenseVector< vector_t > &v1, const number &beta1, const DenseMatrix< matrix_t > &A1, const DenseVector< vector_t > &w1)
#define UG_ASSERT(expr, msg)
#define STATIC_ASSERT(expr, msg)
void SetDirichletRow(TMatrix &mat, const DoFIndex &ind)
number BlockNorm(const number &a)