MBDyn-1.7.3
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups
LinSol Class Reference

#include <linsol.h>

Collaboration diagram for LinSol:

Classes

struct  solver_t
 

Public Types

enum  SolverType {
  EMPTY_SOLVER = 0, HARWELL_SOLVER, LAPACK_SOLVER, MESCHACH_SOLVER,
  NAIVE_SOLVER, SUPERLU_SOLVER, TAUCS_SOLVER, UMFPACK_SOLVER,
  KLU_SOLVER, Y12_SOLVER, LAST_SOLVER
}
 
enum  SolverFlags {
  SOLVER_FLAGS_NONE = 0x00U, SOLVER_FLAGS_ALLOWS_MAP = 0x01U, SOLVER_FLAGS_ALLOWS_CC = 0x02U, SOLVER_FLAGS_ALLOWS_DIR = 0x04U,
  SOLVER_FLAGS_TYPE_MASK = SOLVER_FLAGS_ALLOWS_MAP|SOLVER_FLAGS_ALLOWS_CC|SOLVER_FLAGS_ALLOWS_DIR, SOLVER_FLAGS_ALLOWS_MT_FCT = 0x10U, SOLVER_FLAGS_ALLOWS_MT_ASS = 0x20U, SOLVER_FLAGS_ALLOWS_COLAMD = 0x40U,
  SOLVER_FLAGS_ALLOWS_MMDATA = 0x80U, SOLVER_FLAGS_ALLOWS_REVERSE_CUTHILL_MC_KEE = 0x100U, SOLVER_FLAGS_ALLOWS_KING = 0x200U, SOLVER_FLAGS_ALLOWS_NESTED_DISSECTION = 0x400U,
  SOLVER_FLAGS_ALLOWS_MDAPLUSAT = 0x800U, SOLVER_FLAGS_ALLOWS_SLOAN = 0x1000U, SOLVER_FLAGS_PERM_MASK
}
 

Public Member Functions

 LinSol (void)
 
virtual ~LinSol (void)
 
SolverType GetSolver (void) const
 
const char *const GetSolverName (void) const
 
unsigned GetSolverFlags (void) const
 
unsigned GetNumThreads (void) const
 
integer iGetWorkSpaceSize (void) const
 
const doublerealdGetPivotFactor (void) const
 
const doublerealdGetDropTolerance (void) const
 
unsigned GetBlockSize (void) const
 
const SolutionManager::ScaleOptGetScale (void) const
 
integer GetMaxIterations (void) const
 
const char *const GetSolverName (SolverType t) const
 
unsigned GetSolverFlags (SolverType t) const
 
bool SetSolver (SolverType t, unsigned f=SOLVER_FLAGS_NONE)
 
bool SetSolverFlags (unsigned f)
 
bool AddSolverFlags (unsigned f)
 
bool MaskSolverFlags (unsigned f)
 
bool SetNumThreads (unsigned nt)
 
bool SetWorkSpaceSize (integer)
 
bool SetPivotFactor (const doublereal &d)
 
bool SetDropTolerance (const doublereal &d)
 
bool SetBlockSize (unsigned bs)
 
bool SetScale (const SolutionManager::ScaleOpt &scale)
 
bool SetMaxIterations (integer iMaxIter)
 
SolutionManager *const GetSolutionManager (integer iNLD, integer iLWS=0) const
 

Static Public Attributes

static SolverType defaultSolver
 

Protected Attributes

SolverType currSolver
 
unsigned solverFlags
 
unsigned nThreads
 
integer iWorkSpaceSize
 
unsigned blockSize
 
doublereal dPivotFactor
 
doublereal dDropTolerance
 
SolutionManager::ScaleOpt scale
 
integer iMaxIter
 

Detailed Description

Definition at line 39 of file linsol.h.

Member Enumeration Documentation

Enumerator
SOLVER_FLAGS_NONE 
SOLVER_FLAGS_ALLOWS_MAP 
SOLVER_FLAGS_ALLOWS_CC 
SOLVER_FLAGS_ALLOWS_DIR 
SOLVER_FLAGS_TYPE_MASK 
SOLVER_FLAGS_ALLOWS_MT_FCT 
SOLVER_FLAGS_ALLOWS_MT_ASS 
SOLVER_FLAGS_ALLOWS_COLAMD 
SOLVER_FLAGS_ALLOWS_MMDATA 
SOLVER_FLAGS_ALLOWS_REVERSE_CUTHILL_MC_KEE 
SOLVER_FLAGS_ALLOWS_KING 
SOLVER_FLAGS_ALLOWS_NESTED_DISSECTION 
SOLVER_FLAGS_ALLOWS_MDAPLUSAT 
SOLVER_FLAGS_ALLOWS_SLOAN 
SOLVER_FLAGS_PERM_MASK 

Definition at line 56 of file linsol.h.

56  {
57  SOLVER_FLAGS_NONE = 0x00U,
59  SOLVER_FLAGS_ALLOWS_CC = 0x02U,
61  SOLVER_FLAGS_TYPE_MASK = SOLVER_FLAGS_ALLOWS_MAP|SOLVER_FLAGS_ALLOWS_CC|SOLVER_FLAGS_ALLOWS_DIR,
64  //permutations
68  SOLVER_FLAGS_ALLOWS_KING = 0x200U,
71  SOLVER_FLAGS_ALLOWS_SLOAN = 0x1000U,
72 
74  SOLVER_FLAGS_ALLOWS_COLAMD |
75  SOLVER_FLAGS_ALLOWS_MMDATA |
76  SOLVER_FLAGS_ALLOWS_MDAPLUSAT |
77  SOLVER_FLAGS_ALLOWS_REVERSE_CUTHILL_MC_KEE |
78  SOLVER_FLAGS_ALLOWS_KING |
79  SOLVER_FLAGS_ALLOWS_SLOAN |
80  SOLVER_FLAGS_ALLOWS_NESTED_DISSECTION
81  };
Enumerator
EMPTY_SOLVER 
HARWELL_SOLVER 
LAPACK_SOLVER 
MESCHACH_SOLVER 
NAIVE_SOLVER 
SUPERLU_SOLVER 
TAUCS_SOLVER 
UMFPACK_SOLVER 
KLU_SOLVER 
Y12_SOLVER 
LAST_SOLVER 

Definition at line 41 of file linsol.h.

Constructor & Destructor Documentation

LinSol::LinSol ( void  )

Definition at line 156 of file linsol.cc.

References NO_OP.

158 solverFlags(0),
159 nThreads(1),
160 iWorkSpaceSize(0),
161 blockSize(0),
162 dPivotFactor(-1.),
163 dDropTolerance(0.),
164 iMaxIter(0) // Restore the original behavior by default
165 {
166  NO_OP;
167 }
static SolverType defaultSolver
Definition: linsol.h:151
doublereal dPivotFactor
Definition: linsol.h:127
#define NO_OP
Definition: myassert.h:74
unsigned blockSize
Definition: linsol.h:117
unsigned nThreads
Definition: linsol.h:103
integer iMaxIter
Definition: linsol.h:148
integer iWorkSpaceSize
Definition: linsol.h:110
unsigned solverFlags
Definition: linsol.h:96
SolverType currSolver
Definition: linsol.h:95
doublereal dDropTolerance
Definition: linsol.h:134
LinSol::~LinSol ( void  )
virtual

Definition at line 169 of file linsol.cc.

References NO_OP.

170 {
171  NO_OP;
172 }
#define NO_OP
Definition: myassert.h:74

Member Function Documentation

bool LinSol::AddSolverFlags ( unsigned  f)

Definition at line 292 of file linsol.cc.

References currSolver, and solverFlags.

Referenced by InverseSolver::Prepare(), Solver::Prepare(), and ReadLinSol().

293 {
294  if ((::solver[currSolver].s_flags & f) == f) {
295  solverFlags |= f;
296  return true;
297  }
298 
299  return false;
300 }
const LinSol::solver_t solver[]
Definition: linsol.cc:58
unsigned solverFlags
Definition: linsol.h:96
SolverType currSolver
Definition: linsol.h:95
const doublereal & LinSol::dGetDropTolerance ( void  ) const

Definition at line 351 of file linsol.cc.

References dDropTolerance.

352 {
353  return dDropTolerance;
354 }
doublereal dDropTolerance
Definition: linsol.h:134
const doublereal & LinSol::dGetPivotFactor ( void  ) const

Definition at line 345 of file linsol.cc.

References dPivotFactor.

Referenced by RestartLinSol().

346 {
347  return dPivotFactor;
348 }
doublereal dPivotFactor
Definition: linsol.h:127
unsigned LinSol::GetBlockSize ( void  ) const

Definition at line 396 of file linsol.cc.

References blockSize.

Referenced by RestartLinSol().

397 {
398  return blockSize;
399 }
unsigned blockSize
Definition: linsol.h:117
integer LinSol::GetMaxIterations ( void  ) const

Definition at line 434 of file linsol.cc.

References iMaxIter.

435 {
436  return iMaxIter;
437 }
integer iMaxIter
Definition: linsol.h:148
unsigned LinSol::GetNumThreads ( void  ) const

Definition at line 333 of file linsol.cc.

References nThreads.

Referenced by RestartLinSol().

334 {
335  return nThreads;
336 }
unsigned nThreads
Definition: linsol.h:103
const SolutionManager::ScaleOpt& LinSol::GetScale ( void  ) const
inline

Definition at line 164 of file linsol.h.

References scale.

164 { return scale; }
SolutionManager::ScaleOpt scale
Definition: linsol.h:141
SolutionManager *const LinSol::GetSolutionManager ( integer  iNLD,
integer  iLWS = 0 
) const

Definition at line 455 of file linsol.cc.

References ASSERT, blockSize, currSolver, dDropTolerance, dPivotFactor, EMPTY_SOLVER, HARWELL_SOLVER, iMaxIter, iWorkSpaceSize, KLU_SOLVER, LAPACK_SOLVER, MBDYN_EXCEPT_ARGS, MESCHACH_SOLVER, NAIVE_SOLVER, NO_OP, nThreads, SAFENEWWITHCONSTRUCTOR, scale, SOLVER_FLAGS_ALLOWS_CC, SOLVER_FLAGS_ALLOWS_COLAMD, SOLVER_FLAGS_ALLOWS_DIR, SOLVER_FLAGS_ALLOWS_KING, SOLVER_FLAGS_ALLOWS_MDAPLUSAT, SOLVER_FLAGS_ALLOWS_MMDATA, SOLVER_FLAGS_ALLOWS_MT_FCT, SOLVER_FLAGS_ALLOWS_NESTED_DISSECTION, SOLVER_FLAGS_ALLOWS_REVERSE_CUTHILL_MC_KEE, SOLVER_FLAGS_ALLOWS_SLOAN, SOLVER_FLAGS_PERM_MASK, SOLVER_FLAGS_TYPE_MASK, solverFlags, SUPERLU_SOLVER, TAUCS_SOLVER, UMFPACK_SOLVER, and Y12_SOLVER.

Referenced by Solver::AllocateSolman(), NRTrim::DoTrim(), and DataManager::InitialJointAssembly().

456 {
457  SolutionManager *pCurrSM = NULL;
458  const unsigned type = (solverFlags & LinSol::SOLVER_FLAGS_TYPE_MASK);
459  const unsigned perm = (solverFlags & LinSol::SOLVER_FLAGS_PERM_MASK);
461 
462  /* silence warning */
463  if (mt) {
464  NO_OP;
465  }
466 
467  ASSERT((::solver[currSolver].s_flags & solverFlags) == solverFlags);
468 
469  if (iLWS == 0) {
470  iLWS = iWorkSpaceSize;
471  }
472 
473  switch (currSolver) {
474  case LinSol::Y12_SOLVER:
475 #ifdef USE_Y12
476  switch (type) {
478  typedef Y12SparseCCSolutionManager<DirCColMatrixHandler<1> > CCSM;
479  SAFENEWWITHCONSTRUCTOR(pCurrSM, CCSM,
480  CCSM(iNLD, iLWS, dPivotFactor));
481  break;
482  }
483 
485  typedef Y12SparseCCSolutionManager<CColMatrixHandler<1> > CCSM;
486  SAFENEWWITHCONSTRUCTOR(pCurrSM, CCSM,
487  CCSM(iNLD, iLWS, dPivotFactor));
488  break;
489  }
490 
491  default:
492  SAFENEWWITHCONSTRUCTOR(pCurrSM,
493  Y12SparseSolutionManager,
494  Y12SparseSolutionManager(iNLD, iLWS,
495  dPivotFactor));
496  break;
497  }
498  break;
499 #else /* !USE_Y12 */
500  silent_cerr("Configure with --with-y12 "
501  "to enable Y12 solver" << std::endl);
503 #endif /* !USE_Y12 */
504 
506 #ifdef USE_SUPERLU
507 #ifdef USE_SUPERLU_MT
508  switch (type) {
510  typedef ParSuperLUSparseCCSolutionManager<DirCColMatrixHandler<0> > CCSM;
511  SAFENEWWITHCONSTRUCTOR(pCurrSM, CCSM,
512  CCSM(nThreads, iNLD, dPivotFactor));
513  break;
514  }
515 
517  typedef ParSuperLUSparseCCSolutionManager<CColMatrixHandler<0> > CCSM;
518  SAFENEWWITHCONSTRUCTOR(pCurrSM, CCSM,
519  CCSM(nThreads, iNLD, dPivotFactor));
520  break;
521  }
522 
523  default:
524  SAFENEWWITHCONSTRUCTOR(pCurrSM,
525  ParSuperLUSparseSolutionManager,
526  ParSuperLUSparseSolutionManager(nThreads, iNLD,
527  dPivotFactor));
528  break;
529  }
530  break;
531 
532  if (nThreads == 1) {
533  silent_cerr("warning, using multithread SuperLU with only one thread; "
534  << std::endl);
535  }
536 #else /* !USE_SUPERLU_MT */
537  if (nThreads > 1 && mt) {
538  silent_cerr("multithread SuperLU solver support not compiled; "
539  << std::endl);
541  }
542  switch (type) {
544  typedef SuperLUSparseCCSolutionManager<DirCColMatrixHandler<0> > CCSM;
545  SAFENEWWITHCONSTRUCTOR(pCurrSM, CCSM,
546  CCSM(iNLD, dPivotFactor));
547  break;
548  }
549 
551  typedef SuperLUSparseCCSolutionManager<CColMatrixHandler<0> > CCSM;
552  SAFENEWWITHCONSTRUCTOR(pCurrSM, CCSM,
553  CCSM(iNLD, dPivotFactor));
554  break;
555  }
556 
557  default:
558  SAFENEWWITHCONSTRUCTOR(pCurrSM,
559  SuperLUSparseSolutionManager,
560  SuperLUSparseSolutionManager(iNLD,
561  dPivotFactor));
562  break;
563  }
564 #endif /* !USE_SUPERLU_MT */
565  break;
566 #else /* !USE_SUPERLU */
567  silent_cerr("Configure with --with-superlu "
568  "to enable superlu solver" << std::endl);
570 #endif /* !USE_SUPERLU */
571 
573 #ifdef USE_MESCHACH
574  SAFENEWWITHCONSTRUCTOR(pCurrSM,
575  MeschachSparseSolutionManager,
576  MeschachSparseSolutionManager(iNLD, iLWS,
577  dPivotFactor));
578  break;
579 #else /* !USE_MESCHACH */
580  silent_cerr("Configure with --with-meschach "
581  "to enable Meschach solver" << std::endl);
583 #endif /* !USE_MESCHACH */
584 
586 #ifdef USE_LAPACK
587  SAFENEWWITHCONSTRUCTOR(pCurrSM,
588  LapackSolutionManager,
589  LapackSolutionManager(iNLD, dPivotFactor));
590  break;
591 #else /* !USE_LAPACK */
592  silent_cerr("Configure with --with-lapack "
593  "to enable Lapack dense solver" << std::endl);
595 #endif /* !USE_LAPACK */
596 
597  case LinSol::TAUCS_SOLVER:
598 #ifdef USE_TAUCS
599  switch (type) {
601  typedef TaucsSparseCCSolutionManager<DirCColMatrixHandler<0> > CCSM;
602  SAFENEWWITHCONSTRUCTOR(pCurrSM, CCSM, CCSM(iNLD));
603  break;
604  }
605 
607  typedef TaucsSparseCCSolutionManager<CColMatrixHandler<0> > CCSM;
608  SAFENEWWITHCONSTRUCTOR(pCurrSM, CCSM, CCSM(iNLD));
609  break;
610  }
611 
612  default:
613  SAFENEWWITHCONSTRUCTOR(pCurrSM,
614  TaucsSparseSolutionManager,
615  TaucsSparseSolutionManager(iNLD));
616  break;
617  }
618  break;
619 #else /* !USE_TAUCS */
620  silent_cerr("Configure with --with-taucs "
621  "to enable Taucs sparse solver" << std::endl);
623 #endif /* !USE_TAUCS */
624 
626 #ifdef USE_HARWELL
627  SAFENEWWITHCONSTRUCTOR(pCurrSM,
628  HarwellSparseSolutionManager,
629  HarwellSparseSolutionManager(iNLD, iLWS,
630  dPivotFactor));
631  break;
632 #else /* !USE_HARWELL */
633  silent_cerr("Configure with --with-harwell "
634  "to enable Harwell solver" << std::endl);
636 #endif /* !USE_HARWELL */
637 
639 #ifdef USE_UMFPACK
640  {
641  switch (type) {
643  typedef UmfpackSparseCCSolutionManager<DirCColMatrixHandler<0> > CCSM;
644  SAFENEWWITHCONSTRUCTOR(pCurrSM, CCSM,
646  break;
647  }
648 
650  typedef UmfpackSparseCCSolutionManager<CColMatrixHandler<0> > CCSM;
651  SAFENEWWITHCONSTRUCTOR(pCurrSM, CCSM,
653  break;
654  }
655 
656  default:
657  SAFENEWWITHCONSTRUCTOR(pCurrSM,
658  UmfpackSparseSolutionManager,
659  UmfpackSparseSolutionManager(iNLD,
661  break;
662  }
663  } break;
664 #else /* !USE_UMFPACK */
665  silent_cerr("Configure with --with-umfpack "
666  "to enable Umfpack solver" << std::endl);
668 #endif /* !USE_UMFPACK */
669 
670  case LinSol::KLU_SOLVER:
671 #ifdef USE_KLU
672  {
673  switch (type) {
675  typedef KLUSparseCCSolutionManager<DirCColMatrixHandler<0> > CCSM;
676  SAFENEWWITHCONSTRUCTOR(pCurrSM, CCSM,
677  CCSM(iNLD, dPivotFactor, scale));
678  break;
679  }
680 
682  typedef KLUSparseCCSolutionManager<CColMatrixHandler<0> > CCSM;
683  SAFENEWWITHCONSTRUCTOR(pCurrSM, CCSM,
684  CCSM(iNLD, dPivotFactor, scale));
685  break;
686  }
687 
688  default:
689  SAFENEWWITHCONSTRUCTOR(pCurrSM,
690  KLUSparseSolutionManager,
691  KLUSparseSolutionManager(iNLD,
692  dPivotFactor,
693  scale));
694  break;
695  }
696 
697  } break;
698 #else /* !USE_KLU */
699  silent_cerr("Configure with --with-klu "
700  "to enable KLU solver" << std::endl);
702 #endif /* !USE_KLU */
703 
706  if (nThreads == 1) {
707  SAFENEWWITHCONSTRUCTOR(pCurrSM,
710  } else {
711 #ifdef USE_NAIVE_MULTITHREAD
712  SAFENEWWITHCONSTRUCTOR(pCurrSM,
713  ParNaiveSparsePermSolutionManager,
714  ParNaiveSparsePermSolutionManager(nThreads, iNLD, dPivotFactor));
715 #else /* ! USE_NAIVE_MULTITHREAD */
716  silent_cerr("multithread naive solver support not compiled; "
717  "you can configure --enable-multithread-naive "
718  "on a linux ix86 to get it"
719  << std::endl);
721 #endif /* ! USE_NAIVE_MULTITHREAD */
722  }
723 #ifdef USE_BOOST
724 #ifdef HAVE_BOOST_GRAPH_CUTHILL_MCKEE_ORDERING_HPP
726  if (nThreads == 1) {
727  SAFENEWWITHCONSTRUCTOR(pCurrSM,
730  } else {
731 #ifdef USE_NAIVE_MULTITHREAD
732 #if 0
733  SAFENEWWITHCONSTRUCTOR(pCurrSM,
734  ParNaiveSparsePermSolutionManager,
735  ParNaiveSparsePermSolutionManager(nThreads, iNLD, dPivotFactor));
736 #endif
737  silent_cerr("multithread naive solver with"
738  "reverse Cuthill-McKee permutation not"
739  "available yet. Patches welcome"
740  << std::endl);
742 #else /* ! USE_NAIVE_MULTITHREAD */
743  silent_cerr("multithread naive solver support not compiled; "
744  "you can configure --enable-multithread-naive "
745  "on a linux ix86 to get it"
746  << std::endl);
748 #endif /* ! USE_NAIVE_MULTITHREAD */
749  }
750 #endif /* HAVE_BOOST_GRAPH_CUTHILL_MCKEE_ORDERING_HPP */
751 
752 #if 0 /* ?!? */
753  } else if (perm == LinSol::SOLVER_FLAGS_ALLOWS_MMDATA) {
754  if (nThreads == 1) {
755  SAFENEWWITHCONSTRUCTOR(pCurrSM,
758  } else {
759 #ifdef USE_NAIVE_MULTITHREAD
760 #if 0
761  SAFENEWWITHCONSTRUCTOR(pCurrSM,
762  ParNaiveSparsePermSolutionManager,
763  ParNaiveSparsePermSolutionManager(nThreads, iNLD, dPivotFactor));
764 #endif
765  silent_cerr("multithread naive solver with"
766  "approximate minimum degree permutation not"
767  "available yet. Patches welcome"
768  << std::endl);
770 #else
771  silent_cerr("multithread naive solver support not compiled; "
772  "you can configure --enable-multithread-naive "
773  "on a linux ix86 to get it"
774  << std::endl);
776 #endif /* USE_NAIVE_MULTITHREAD */
777  }
778 #endif
779 
780 #ifdef HAVE_BOOST_GRAPH_KING_ORDERING_HPP
781  } else if (perm == LinSol::SOLVER_FLAGS_ALLOWS_KING) {
782  if (nThreads == 1) {
783  SAFENEWWITHCONSTRUCTOR(pCurrSM,
786  } else {
787 #ifdef USE_NAIVE_MULTITHREAD
788 #if 0
789  SAFENEWWITHCONSTRUCTOR(pCurrSM,
790  ParNaiveSparsePermSolutionManager,
791  ParNaiveSparsePermSolutionManager(nThreads, iNLD, dPivotFactor));
792 #endif
793  silent_cerr("multithread naive solver with"
794  "king permutation not"
795  "available yet. Patches welcome"
796  << std::endl);
798 #else /* ! USE_NAIVE_MULTITHREAD */
799  silent_cerr("multithread naive solver support not compiled; "
800  "you can configure --enable-multithread-naive "
801  "on a linux ix86 to get it"
802  << std::endl);
804 #endif /* ! USE_NAIVE_MULTITHREAD */
805  }
806 #endif /* HAVE_BOOST_GRAPH_KING_ORDERING_HPP */
807 
808 #ifdef HAVE_BOOST_GRAPH_SLOAN_ORDERING_HPP
809  } else if (perm == LinSol::SOLVER_FLAGS_ALLOWS_SLOAN) {
810  if (nThreads == 1) {
811  SAFENEWWITHCONSTRUCTOR(pCurrSM,
814  } else {
815 #ifdef USE_NAIVE_MULTITHREAD
816 #if 0
817  SAFENEWWITHCONSTRUCTOR(pCurrSM,
818  ParNaiveSparsePermSolutionManager,
819  ParNaiveSparsePermSolutionManager(nThreads, iNLD, dPivotFactor));
820 #endif
821  silent_cerr("multithread naive solver with"
822  "sloan permutation not"
823  "available yet. Patches welcome"
824  << std::endl);
826 #else /* ! USE_NAIVE_MULTITHREAD */
827  silent_cerr("multithread naive solver support not compiled; "
828  "you can configure --enable-multithread-naive "
829  "on a linux ix86 to get it"
830  << std::endl);
832 #endif /* ! USE_NAIVE_MULTITHREAD */
833  }
834 #endif /* HAVE_BOOST_GRAPH_SLOAN_ORDERING_HPP */
835 
836 #ifdef HAVE_BOOST_GRAPH_MINIMUM_DEGREE_ORDERING_HPP
837  } else if (perm == LinSol::SOLVER_FLAGS_ALLOWS_MDAPLUSAT) {
838  if (nThreads == 1) {
839  SAFENEWWITHCONSTRUCTOR(pCurrSM,
842  } else {
843 #ifdef USE_NAIVE_MULTITHREAD
844 #if 0
845  SAFENEWWITHCONSTRUCTOR(pCurrSM,
846  ParNaiveSparsePermSolutionManager,
847  ParNaiveSparsePermSolutionManager(nThreads, iNLD, dPivotFactor));
848 #endif
849  silent_cerr("multithread naive solver with"
850  "md permutation not"
851  "available yet. Patches welcome"
852  << std::endl);
854 #else /* ! USE_NAIVE_MULTITHREAD */
855  silent_cerr("multithread naive solver support not compiled; "
856  "you can configure --enable-multithread-naive "
857  "on a linux ix86 to get it"
858  << std::endl);
860 #endif /* ! USE_NAIVE_MULTITHREAD */
861  }
862 #endif /* HAVE_BOOST_GRAPH_MINIMUM_DEGREE_ORDERING_HPP */
863 #endif /* USE_BOOST */
864 
866 #ifdef USE_METIS
867  if (nThreads == 1) {
868  SAFENEWWITHCONSTRUCTOR(pCurrSM,
871  } else {
872 #ifdef USE_NAIVE_MULTITHREAD
873 #if 0
874  SAFENEWWITHCONSTRUCTOR(pCurrSM,
875  ParNaiveSparsePermSolutionManager,
876  ParNaiveSparsePermSolutionManager(nThreads, iNLD, dPivotFactor));
877 #endif
878  silent_cerr("multithread naive solver with"
879  "nested dissection permutation not"
880  "available yet. Patches welcome"
881  << std::endl);
883 #else /* ! USE_NAIVE_MULTITHREAD */
884  silent_cerr("multithread naive solver support not compiled; "
885  "you can configure --enable-multithread-naive "
886  "on a linux ix86 to get it"
887  << std::endl);
889 #endif /* ! USE_NAIVE_MULTITHREAD */
890  }
891 #else /* ! USE_METIS */
892  silent_cerr("you should not get here("<< __FILE__ << ":" <<
893  __LINE__ << ")" << std::endl);
895 #endif /* ! USE_METIS */
896 
897  } else {
898  if (nThreads == 1) {
899  SAFENEWWITHCONSTRUCTOR(pCurrSM,
902  } else {
903 #ifdef USE_NAIVE_MULTITHREAD
904  SAFENEWWITHCONSTRUCTOR(pCurrSM,
905  ParNaiveSparseSolutionManager,
906  ParNaiveSparseSolutionManager(nThreads, iNLD, dPivotFactor));
907 #else
908  silent_cerr("multithread naive solver support not compiled; "
909  "you can configure --enable-multithread-naive "
910  "on a linux ix86 to get it"
911  << std::endl);
913 #endif /* USE_NAIVE_MULTITHREAD */
914  }
915  }
916  break;
917 
919  break;
920 
921  default:
922  ASSERT(0);
924 
925  }
926 
927  return pCurrSM;
928 }
#define MBDYN_EXCEPT_ARGS
Definition: except.h:63
doublereal dPivotFactor
Definition: linsol.h:127
#define NO_OP
Definition: myassert.h:74
unsigned blockSize
Definition: linsol.h:117
unsigned nThreads
Definition: linsol.h:103
const LinSol::solver_t solver[]
Definition: linsol.cc:58
integer iMaxIter
Definition: linsol.h:148
SolutionManager::ScaleOpt scale
Definition: linsol.h:141
integer iWorkSpaceSize
Definition: linsol.h:110
unsigned solverFlags
Definition: linsol.h:96
#define ASSERT(expression)
Definition: colamd.c:977
#define SAFENEWWITHCONSTRUCTOR(pnt, item, constructor)
Definition: mynewmem.h:698
SolverType currSolver
Definition: linsol.h:95
doublereal dDropTolerance
Definition: linsol.h:134
LinSol::SolverType LinSol::GetSolver ( void  ) const

Definition at line 175 of file linsol.cc.

References currSolver.

Referenced by Solver::AllocateSchurSolman(), Solver::AllocateSolman(), ReadLinSol(), and RestartLinSol().

176 {
177  return currSolver;
178 }
SolverType currSolver
Definition: linsol.h:95
unsigned LinSol::GetSolverFlags ( void  ) const

Definition at line 257 of file linsol.cc.

References solverFlags.

Referenced by InverseSolver::Prepare(), Solver::Prepare(), ReadLinSol(), RestartLinSol(), and SetNumThreads().

258 {
259  return solverFlags;
260 }
unsigned solverFlags
Definition: linsol.h:96
unsigned LinSol::GetSolverFlags ( SolverType  t) const

Definition at line 263 of file linsol.cc.

References LinSol::solver_t::s_flags, and solver.

264 {
266 }
const LinSol::solver_t solver[]
Definition: linsol.cc:58
unsigned s_flags
Definition: linsol.h:88
const char *const LinSol::GetSolverName ( void  ) const

Definition at line 269 of file linsol.cc.

References currSolver, LinSol::solver_t::s_name, and solver.

Referenced by Solver::AllocateSchurSolman(), InverseSolver::Prepare(), Solver::Prepare(), InverseSolver::ReadData(), Solver::ReadData(), ReadLinSol(), and RestartLinSol().

270 {
272 }
const char *const s_name
Definition: linsol.h:85
const LinSol::solver_t solver[]
Definition: linsol.cc:58
SolverType currSolver
Definition: linsol.h:95
const char *const LinSol::GetSolverName ( SolverType  t) const

Definition at line 275 of file linsol.cc.

References LinSol::solver_t::s_name, and solver.

276 {
278 }
const char *const s_name
Definition: linsol.h:85
const LinSol::solver_t solver[]
Definition: linsol.cc:58
integer LinSol::iGetWorkSpaceSize ( void  ) const

Definition at line 339 of file linsol.cc.

References iWorkSpaceSize.

Referenced by RestartLinSol(), and Solver::SetupSolmans().

340 {
341  return iWorkSpaceSize;
342 }
integer iWorkSpaceSize
Definition: linsol.h:110
bool LinSol::MaskSolverFlags ( unsigned  f)

Definition at line 303 of file linsol.cc.

References currSolver, and solverFlags.

Referenced by ReadLinSol().

304 {
305  if ((::solver[currSolver].s_flags & f) == f) {
306  solverFlags &= ~f;
307  return true;
308  }
309 
310  return false;
311 }
const LinSol::solver_t solver[]
Definition: linsol.cc:58
unsigned solverFlags
Definition: linsol.h:96
SolverType currSolver
Definition: linsol.h:95
bool LinSol::SetBlockSize ( unsigned  bs)

Definition at line 402 of file linsol.cc.

References blockSize, currSolver, and UMFPACK_SOLVER.

Referenced by ReadLinSol().

403 {
404  switch (currSolver) {
406  blockSize = bs;
407  break;
408 
409  default:
410  return false;
411  }
412 
413  return true;
414 }
unsigned blockSize
Definition: linsol.h:117
SolverType currSolver
Definition: linsol.h:95
bool LinSol::SetDropTolerance ( const doublereal d)

Definition at line 384 of file linsol.cc.

References currSolver, and dDropTolerance.

Referenced by ReadLinSol().

385 {
386  if (::solver[currSolver].s_drop_tolerance == -1.) {
387  return false;
388  }
389 
390  dDropTolerance = d;
391 
392  return true;
393 }
const LinSol::solver_t solver[]
Definition: linsol.cc:58
SolverType currSolver
Definition: linsol.h:95
doublereal dDropTolerance
Definition: linsol.h:134
bool LinSol::SetMaxIterations ( integer  iMaxIter)

Definition at line 440 of file linsol.cc.

References currSolver, iMaxIter, and UMFPACK_SOLVER.

Referenced by ReadLinSol().

441 {
442  switch (currSolver) {
444  iMaxIter = iMaxIterations;
445  break;
446 
447  default:
448  return false;
449  }
450 
451  return true;
452 }
integer iMaxIter
Definition: linsol.h:148
SolverType currSolver
Definition: linsol.h:95
bool LinSol::SetNumThreads ( unsigned  nt)

Definition at line 314 of file linsol.cc.

References currSolver, GetSolverFlags(), nThreads, SOLVER_FLAGS_ALLOWS_MT_FCT, and solverFlags.

Referenced by InverseSolver::ReadData(), Solver::ReadData(), and ReadLinSol().

315 {
317  if (nt == 0) {
319 
320  } else {
322  }
323 
324  nThreads = nt;
325 
326  return true;
327  }
328 
329  return false;
330 }
unsigned GetSolverFlags(void) const
Definition: linsol.cc:257
Definition: linsol.h:39
unsigned nThreads
Definition: linsol.h:103
unsigned solverFlags
Definition: linsol.h:96
SolverType currSolver
Definition: linsol.h:95

Here is the call graph for this function:

bool LinSol::SetPivotFactor ( const doublereal d)

Definition at line 372 of file linsol.cc.

References currSolver, and dPivotFactor.

Referenced by ReadLinSol().

373 {
374  if (::solver[currSolver].s_pivot_factor == -1.) {
375  return false;
376  }
377 
378  dPivotFactor = d;
379 
380  return true;
381 }
doublereal dPivotFactor
Definition: linsol.h:127
const LinSol::solver_t solver[]
Definition: linsol.cc:58
SolverType currSolver
Definition: linsol.h:95
bool LinSol::SetScale ( const SolutionManager::ScaleOpt scale)

Definition at line 417 of file linsol.cc.

References currSolver, KLU_SOLVER, NAIVE_SOLVER, scale, and UMFPACK_SOLVER.

Referenced by ReadLinSol().

418 {
419  switch (currSolver) {
422  case LinSol::KLU_SOLVER:
423  scale = s;
424  break;
425 
426  default:
427  return false;
428  }
429 
430  return true;
431 }
SolverType currSolver
Definition: linsol.h:95
bool LinSol::SetSolver ( LinSol::SolverType  t,
unsigned  f = SOLVER_FLAGS_NONE 
)

Definition at line 181 of file linsol.cc.

References currSolver, EMPTY_SOLVER, HARWELL_SOLVER, KLU_SOLVER, LAPACK_SOLVER, MESCHACH_SOLVER, NAIVE_SOLVER, solverFlags, SUPERLU_SOLVER, TAUCS_SOLVER, UMFPACK_SOLVER, and Y12_SOLVER.

Referenced by InverseSolver::Prepare(), Solver::Prepare(), and ReadLinSol().

182 {
183  /* check flags */
184  if ((::solver[t].s_flags & f) != f) {
185  return false;
186  }
187 
188  solverFlags = f;
189 
190  switch (t) {
192 #ifdef USE_UMFPACK
193  currSolver = t;
194  return true;
195 #endif /* USE_UMFPACK */
196 
197  case LinSol::KLU_SOLVER:
198 #ifdef USE_KLU
199  currSolver = t;
200  return true;
201 #endif /* USE_UMFPACK */
202 
204 #ifdef USE_SUPERLU
205  currSolver = t;
206  return true;
207 #endif /* USE_SUPERLU */
208 
210 #ifdef USE_LAPACK
211  currSolver = t;
212  return true;
213 #endif /* USE_LAPACK */
214 
216 #ifdef USE_TAUCS
217  currSolver = t;
218  return true;
219 #endif /* USE_TAUCS */
220 
221  case LinSol::Y12_SOLVER:
222 #ifdef USE_Y12
223  currSolver = t;
224  return true;
225 #endif /* USE_Y12 */
226 
228 #ifdef USE_HARWELL
229  currSolver = t;
230  return true;
231 #endif /* USE_HARWELL */
232 
234 #ifdef USE_MESCHACH
235  currSolver = t;
236  return true;
237 #endif /* USE_MESCHACH */
238 
239  /* else */
240  silent_cerr(::solver[t].s_name << " unavailable" << std::endl);
241  return false;
242 
244  currSolver = t;
245  return true;
246 
248  currSolver = t;
249  return true;
250 
251  default:
252  return false;
253  }
254 }
const LinSol::solver_t solver[]
Definition: linsol.cc:58
unsigned solverFlags
Definition: linsol.h:96
SolverType currSolver
Definition: linsol.h:95
bool LinSol::SetSolverFlags ( unsigned  f)

Definition at line 281 of file linsol.cc.

References currSolver, and solverFlags.

Referenced by ReadLinSol().

282 {
283  if ((::solver[currSolver].s_flags & f) == f) {
284  solverFlags = f;
285  return true;
286  }
287 
288  return false;
289 }
const LinSol::solver_t solver[]
Definition: linsol.cc:58
unsigned solverFlags
Definition: linsol.h:96
SolverType currSolver
Definition: linsol.h:95
bool LinSol::SetWorkSpaceSize ( integer  i)

Definition at line 357 of file linsol.cc.

References currSolver, iWorkSpaceSize, and Y12_SOLVER.

Referenced by ReadLinSol().

358 {
359  switch (currSolver) {
360  case LinSol::Y12_SOLVER:
361  iWorkSpaceSize = i;
362  break;
363 
364  default:
365  return false;
366  }
367 
368  return true;
369 }
integer iWorkSpaceSize
Definition: linsol.h:110
SolverType currSolver
Definition: linsol.h:95

Member Data Documentation

unsigned LinSol::blockSize
protected

Definition at line 117 of file linsol.h.

Referenced by GetBlockSize(), GetSolutionManager(), and SetBlockSize().

doublereal LinSol::dDropTolerance
protected

Definition at line 134 of file linsol.h.

Referenced by dGetDropTolerance(), GetSolutionManager(), and SetDropTolerance().

LinSol::SolverType LinSol::defaultSolver
static
Initial value:

Definition at line 151 of file linsol.h.

doublereal LinSol::dPivotFactor
protected

Definition at line 127 of file linsol.h.

Referenced by dGetPivotFactor(), GetSolutionManager(), and SetPivotFactor().

integer LinSol::iMaxIter
protected

Definition at line 148 of file linsol.h.

Referenced by GetMaxIterations(), GetSolutionManager(), and SetMaxIterations().

integer LinSol::iWorkSpaceSize
protected

Definition at line 110 of file linsol.h.

Referenced by GetSolutionManager(), iGetWorkSpaceSize(), and SetWorkSpaceSize().

unsigned LinSol::nThreads
protected

Definition at line 103 of file linsol.h.

Referenced by GetNumThreads(), GetSolutionManager(), and SetNumThreads().

SolutionManager::ScaleOpt LinSol::scale
protected

Definition at line 141 of file linsol.h.

Referenced by GetScale(), GetSolutionManager(), and SetScale().

unsigned LinSol::solverFlags
protected

The documentation for this class was generated from the following files: