MBDyn-1.7.3
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups
readlinsol.cc File Reference
#include "mbconfig.h"
#include "ac/sys_sysinfo.h"
#include "dataman.h"
#include "readlinsol.h"
Include dependency graph for readlinsol.cc:

Go to the source code of this file.

Functions

void ReadLinSol (LinSol &cs, HighParser &HP, bool bAllowEmpty)
 
std::ostream & RestartLinSol (std::ostream &out, const LinSol &cs)
 

Function Documentation

void ReadLinSol ( LinSol cs,
HighParser HP,
bool  bAllowEmpty 
)

USE_METIS

Definition at line 39 of file readlinsol.cc.

References LinSol::AddSolverFlags(), SolutionManager::ScaleOpt::algorithm, DEBUGLCOUT, SolutionManager::ScaleOpt::dTol, EMPTY, LinSol::EMPTY_SOLVER, get_nprocs(), HighParser::GetInt(), HighParser::GetLineData(), HighParser::GetReal(), LinSol::GetSolver(), LinSol::GetSolverFlags(), LinSol::GetSolverName(), HighParser::GetWord(), HighParser::GetYesNoOrBool(), LinSol::HARWELL_SOLVER, SolutionManager::ScaleOpt::iMaxIter, HighParser::IsKeyWord(), LinSol::KLU_SOLVER, LinSol::LAPACK_SOLVER, LASTKEYWORD, LinSol::MaskSolverFlags(), MBDYN_EXCEPT_ARGS, LinSol::MESCHACH_SOLVER, MYDEBUG_INPUT, LinSol::NAIVE_SOLVER, LinSol::solver_t::s_alias, LinSol::solver_t::s_default_flags, LinSol::solver_t::s_drop_tolerance, LinSol::solver_t::s_flags, LinSol::solver_t::s_name, LinSol::solver_t::s_pivot_factor, SolutionManager::SCALEA_COL_MAX, SolutionManager::SCALEA_COL_SUM, SolutionManager::SCALEA_ITERATIVE, SolutionManager::SCALEA_LAPACK, SolutionManager::SCALEA_NONE, SolutionManager::SCALEA_ROW_MAX, SolutionManager::SCALEA_ROW_MAX_COL_MAX, SolutionManager::SCALEA_ROW_SUM, SolutionManager::SCALEA_UNDEF, SolutionManager::SCALEF_COND_NUM_1, SolutionManager::SCALEF_COND_NUM_INF, SolutionManager::SCALEF_VERBOSE, SolutionManager::SCALEF_WARN, SolutionManager::SCALEW_ALWAYS, SolutionManager::SCALEW_NEVER, SolutionManager::SCALEW_ONCE, LinSol::SetBlockSize(), LinSol::SetDropTolerance(), LinSol::SetMaxIterations(), LinSol::SetNumThreads(), LinSol::SetPivotFactor(), LinSol::SetScale(), LinSol::SetSolver(), LinSol::SetSolverFlags(), LinSol::SetWorkSpaceSize(), solver, LinSol::SOLVER_FLAGS_ALLOWS_CC, LinSol::SOLVER_FLAGS_ALLOWS_COLAMD, LinSol::SOLVER_FLAGS_ALLOWS_DIR, LinSol::SOLVER_FLAGS_ALLOWS_KING, LinSol::SOLVER_FLAGS_ALLOWS_MAP, LinSol::SOLVER_FLAGS_ALLOWS_MDAPLUSAT, LinSol::SOLVER_FLAGS_ALLOWS_MMDATA, LinSol::SOLVER_FLAGS_ALLOWS_MT_FCT, LinSol::SOLVER_FLAGS_ALLOWS_NESTED_DISSECTION, LinSol::SOLVER_FLAGS_ALLOWS_REVERSE_CUTHILL_MC_KEE, LinSol::SOLVER_FLAGS_TYPE_MASK, LinSol::SUPERLU_SOLVER, LinSol::TAUCS_SOLVER, SolutionManager::ScaleOpt::uFlags, LinSol::UMFPACK_SOLVER, SolutionManager::ScaleOpt::when, and LinSol::Y12_SOLVER.

Referenced by DataManager::ReadControl(), InverseSolver::ReadData(), and Solver::ReadData().

40 {
41  /* parole chiave */
42  const char* sKeyWords[] = {
54  NULL
55  };
56 
57  enum KeyWords {
58  EMPTY,
59  HARWELL,
60  LAPACK,
61  MESCHACH,
62  NAIVE,
63  SUPERLU,
64  TAUCS,
65  UMFPACK,
66  UMFPACK3,
67  KLU,
68  Y12,
69 
71  };
72 
73  /* tabella delle parole chiave */
74  KeyTable K(HP, sKeyWords);
75 
76  bool bGotIt = false;
77  switch (KeyWords(HP.GetWord())) {
78  case EMPTY:
79  if (!bAllowEmpty) {
80  silent_cerr("empty solver not allowed at line "
81  << HP.GetLineData() << std::endl);
83  }
84 
87  "No LU solver" << std::endl);
88  bGotIt = true;
89  break;
90 
91  case HARWELL:
94  "Using harwell sparse LU solver" << std::endl);
95 #ifdef USE_HARWELL
96  bGotIt = true;
97 #endif /* USE_HARWELL */
98  break;
99 
100  case LAPACK:
103  "Using Lapack dense LU solver" << std::endl);
104 #ifdef USE_LAPACK
105  bGotIt = true;
106 #endif /* USE_LAPACK */
107  break;
108 
109  case MESCHACH:
112  "Using meschach sparse LU solver"
113  << std::endl);
114 #ifdef USE_MESCHACH
115  bGotIt = true;
116 #endif /* USE_MESCHACH */
117  break;
118 
119  case NAIVE:
121  bGotIt = true;
122  break;
123 
124  case SUPERLU:
125  /*
126  * FIXME: use CC as default???
127  */
130  "Using SuperLU sparse LU solver" << std::endl);
131 #ifdef USE_SUPERLU
132  bGotIt = true;
133 #endif /* USE_SUPERLU */
134  break;
135 
136  case TAUCS:
137  /*
138  * FIXME: use CC as default???
139  */
142  "Using Taucs sparse solver" << std::endl);
143 #ifdef USE_TAUCS
144  bGotIt = true;
145 #endif /* USE_TAUCS */
146  break;
147 
148  case UMFPACK3:
149  pedantic_cerr("\"umfpack3\" is deprecated; "
150  "use \"umfpack\" instead" << std::endl);
151  case UMFPACK:
152  /*
153  * FIXME: use CC as default???
154  */
157  "Using umfpack sparse LU solver" << std::endl);
158 #ifdef USE_UMFPACK
159  bGotIt = true;
160 #endif /* USE_UMFPACK */
161  break;
162 
163  case KLU:
164  /*
165  * FIXME: use CC as default???
166  */
169  "Using KLU sparse LU solver" << std::endl);
170 #ifdef USE_KLU
171  bGotIt = true;
172 #endif /* USE_KLU */
173  break;
174 
175  case Y12:
176  /*
177  * FIXME: use CC as default???
178  */
181  "Using y12 sparse LU solver" << std::endl);
182 #ifdef USE_Y12
183  bGotIt = true;
184 #endif /* USE_Y12 */
185  break;
186 
187  default:
188  silent_cerr("unknown solver" << std::endl);
190  }
191 
192  const LinSol::solver_t currSolver = ::solver[cs.GetSolver()];
193 
194  if (!bGotIt) {
195  silent_cerr(currSolver.s_name << " solver "
196  "not available; requested at line " << HP.GetLineData()
197  << std::endl);
199  }
200 
201  cs.SetSolverFlags(currSolver.s_default_flags);
202 
203  /* map? */
204  if (HP.IsKeyWord("map")) {
205  if (currSolver.s_flags & LinSol::SOLVER_FLAGS_ALLOWS_MAP) {
208  pedantic_cout("using map matrix handling for "
209  << currSolver.s_name
210  << " solver" << std::endl);
211  } else {
212  pedantic_cerr("map is meaningless for "
213  << currSolver.s_name
214  << " solver" << std::endl);
215  }
216 
217  /* CC? */
218  } else if (HP.IsKeyWord("column" "compressed") || HP.IsKeyWord("cc")) {
219  if (currSolver.s_flags & LinSol::SOLVER_FLAGS_ALLOWS_CC) {
222  pedantic_cout("using column compressed matrix handling for "
223  << currSolver.s_name
224  << " solver" << std::endl);
225 
226  } else {
227  pedantic_cerr("column compressed is meaningless for "
228  << currSolver.s_name
229  << " solver" << std::endl);
230  }
231 
232  /* direct? */
233  } else if (HP.IsKeyWord("direct" "access") || HP.IsKeyWord("dir")) {
234  if (currSolver.s_flags & LinSol::SOLVER_FLAGS_ALLOWS_DIR) {
237  pedantic_cout("using direct access matrix handling for "
238  << currSolver.s_name
239  << " solver" << std::endl);
240  } else {
241  pedantic_cerr("direct is meaningless for "
242  << currSolver.s_name
243  << " solver" << std::endl);
244  }
245  }
246 
247  /* colamd? */
248  if (HP.IsKeyWord("colamd")) {
249  if (currSolver.s_flags & LinSol::SOLVER_FLAGS_ALLOWS_COLAMD) {
251  pedantic_cout("using colamd symmetric preordering for "
252  << currSolver.s_name
253  << " solver" << std::endl);
254  } else {
255  pedantic_cerr("colamd preordering is meaningless for "
256  << currSolver.s_name
257  << " solver" << std::endl);
258  }
259  /* amd ata? */
260  } else if (HP.IsKeyWord("mmdata")) {
261  silent_cerr("approximate minimum degree solver support is still TODO"
262  "task: detect (or import) the MD library;"
263  "uncomment the relevant bits in naivewrap;"
264  "remove this check (readlinsol.cc)."
265  "Patches welcome"
266  << std::endl);
268  if (currSolver.s_flags & LinSol::SOLVER_FLAGS_ALLOWS_MMDATA) {
270  pedantic_cout("using mmd symmetric preordering for "
271  << currSolver.s_name
272  << " solver" << std::endl);
273  } else {
274  pedantic_cerr("mmdata preordering is meaningless for "
275  << currSolver.s_name
276  << " solver" << std::endl);
277  }
278  /* minimum degree ?*/
279  } else if (HP.IsKeyWord("minimum" "degree")) {
282  pedantic_cout("using minimum degree symmetric preordering of A+A^T for "
283  << currSolver.s_name
284  << " solver" << std::endl);
285  } else {
286  pedantic_cerr("md preordering is meaningless for "
287  << currSolver.s_name
288  << " solver" << std::endl);
289  }
290  /* Reverse Kuthill McKee? */
291  } else if (HP.IsKeyWord("rcmk")) {
294  pedantic_cout("using rcmk symmetric preordering for "
295  << currSolver.s_name
296  << " solver" << std::endl);
297 
298  } else {
299  pedantic_cerr("rcmk preordering is meaningless for "
300  << currSolver.s_name
301  << " solver" << std::endl);
302  }
303  /* king ?*/
304  } else if (HP.IsKeyWord("king")) {
305  if (currSolver.s_flags & LinSol::SOLVER_FLAGS_ALLOWS_KING) {
307  pedantic_cout("using king symmetric preordering for "
308  << currSolver.s_name
309  << " solver" << std::endl);
310 
311  } else {
312  pedantic_cerr("king preordering is meaningless for "
313  << currSolver.s_name
314  << " solver" << std::endl);
315  }
316  /* sloan ? */
317  } else if (HP.IsKeyWord("sloan")) {
318  if (currSolver.s_flags & LinSol::SOLVER_FLAGS_ALLOWS_KING) {
320  pedantic_cout("using sloan symmetric preordering for "
321  << currSolver.s_name
322  << " solver" << std::endl);
323 
324  } else {
325  pedantic_cerr("sloan preordering is meaningless for "
326  << currSolver.s_name
327  << " solver" << std::endl);
328  }
329  /* nested dissection ? */
330  } else if (HP.IsKeyWord("nested" "dissection")) {
331 #ifdef USE_METIS
334  pedantic_cout("using nested dissection symmetric preordering for "
335  << currSolver.s_name
336  << " solver" << std::endl);
337 
338  } else {
339  pedantic_cerr("nested dissection preordering is meaningless for "
340  << currSolver.s_name
341  << " solver" << std::endl);
342  }
343 #else
344  silent_cerr("nested dissection permutation not built in;"
345  "please configure --with-metis to get it"
346  << std::endl);
348 #endif //USE_METIS
349  }
350 
351  /* multithread? */
352  if (HP.IsKeyWord("multi" "thread") || HP.IsKeyWord("mt")) {
353  int nThreads = HP.GetInt();
354 
355  if (currSolver.s_flags & LinSol::SOLVER_FLAGS_ALLOWS_MT_FCT) {
357  if (nThreads < 1) {
358  silent_cerr("illegal thread number, using 1" << std::endl);
359  nThreads = 1;
360  }
361  cs.SetNumThreads(nThreads);
362 
363  } else if (nThreads != 1) {
364  pedantic_cerr("multithread is meaningless for "
365  << currSolver.s_name
366  << " solver" << std::endl);
367  }
368 
369 
370  } else {
371  if (currSolver.s_flags & LinSol::SOLVER_FLAGS_ALLOWS_MT_FCT) {
372  int nThreads = get_nprocs();
373 
374  if (nThreads > 1) {
375  silent_cout("no multithread requested "
376  "with a potential "
377  "of " << nThreads << " CPUs"
378  << std::endl);
379  }
380 
381  cs.SetNumThreads(nThreads);
382  }
383  }
384 
385  if (HP.IsKeyWord("workspace" "size")) {
386  integer iWorkSpaceSize = HP.GetInt();
387  if (iWorkSpaceSize < 0) {
388  iWorkSpaceSize = 0;
389  }
390 
391  switch (cs.GetSolver()) {
392  case LinSol::Y12_SOLVER:
393  cs.SetWorkSpaceSize(iWorkSpaceSize);
394  break;
395 
396  default:
397  pedantic_cerr("workspace size is meaningless for "
398  << currSolver.s_name
399  << " solver" << std::endl);
400  break;
401  }
402  }
403 
404  if (HP.IsKeyWord("pivot" "factor")) {
405  doublereal dPivotFactor = HP.GetReal();
406 
407  if (currSolver.s_pivot_factor == -1.) {
408  pedantic_cerr("pivot factor is meaningless for "
409  << currSolver.s_name
410  << " solver" << std::endl);
411 
412  } else {
413  if (dPivotFactor <= 0. || dPivotFactor > 1.) {
414  silent_cerr("pivot factor " << dPivotFactor
415  << " is out of bounds; "
416  "using default "
417  "(" << currSolver.s_pivot_factor << ")"
418  << std::endl);
419  dPivotFactor = currSolver.s_pivot_factor;
420  }
421  cs.SetPivotFactor(dPivotFactor);
422  }
423 
424  } else {
425  if (currSolver.s_pivot_factor != -1.) {
426  cs.SetPivotFactor(currSolver.s_pivot_factor);
427  }
428  }
429 
430  if (HP.IsKeyWord("drop" "tolerance")) {
431  doublereal dDropTolerance = HP.GetReal();
432 
433  if (currSolver.s_drop_tolerance == -1.) {
434  pedantic_cerr("\"drop tolerance\" is meaningless for "
435  << currSolver.s_name
436  << " solver" << std::endl);
437 
438  } else {
439  if (dDropTolerance < 0.) {
440  silent_cerr("drop tolerance " << dDropTolerance
441  << " is out of bounds; "
442  "using default "
443  "(" << currSolver.s_drop_tolerance << ")"
444  << std::endl);
445  dDropTolerance = currSolver.s_drop_tolerance;
446  }
447  cs.SetDropTolerance(dDropTolerance);
448  }
449 
450  } else {
451  if (currSolver.s_drop_tolerance != -1.) {
452  cs.SetDropTolerance(currSolver.s_drop_tolerance);
453  }
454  }
455 
456  if (HP.IsKeyWord("block" "size")) {
457  integer blockSize = HP.GetInt();
458  if (blockSize < 1) {
459  silent_cerr("illegal negative block size; "
460  "using default" << std::endl);
461  blockSize = 0;
462  }
463 
464  switch (cs.GetSolver()) {
466  cs.SetBlockSize(blockSize);
467  break;
468 
469  default:
470  pedantic_cerr("block size is meaningless for "
471  << currSolver.s_name
472  << " solver" << std::endl);
473  break;
474  }
475  }
476 
477  if (HP.IsKeyWord("scale")) {
478  switch (cs.GetSolver()) {
480  case LinSol::KLU_SOLVER:
481  case LinSol::UMFPACK_SOLVER: {
483 
484  if (HP.IsKeyWord("no")) {
486 
487  } else if (HP.IsKeyWord("max") || HP.IsKeyWord("row" "max")) {
489 
490  } else if (HP.IsKeyWord("sum") || HP.IsKeyWord("row" "sum")) {
492  } else if (HP.IsKeyWord("column" "max")) {
494  } else if (HP.IsKeyWord("column" "sum")) {
496  } else if (HP.IsKeyWord("lapack")) {
498  } else if (HP.IsKeyWord("iterative")) {
500  } else if (HP.IsKeyWord("row" "max" "column" "max") || HP.IsKeyWord("sinkhorn" "knopp")) {
502  }
503 
504  if (HP.IsKeyWord("scale" "tolerance")) {
505  scale.dTol = HP.GetReal();
506  }
507 
508  if (HP.IsKeyWord("scale" "iterations")) {
509  scale.iMaxIter = HP.GetInt();
510  }
511 
512  if (HP.IsKeyWord("once")) {
514 
515  } else if (HP.IsKeyWord("always")) {
517 
518  } else if (HP.IsKeyWord("never")) {
520  }
521 
522  switch (scale.when) {
525  switch (scale.algorithm) {
527  scale.algorithm = SolutionManager::SCALEA_LAPACK; // Restore the original behavior for Naive
528  break;
529 
530  default:
531  // Use the value provided by the input file
532  ;
533  }
534  break;
535 
538  break;
539  }
540 
541  if (HP.IsKeyWord("verbose") && HP.GetYesNoOrBool()) {
543  }
544 
545  if (HP.IsKeyWord("warnings") && HP.GetYesNoOrBool()) {
547  }
548 
549  unsigned uCondFlag = 0;
550 
551  if (HP.IsKeyWord("print" "condition" "number")) {
552  if (HP.IsKeyWord("norm")) {
553  if (HP.IsKeyWord("inf")) {
555  } else {
556  const doublereal dNorm = HP.GetReal();
557 
558  if (dNorm != 1.) {
559  silent_cerr("Only one norm or infinity norm are supported for condition numbers at line " << HP.GetLineData() << std::endl);
561  }
562 
564  }
565  } else {
567  }
568  }
569 
570  if (uCondFlag != 0 && HP.GetYesNoOrBool()) {
571  scale.uFlags |= uCondFlag;
572  }
573 
574  if (!cs.SetScale(scale)) {
575  silent_cerr("Warning: Scale options are not available for "
576  << cs.GetSolverName()
577  << " at line "
578  << HP.GetLineData() << std::endl);
579  }
580 
581  } break;
582 
583  default:
584  pedantic_cerr("scale is meaningless for "
585  << currSolver.s_name
586  << " solver" << std::endl);
587  break;
588  }
589  }
590 
591  if (HP.IsKeyWord("max" "iterations")) {
592  if (!cs.SetMaxIterations(HP.GetInt())) {
593  silent_cerr("Warning: iterative refinement is not supported by " << cs.GetSolverName() << " at line " << HP.GetLineData() << std::endl);
594  }
595  }
596 
597  switch (cs.GetSolver()) {
600  silent_cout("warning: \"naive\" solver should be used with \"colamd\"" << std::endl);
601  }
602  break;
603 
604  // add more warnings...
605 
606  default:
607  break;
608  }
609 }
bool SetScale(const SolutionManager::ScaleOpt &scale)
Definition: linsol.cc:417
bool SetPivotFactor(const doublereal &d)
Definition: linsol.cc:372
unsigned GetSolverFlags(void) const
Definition: linsol.cc:257
bool MaskSolverFlags(unsigned f)
Definition: linsol.cc:303
#define MBDYN_EXCEPT_ARGS
Definition: except.h:63
virtual integer GetInt(integer iDefval=0)
Definition: parser.cc:1050
virtual HighParser::ErrOut GetLineData(void) const
Definition: parser.cc:681
bool SetNumThreads(unsigned nt)
Definition: linsol.cc:314
virtual bool GetYesNoOrBool(bool bDefval=false)
Definition: parser.cc:1038
const char *const s_name
Definition: linsol.h:85
bool SetSolver(SolverType t, unsigned f=SOLVER_FLAGS_NONE)
Definition: linsol.cc:181
const LinSol::solver_t solver[]
Definition: linsol.cc:58
unsigned s_default_flags
Definition: linsol.h:89
bool SetSolverFlags(unsigned f)
Definition: linsol.cc:281
bool SetMaxIterations(integer iMaxIter)
Definition: linsol.cc:440
bool SetWorkSpaceSize(integer)
Definition: linsol.cc:357
bool AddSolverFlags(unsigned f)
Definition: linsol.cc:292
virtual bool IsKeyWord(const char *sKeyWord)
Definition: parser.cc:910
doublereal s_drop_tolerance
Definition: linsol.h:91
doublereal s_pivot_factor
Definition: linsol.h:90
KeyWords
Definition: dataman4.cc:94
const char *const s_alias
Definition: linsol.h:86
bool SetBlockSize(unsigned bs)
Definition: linsol.cc:402
SolverType GetSolver(void) const
Definition: linsol.cc:175
virtual int GetWord(void)
Definition: parser.cc:1083
ScaleAlgorithm algorithm
Definition: solman.h:145
int get_nprocs(void)
Definition: get_nprocs.c:44
#define EMPTY
Definition: colamd.c:768
bool SetDropTolerance(const doublereal &d)
Definition: linsol.cc:384
const char *const GetSolverName(void) const
Definition: linsol.cc:269
double doublereal
Definition: colamd.c:52
long int integer
Definition: colamd.c:51
#define DEBUGLCOUT(level, msg)
Definition: myassert.h:244
unsigned s_flags
Definition: linsol.h:88
virtual doublereal GetReal(const doublereal &dDefval=0.0)
Definition: parser.cc:1056

Here is the call graph for this function:

std::ostream& RestartLinSol ( std::ostream &  out,
const LinSol cs 
)

Definition at line 611 of file readlinsol.cc.

References LinSol::dGetPivotFactor(), LinSol::GetBlockSize(), LinSol::GetNumThreads(), LinSol::GetSolver(), LinSol::GetSolverFlags(), LinSol::GetSolverName(), LinSol::iGetWorkSpaceSize(), LinSol::solver_t::s_default_flags, solver, LinSol::SOLVER_FLAGS_ALLOWS_CC, LinSol::SOLVER_FLAGS_ALLOWS_COLAMD, LinSol::SOLVER_FLAGS_ALLOWS_DIR, LinSol::SOLVER_FLAGS_ALLOWS_MAP, and LinSol::SOLVER_FLAGS_ALLOWS_MT_FCT.

Referenced by InverseSolver::Restart(), and Solver::Restart().

612 {
613  out << cs.GetSolverName();
614 
615  const LinSol::solver_t currSolver = ::solver[cs.GetSolver()];
616 
617  if (cs.GetSolverFlags() != currSolver.s_default_flags) {
618  unsigned f = cs.GetSolverFlags();
620  LinSol::SOLVER_FLAGS_ALLOWS_MAP) {
621  /*Map*/
622  out << ", map ";
623  }
624  if((f & LinSol::SOLVER_FLAGS_ALLOWS_CC) ==
625  LinSol::SOLVER_FLAGS_ALLOWS_CC) {
626  /*column compressed*/
627  out << ", cc ";
628  }
630  LinSol::SOLVER_FLAGS_ALLOWS_DIR) {
631  /*direct access*/
632  out << ", dir ";
633  }
635  LinSol::SOLVER_FLAGS_ALLOWS_COLAMD) {
636  /*colamd*/
637  out << ", colamd ";
638  }
639  unsigned nt = cs.GetNumThreads();
641  LinSol::SOLVER_FLAGS_ALLOWS_MT_FCT) &&
642  (nt != 1)) {
643  /*multi thread*/
644  out << ", mt , " << nt;
645  }
646  }
647  integer ws = cs.iGetWorkSpaceSize();
648  if(ws > 0) {
649  /*workspace size*/
650  out << ", workspace size, " << ws;
651  }
652  doublereal pf = cs.dGetPivotFactor();
653  if(pf != -1.) {
654  /*pivot factor*/
655  out << ", pivot factor, " << pf;
656  }
657  unsigned bs = cs.GetBlockSize();
658  if(bs != 0) {
659  /*block size*/
660  out << ", block size, " << bs ;
661  }
662  out << ";" << std::endl;
663  return out;
664 }
unsigned GetSolverFlags(void) const
Definition: linsol.cc:257
unsigned GetNumThreads(void) const
Definition: linsol.cc:333
const LinSol::solver_t solver[]
Definition: linsol.cc:58
unsigned s_default_flags
Definition: linsol.h:89
const doublereal & dGetPivotFactor(void) const
Definition: linsol.cc:345
SolverType GetSolver(void) const
Definition: linsol.cc:175
integer iGetWorkSpaceSize(void) const
Definition: linsol.cc:339
unsigned GetBlockSize(void) const
Definition: linsol.cc:396
const char *const GetSolverName(void) const
Definition: linsol.cc:269
double doublereal
Definition: colamd.c:52
long int integer
Definition: colamd.c:51

Here is the call graph for this function: