42 #ifdef SPARSE_SOLVER_PACK_USE_SUPERLU
61 int local_panel_size = 0;
64 class StaticSuperLUInit {
68 local_panel_size = sp_ienv(1);
69 local_relax = sp_ienv(2);
70 StatInit(local_panel_size,local_relax);
78 StaticSuperLUInit static_super_lu_init;
87 std::valarray<T>& cva(
const std::valarray<T>& va )
89 return const_cast<std::valarray<T>&
>(va);
94 namespace SuperLUPack {
96 class SuperLUSolverImpl;
102 class SuperLUSolverImpl :
public SuperLUSolver {
109 class FactorizationStructureImpl :
public FactorizationStructure {
111 friend class SuperLUSolverImpl;
115 std::valarray<int> perm_r_;
116 std::valarray<int> perm_c_;
117 std::valarray<int> etree_;
121 std::valarray<int> perm_r_orig_;
122 std::valarray<int> perm_c_orig_;
126 class FactorizationNonzerosImpl :
public FactorizationNonzeros {
128 friend class SuperLUSolverImpl;
140 void analyze_and_factor(
144 ,
const double a_val[]
146 ,
const int a_col_ptr[]
147 ,FactorizationStructure *fact_struct
148 ,FactorizationNonzeros *fact_nonzeros
158 ,
const double a_val[]
160 ,
const int a_col_ptr[]
161 ,
const FactorizationStructure &fact_struct
162 ,FactorizationNonzeros *fact_nonzeros
166 const FactorizationStructure &fact_struct
167 ,
const FactorizationNonzeros &fact_nonzeros
180 void copy_basis_nonzeros(
184 ,
const double a_orig_val[]
185 ,
const int a_orig_row_i[]
186 ,
const int a_orig_col_ptr[]
187 ,
const int a_orig_perm_r[]
188 ,
const int a_orig_perm_c[]
203 SuperLUSolver::create_solver()
209 SuperLUSolver::create_fact_struct()
211 return Teuchos::rcp(
new SuperLUSolverImpl::FactorizationStructureImpl());
215 SuperLUSolver::create_fact_nonzeros()
217 return Teuchos::rcp(
new SuperLUSolverImpl::FactorizationNonzerosImpl());
226 void SuperLUSolverImpl::analyze_and_factor(
230 ,
const double a_val[]
232 ,
const int a_col_ptr[]
233 ,FactorizationStructure *fact_struct
234 ,FactorizationNonzeros *fact_nonzeros
244 FactorizationStructureImpl
245 &fs =
dyn_cast<FactorizationStructureImpl>(*fact_struct);
246 FactorizationNonzerosImpl
247 &fn =
dyn_cast<FactorizationNonzerosImpl>(*fact_nonzeros);
257 fs.perm_r_.resize(m);
258 fs.perm_c_.resize(n);
263 dCreate_CompCol_Matrix(
265 ,const_cast<double*>(a_val)
266 ,const_cast<int*>(a_row_i)
267 ,const_cast<int*>(a_col_ptr)
273 get_perm_c(permc_spec, &A, &fs.perm_c_[0]);
277 sp_preorder(refact,&A,&fs.perm_c_[0],&fs.etree_[0],&AC);
298 info != 0, std::runtime_error
299 ,
"SuperLUSolverImpl::analyze_and_factor(...): Error, dgstrf(...) returned info = " << info
302 std::copy( &fs.perm_r_[0], &fs.perm_r_[0] + m, perm_r );
303 std::copy( &fs.perm_c_[0], &fs.perm_c_[0] + n, perm_c );
314 fs.perm_r_orig_ = fs.perm_r_;
315 fs.perm_c_orig_ = fs.perm_c_;
317 Workspace<double> b_val(wss,nz);
318 Workspace<int> b_row_i(wss,nz);
319 Workspace<int> b_col_ptr(wss,n+1);
321 m,n,nz,a_val,a_row_i,a_col_ptr
322 ,&fs.perm_r_orig_[0],&fs.perm_c_orig_[0],fs.rank_
323 ,&b_val[0],&b_row_i[0],&b_col_ptr[0]
329 fs.rank_, fs.rank_, fs.nz_, &b_val[0], &b_row_i[0], &b_col_ptr[0]
330 ,fact_struct, fact_nonzeros
331 ,&fs.perm_r_[0], &fs.perm_c_[0], &b_rank
334 (b_rank != *rank), std::runtime_error
335 ,
"SuperLUSolverImpl::analyze_and_factor(...): Error, the rank determined by "
336 "the factorization of the rectangular " << m <<
" x " << n <<
" matrix of "
337 << (*rank) <<
" is not the same as the refactorization of the basis returned as "
348 void SuperLUSolverImpl::factor(
352 ,
const double a_val[]
354 ,
const int a_col_ptr[]
355 ,
const FactorizationStructure &fact_struct
356 ,FactorizationNonzeros *fact_nonzeros
363 const FactorizationStructureImpl
364 &fs =
dyn_cast<
const FactorizationStructureImpl>(fact_struct);
365 FactorizationNonzerosImpl
366 &fn =
dyn_cast<FactorizationNonzerosImpl>(*fact_nonzeros);
371 Workspace<double> b_val(wss,fs.nz_);
372 Workspace<int> b_row_i(wss,fs.nz_);
373 Workspace<int> b_col_ptr(wss,fs.rank_+1);
374 if(fs.m_orig_ > fs.n_orig_) {
377 m,n,nz,a_val,a_row_i,a_col_ptr
378 ,&cva(fs.perm_r_orig_)[0],&cva(fs.perm_c_orig_)[0],fs.rank_
379 ,&b_val[0],&b_row_i[0],&b_col_ptr[0]
383 (b_nz != fs.nz_), std::runtime_error
384 ,
"SuperLUSolverImpl::factor(...): Error!"
388 std::copy( a_val, a_val + nz, &b_val[0] );
389 std::copy( a_row_i, a_row_i + nz, &b_row_i[0] );
390 std::copy( a_col_ptr, a_col_ptr + n+1, &b_col_ptr[0] );
395 dCreate_CompCol_Matrix(
396 &A, fs.rank_, fs.rank_, fs.nz_
420 ,const_cast<int*>(&cva(fs.etree_)[0])
431 info != 0, std::runtime_error
432 ,
"SuperLUSolverImpl::factor(...): Error, dgstrf(...) returned info = " << info
437 void SuperLUSolverImpl::solve(
438 const FactorizationStructure &fact_struct
439 ,
const FactorizationNonzeros &fact_nonzeros
450 const FactorizationStructureImpl
451 &fs =
dyn_cast<
const FactorizationStructureImpl>(fact_struct);
452 const FactorizationNonzerosImpl
453 &fn =
dyn_cast<
const FactorizationNonzerosImpl>(fact_nonzeros);
456 n != fs.rank_, std::runtime_error
457 ,
"SuperLUSolverImpl::solve(...): Error, the dimmensions n = " << n <<
" and fs.rank = " << fs.rank_
458 <<
" do not match up!"
462 dCreate_Dense_Matrix(&B, n, nrhs, rhs, ldrhs, DN, D_, GE);
465 transc[0] = ( transp ?
'T' :
'N' );
470 ,const_cast<SuperMatrix*>(&fn.L_)
471 ,const_cast<SuperMatrix*>(&fn.U_)
478 info != 0, std::runtime_error
479 ,
"SuperLUSolverImpl::solve(...): Error, dgssv(...) returned info = " << info
486 void SuperLUSolverImpl::copy_basis_nonzeros(
490 ,
const double a_orig_val[]
491 ,
const int a_orig_row_i[]
492 ,
const int a_orig_col_ptr[]
493 ,
const int a_orig_perm_r[]
494 ,
const int a_orig_perm_c[]
503 b_col_ptr[0] = *b_nz;
504 for(
int j = 0; j < rank; ++j ) {
505 const int col_start_k = a_orig_col_ptr[j];
506 const int col_end_k = a_orig_col_ptr[j+1];
507 for(
int k = col_start_k; k < col_end_k; ++k ) {
508 const int i_orig = a_orig_row_i[k];
510 b_val[*b_nz] = a_orig_val[k];
511 b_row_i[*b_nz] = a_orig_row_i[k];
515 b_col_ptr[j+1] = *b_nz;
521 #endif // SPARSE_SOLVER_PACK_USE_SUPERLU
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
const LAPACK_C_Decl::f_int LAPACK_C_Decl::f_dbl_prec * A
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
void copy(const f_int &N, const f_dbl_prec *X, const f_int &INCX, f_dbl_prec *Y, const f_int &INCY)
T_To & dyn_cast(T_From &from)
const LAPACK_C_Decl::f_int const LAPACK_C_Decl::f_int const LAPACK_C_Decl::f_dbl_prec const LAPACK_C_Decl::f_int const LAPACK_C_Decl::f_int LAPACK_C_Decl::f_dbl_prec B[]
const f_dbl_prec const f_int const f_int const f_int f_dbl_prec rhs[]
TEUCHOSCORE_LIB_DLL_EXPORT Teuchos::RCP< WorkspaceStore > get_default_workspace_store()