47 #include "Moocho_ConfigDefs.hpp"
50 #ifdef HAVE_MOOCHO_MA28
53 #include "AbstractLinAlgPack_DirectSparseSolverMA28.hpp"
54 #include "AbstractLinAlgPack_MatrixScaling_Strategy.hpp"
55 #include "AbstractLinAlgPack_VectorDenseEncap.hpp"
56 #include "DenseLinAlgPack_PermVecMat.hpp"
57 #include "Teuchos_AbstractFactoryStd.hpp"
58 #include "Teuchos_Assert.hpp"
59 #include "Teuchos_Workspace.hpp"
60 #include "Teuchos_dyn_cast.hpp"
61 #include "FortranTypes_f_open_file.hpp"
67 T my_max(
const T& v1,
const T& v2 ) {
return v1 > v2 ? v1 : v2; }
71 std::valarray<T>& cva(
const std::valarray<T>& va )
73 return const_cast<std::valarray<T>&
>(va);
77 namespace AbstractLinAlgPack {
109 DirectSparseSolverMA28::FactorizationStructureMA28::FactorizationStructureMA28()
110 :m_(0),n_(0),max_n_(0),nz_(0),licn_(0),lirn_(0)
111 ,u_(0.1),scaling_(NO_SCALING)
135 const FactorizationStructureMA28
137 const FactorizationNonzerosMA28
138 &fn =
dyn_cast<
const FactorizationNonzerosMA28>(*this->get_fact_nonzeros());
143 y == NULL, std::invalid_argument
144 ,
"DirectSparseSolverMA28::BasisMatrixMA28::V_InvMtV(...) : Error! " );
146 const size_type y_dim = y->dim(), x_dim = x.dim();
149 fs.rank_ != y_dim, std::invalid_argument
150 ,
"DirectSparseSolverMA28::BasisMatrixMA28::V_InvMtV(...) : Error! " );
152 fs.rank_ != x_dim, std::invalid_argument
153 ,
"DirectSparseSolverMA28::BasisMatrixMA28::V_InvMtV(...) : Error! " );
156 VectorDenseMutableEncap yd(*y);
157 VectorDenseEncap xd(x);
160 Workspace<value_type> xfull_s(wss,fs.max_n_,
false);
161 DVectorSlice xfull(&xfull_s[0],xfull_s.size());
162 Workspace<value_type> ws(wss,fs.max_n_,
false);
163 DVectorSlice w(&ws[0],ws.size());
176 for( k = 1; k <= x_dim; ++k )
177 xfull(row_perm(k)) = xd()(k);
180 if( fs.matrix_scaling_.get() )
181 fs.matrix_scaling_->scale_rhs( M_trans, xfull.raw_ptr() );
186 fs.max_n_, &cva(fn.a_)[0], fs.licn_, &cva(fs.icn_)[0], &cva(fs.ikeep_)[0]
187 ,xfull.raw_ptr(), w.raw_ptr(), mtype );
190 if( fs.matrix_scaling_.get() )
191 fs.matrix_scaling_->scale_rhs( M_trans, xfull.raw_ptr() );
198 for( k = 1; k <= y_dim; ++k )
199 yd()(k) = xfull(col_perm(k));
215 ,
bool print_ma28_outputs
216 ,
const std::string& output_file_name
218 :estimated_fillin_ratio_(estimated_fillin_ratio)
224 ,print_ma28_outputs_(print_ma28_outputs)
225 ,output_file_name_(output_file_name)
234 namespace mmp = MemMngPack;
261 ,FactorizationStructure *fact_struc
262 ,FactorizationNonzeros *fact_nonzeros
270 typedef MatrixConvertToSparse MCTS;
275 *out <<
"\nUsing MA28 to analyze and factor a new matrix ...\n";
278 FactorizationStructureMA28
279 &fs =
dyn_cast<FactorizationStructureMA28>(*fact_struc);
280 FactorizationNonzerosMA28
281 &fn =
dyn_cast<FactorizationNonzerosMA28>(*fact_nonzeros);
284 set_ma28_parameters(&fs);
290 nz = A.
num_nonzeros( MCTS::EXTRACT_FULL_MATRIX ,MCTS::ELEMENTS_ALLOW_DUPLICATES_SUM );
294 n <= 0 || m <= 0 || m > n, std::invalid_argument
295 ,
"DirectSparseSolverMA28::imp_analyze_and_factor(...) : Error!" );
298 fs.m_ = m; fs.n_ = n; fs.nz_ = nz;
299 fs.max_n_ = my_max(fs.m_,fs.n_);
302 if( estimated_fillin_ratio_ < 1.0 ) {
303 if( out ) *out <<
"Warning, client set estimated_fillin_ratio = " << estimated_fillin_ratio_
304 <<
" < 1.0.\nSetting estimated_fillin_ratio = 10.0 ...\n";
305 estimated_fillin_ratio_ = 10.0;
307 if( fs.licn_ < fs.nz_ ) fs.licn_ =
static_cast<index_type
>(estimated_fillin_ratio_ * fs.nz_);
308 if( fs.lirn_ < fs.nz_ ) fs.lirn_ =
static_cast<index_type
>(estimated_fillin_ratio_ * fs.nz_);
311 fs.ivect_.resize(fs.nz_);
312 fs.jvect_.resize(fs.nz_);
314 index_type iflag = 0;
315 for(
int num_fac = 0; num_fac < 5; ++num_fac ) {
318 fs.icn_.resize(fs.licn_);
319 fn.a_.resize(fs.licn_);
320 fs.ikeep_.resize( fs.ma28_.lblock() ? 5*fs.max_n_ : 4*fs.max_n_ + 1 );
321 Workspace<index_type> irn_tmp_(wss,fs.lirn_), iw(wss,8*fs.max_n_);
322 Workspace<value_type> w(wss,fs.max_n_);
326 MCTS::EXTRACT_FULL_MATRIX
327 ,MCTS::ELEMENTS_ALLOW_DUPLICATES_SUM
334 std::copy( &fs.ivect_[0], &fs.ivect_[0] + fs.nz_, &irn_tmp_[0] );
335 std::copy( &fs.jvect_[0], &fs.jvect_[0] + fs.nz_, &fs.icn_[0] );
338 if( fs.matrix_scaling_.get() )
339 fs.matrix_scaling_->scale_matrix(
340 fs.m_, fs.n_, fs.nz_, &fs.ivect_[0], &fs.jvect_[0],
true
346 *out <<
"\nCalling ma28ad(...) ...\n";
348 fs.max_n_, fs.nz_, &fn.a_[0], fs.licn_, &irn_tmp_[0], fs.lirn_, &fs.icn_[0], fs.u_
349 ,&fs.ikeep_[0], &iw[0], &w[0], &iflag
352 if(iflag != 0 && out)
353 *out <<
"\nMA28AD returned iflag = " << iflag <<
" != 0!\n";
356 print_ma28_outputs(
true,iflag,fs,&w[0],out);
358 if( iflag >= 0 )
break;
361 case LICN_AND_LIRN_TOO_SMALL:
363 case LICN_FAR_TOO_SMALL:
366 *out <<
"\nWarning! iflag = " << iflag <<
", LIRN and/or LICN are too small!\n"
367 <<
"Increasing lirn = " << fs.lirn_ <<
" amd licn = " << fs.licn_
368 <<
" by a factor of 10\n"
369 <<
"(try increasing estimated_fillin_ratio = " << estimated_fillin_ratio_
370 <<
" to a larger value next time)...\n";
371 fs.lirn_ = 10 * fs.lirn_;
372 fs.licn_ = 10 * fs.licn_;
378 ThrowIFlagException(iflag);
382 *rank = fs.ma28_.irank();
384 row_perm->resize(fs.m_);
385 if( *rank < fs.m_ ) {
387 *row_perm_ikeep = &fs.ikeep_[fs.max_n_],
388 *row_perm_itr = &(*row_perm)(1),
389 *row_perm_last = row_perm_itr + fs.m_;
390 for(; row_perm_itr != row_perm_last;)
391 *row_perm_itr++ = abs(*row_perm_ikeep++);
393 std::sort(&(*row_perm)[0] , &(*row_perm)[0] + (*rank) );
394 std::sort(&(*row_perm)[0] + (*rank) , &(*row_perm)[0] + m );
397 DenseLinAlgPack::identity_perm( row_perm );
400 col_perm->resize(fs.n_);
401 if( *rank < fs.n_ ) {
403 *col_perm_ikeep = &fs.ikeep_[2*fs.max_n_],
404 *col_perm_itr = &(*col_perm)(1),
405 *col_perm_last = col_perm_itr + fs.n_;
406 for(; col_perm_itr != col_perm_last;)
407 *col_perm_itr++ = abs(*col_perm_ikeep++);
409 std::sort(&(*col_perm)[0] , &(*col_perm)[0] + (*rank) );
410 std::sort(&(*col_perm)[0] + (*rank) , &(*col_perm)[0] + n );
413 DenseLinAlgPack::identity_perm( col_perm );
417 fs.row_perm_ = *row_perm;
418 fs.col_perm_ = *col_perm;
425 ,
const FactorizationStructure &fact_struc
426 ,FactorizationNonzeros *fact_nonzeros
431 typedef MatrixConvertToSparse MCTS;
436 *out <<
"\nUsing MA28 to factor a new matrix ...\n";
439 const FactorizationStructureMA28
440 &fs =
dyn_cast<
const FactorizationStructureMA28>(fact_struc);
441 FactorizationNonzerosMA28
442 &fn =
dyn_cast<FactorizationNonzerosMA28>(*fact_nonzeros);
448 nz = A.
num_nonzeros( MCTS::EXTRACT_FULL_MATRIX ,MCTS::ELEMENTS_ALLOW_DUPLICATES_SUM );
453 m != fs.m_ || n != fs.n_ || fs.nz_ != nz, std::invalid_argument
454 ,
"DirectSparseSolverMA28::imp_factor(...) : Error, "
455 "A is not compatible with matrix passed to imp_analyze_and_factor()!" );
459 if(fn.a_.size() < fs.licn_) fn.a_.resize(fs.licn_);
460 Workspace<index_type> iw(wss,5*fs.max_n_);
461 Workspace<value_type> w(wss,fs.max_n_);
465 MCTS::EXTRACT_FULL_MATRIX
466 ,MCTS::ELEMENTS_ALLOW_DUPLICATES_SUM
475 if( fs.matrix_scaling_.get() )
476 fs.matrix_scaling_->scale_matrix(
477 fs.m_, fs.n_, fs.nz_, &cva(fs.ivect_)[0], &cva(fs.jvect_)[0],
false
482 index_type iflag = 0;
484 fs.max_n_, fs.nz_, &fn.a_[0], fs.licn_, &cva(fs.ivect_)[0], &cva(fs.jvect_)[0], &cva(fs.icn_)[0]
485 ,&cva(fs.ikeep_)[0], &iw[0], &w[0], &iflag
488 if(iflag != 0 && out)
489 *out <<
"\nMA28BD returned iflag = " << iflag <<
" != 0!\n";
492 print_ma28_outputs(
false,iflag,fs,&w[0],out);
495 ThrowIFlagException(iflag);
501 void DirectSparseSolverMA28::set_ma28_parameters( FactorizationStructureMA28* fs )
504 fs->ma28_.lblock( FortranTypes::F_FALSE );
506 fs->ma28_.grow( grow_ ? FortranTypes::F_TRUE : FortranTypes::F_FALSE );
508 fs->ma28_.nsrch(nsrch_);
509 fs->ma28_.lbig( lbig_ ? FortranTypes::F_TRUE : FortranTypes::F_FALSE );
511 if( output_file_name_.length() > 0 && fs->ma28_.lp() == 0 ) {
513 index_type iout = 25;
514 FortranTypes::f_open_file( iout, output_file_name_.c_str() );
518 else if( output_file_name_.length() == 0 && fs->ma28_.lp() ) {
524 void DirectSparseSolverMA28::print_ma28_outputs(
527 ,
const FactorizationStructureMA28 &fs
528 ,
const value_type w[]
532 if( print_ma28_outputs_ && out ) {
533 *out <<
"\nReturn parameters from MA28 (call number = " << ++file_output_num_ <<
")\n";
534 if( fs.ma28_.grow() == FortranTypes::F_TRUE )
535 *out <<
"w(1) = " << w[0] << std::endl;
536 *out <<
"rmin = " << fs.ma28_.rmin() << std::endl;
537 *out <<
"irncp = " << fs.ma28_.irncp() << std::endl;
538 *out <<
"icncp = " << fs.ma28_.icncp() << std::endl;
539 *out <<
"minirn = " << fs.ma28_.minirn() << std::endl;
540 *out <<
"minicn = " << fs.ma28_.minicn() << std::endl;
541 *out <<
"irank = " << fs.ma28_.irank() << std::endl;
542 *out <<
"themax = " << fs.ma28_.themax() << std::endl;
543 if( fs.ma28_.lbig() == FortranTypes::F_TRUE )
544 *out <<
"big = " << fs.ma28_.big() << std::endl;
545 *out <<
"ndrop = " << fs.ma28_.ndrop() << std::endl;
547 *out <<
"\nAnalysis:\n"
548 <<
"estimated_fillin_ratio can be reduced to max(minirn,minicn)/nz = "
549 <<
"max(" << fs.ma28_.minirn() <<
"," << fs.ma28_.minicn() <<
")/" << fs.nz_
550 <<
" = " << my_max( fs.ma28_.minirn(), fs.ma28_.minicn() ) / (
double)fs.nz_
556 void DirectSparseSolverMA28::ThrowIFlagException(index_type iflag)
558 E_IFlag e_iflag =
static_cast<E_IFlag
>(iflag);
559 const char msg_err_head[] =
"DirectSparseSolverMA28::ThrowIFlagException(iflag) : Error";
561 case SLOW_ITER_CONV :
563 true, std::runtime_error
564 ,msg_err_head <<
", Convergence to slow" );
567 true, std::runtime_error
568 ,msg_err_head <<
", Maximum iterations exceeded");
569 case MA28BD_CALLED_WITH_DROPPED :
571 true, std::logic_error
572 ,msg_err_head <<
", ma28bd called with elements dropped in ma28ad");
573 case DUPLICATE_ELEMENTS :
575 true, FactorizationFailure
576 ,msg_err_head <<
", Duplicate elements have been detected");
577 case NEW_NONZERO_ELEMENT :
579 true, FactorizationFailure
580 ,msg_err_head <<
", A new non-zero element has be passed to ma28bd that was not ot ma28ad");
581 case N_OUT_OF_RANGE :
583 true, FactorizationFailure
584 ,msg_err_head <<
", 1 <=max(n,m) <= 32767 has been violated");
587 true, std::logic_error
588 ,msg_err_head <<
", nz <= 0 has been violated");
591 true, std::logic_error
592 ,msg_err_head <<
", licn <= nz has been violated");
595 true, std::logic_error
596 ,msg_err_head <<
", lirn <= nz has been violated");
597 case ERROR_DURRING_BLOCK_TRI :
599 true, FactorizationFailure
600 ,msg_err_head <<
", An error has occured durring block triangularization");
601 case LICN_AND_LIRN_TOO_SMALL :
603 true, FactorizationFailure
604 ,msg_err_head <<
", licn and lirn are to small to hold matrix factorization");
605 case LICN_TOO_SMALL :
607 true, FactorizationFailure
608 ,msg_err_head <<
", licn is to small to hold matrix factorization");
609 case LICN_FAR_TOO_SMALL :
611 true, FactorizationFailure
612 ,msg_err_head <<
", licn is to far small to hold matrix factorization");
613 case LIRN_TOO_SMALL :
615 true, FactorizationFailure
616 ,msg_err_head <<
", lirn is to small to hold matrix factorization");
617 case NUMERICALLY_SINGULAR :
619 true, FactorizationFailure
620 ,msg_err_head <<
", matrix is numerically singular, see \'abort2\'");
621 case STRUCTURALLY_SINGULAR :
623 true, FactorizationFailure
624 ,msg_err_head <<
", matrix is structurally singular, see \'abort1\'");
633 #endif // HAVE_MOOCHO_MA28
DirectSparseSolverMA28(value_type estimated_fillin_ratio=10.0, value_type u=0.1, bool grow=false, value_type tol=0.0, index_type nsrch=4, bool lbig=false, bool print_ma28_outputs=false, const std::string &output_file_name="")
Default constructor.
Teuchos::RCP< const Teuchos::AbstractFactory< BasisMatrix > > basis_matrix_factory_ptr_t
virtual void coor_extract_nonzeros(EExtractRegion extract_region, EElementUniqueness element_uniqueness, const index_type len_Aval, value_type Aval[], const index_type len_Aij, index_type Arow[], index_type Acol[], const index_type row_offset=0, const index_type col_offset=0) const =0
Extract sparse elements in a coordinate data structure.
const Teuchos::RCP< FactorizationNonzeros > create_fact_nonzeros() const
const BasisMatrix::fact_struc_ptr_t & get_fact_struc() const
virtual index_type num_nonzeros(EExtractRegion extract_region, EElementUniqueness element_uniqueness) const =0
Returns the number of nonzeros in the matrix.
Teuchos::RCP< BasisMatrixImp > create_matrix() const
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
T_To & dyn_cast(T_From &from)
const Teuchos::RCP< FactorizationStructure > create_fact_struc() const
void V_InvMtV(VectorMutable *v_lhs, BLAS_Cpp::Transp trans_rhs1, const Vector &v_rhs2) const
const basis_matrix_factory_ptr_t basis_matrix_factory() const
virtual void estimated_fillin_ratio(value_type estimated_fillin_ratio)=0
Set the estimate of the fill-in ratio of the factorization.
virtual size_type cols() const
Return the number of columns in the matrix.
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
void imp_factor(const AbstractLinAlgPack::MatrixConvertToSparse &A, const FactorizationStructure &fact_struc, FactorizationNonzeros *fact_nonzeros, std::ostream *out)
Mix-in interface for extracing explicit elements from a sparse matrix in one of several Fortran compa...
void estimated_fillin_ratio(value_type estimated_fillin_ratio)
virtual size_type rows() const
Return the number of rows in the matrix.
TEUCHOSCORE_LIB_DLL_EXPORT Teuchos::RCP< WorkspaceStore > get_default_workspace_store()
void imp_analyze_and_factor(const AbstractLinAlgPack::MatrixConvertToSparse &A, FactorizationStructure *fact_struc, FactorizationNonzeros *fact_nonzeros, DenseLinAlgPack::IVector *row_perm, DenseLinAlgPack::IVector *col_perm, size_type *rank, std::ostream *out)