MBDyn-1.7.3
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups
kluwrap.cc
Go to the documentation of this file.
1 /* $Header: /var/cvs/mbdyn/mbdyn/mbdyn-1.0/libraries/libmbwrap/kluwrap.cc,v 1.8 2014/01/13 08:59:07 masarati Exp $ */
2 /*
3  * HmFe (C) is a FEM analysis code.
4  *
5  * Copyright (C) 1996-2006
6  *
7  * Marco Morandini <morandini@aero.polimi.it>
8  *
9  * Dipartimento di Ingegneria Aerospaziale - Politecnico di Milano
10  * via La Masa, 34 - 20156 Milano, Italy
11  * http://www.aero.polimi.it
12  *
13  * Changing this copyright notice is forbidden.
14  * This program is distributed in the hope that it will be useful,
15  * but WITHOUT ANY WARRANTY; without even the implied warranty of
16  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.
17  *
18  */
19 /* December 2001
20  * Modified to add a Sparse matrix in row form and to implement methods
21  * to be used in the parallel MBDyn Solver.
22  *
23  * Copyright (C) 2001-2005
24  *
25  * Giuseppe Quaranta <quaranta@aero.polimi.it>
26  *
27  * Dipartimento di Ingegneria Aerospaziale - Politecnico di Milano
28  * via La Masa, 34 - 20156 Milano, Italy
29  * http://www.aero.polimi.it
30  *
31  */
32 /*
33  * MBDyn (C) is a multibody analysis code.
34  * http://www.mbdyn.org
35  *
36  * Copyright (C) 1996-2005
37  *
38  * Pierangelo Masarati <masarati@aero.polimi.it>
39  * Paolo Mantegazza <mantegazza@aero.polimi.it>
40  *
41  * Dipartimento di Ingegneria Aerospaziale - Politecnico di Milano
42  * via La Masa, 34 - 20156 Milano, Italy
43  * http://www.aero.polimi.it
44  *
45  * Changing this copyright notice is forbidden.
46  *
47  * This program is free software; you can redistribute it and/or modify
48  * it under the terms of the GNU General Public License as published by
49  * the Free Software Foundation (version 2 of the License).
50  *
51  *
52  * This program is distributed in the hope that it will be useful,
53  * but WITHOUT ANY WARRANTY; without even the implied warranty of
54  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
55  * GNU General Public License for more details.
56  *
57  * You should have received a copy of the GNU General Public License
58  * along with this program; if not, write to the Free Software
59  * Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
60  */
61 
62 /*
63  * Umfpack is used by permission; please read its Copyright,
64  * License and Availability note.
65  */
66 
67 #include "mbconfig.h" /* This goes first in every *.c,*.cc file */
68 
69 #ifdef USE_KLU
70 #include <cassert>
71 #include "solman.h"
72 #include "spmapmh.h"
73 #include "ccmh.h"
74 #include "dirccmh.h"
75 #include "kluwrap.h"
76 #include "dgeequ.h"
77 
78 /* KLUSolver - begin */
79 
80 KLUSolver::KLUSolver(const integer &size, const doublereal &dPivot, Scale scale)
81 : LinearSolver(0),
82 iSize(size),
83 Axp(0),
84 Aip(0),
85 App(0),
86 Symbolic(0),
87 Numeric(0)
88 {
89  klu_defaults(&Control);
90 
91  if (dPivot != -1. && (dPivot >= 0. && dPivot <= 1.)) {
92  /*
93  * 1.0: true partial pivoting
94  * 0.0: treated as 1.0
95  *
96  * default: 0.1
97  */
98  Control.tol = dPivot;
99  }
100  // Control.ordering = 1; /* colamd */
101  Control.ordering = 0; /* amd */
102 
103  switch (scale) {
104  case SCALE_UNDEF:
105  // use the default value provided by KLU
106  break;
107  default:
108  Control.scale = scale;
109  }
110 }
111 
112 KLUSolver::~KLUSolver(void)
113 {
114  if (Symbolic) {
115  klu_free_symbolic(&Symbolic, &Control);
116  }
117  ASSERT(Symbolic == 0);
118 
119  if (Numeric) {
120  klu_free_numeric(&Numeric, &Control);
121  }
122  ASSERT(Numeric == 0);
123 }
124 
125 void
126 KLUSolver::Reset(void)
127 {
128  if (Symbolic) {
129  // FIXME: This code might be an performance issue!
130  // However if we do not reset the symbolic object
131  // the code will fail if zero entries become nonzero
132  // during simulation.
133  klu_free_symbolic(&Symbolic, &Control);
134  ASSERT(Symbolic == 0);
135  }
136 
137  if (Numeric) {
138  klu_free_numeric(&Numeric, &Control);
139  ASSERT(Numeric == 0);
140  }
141 
142  bHasBeenReset = true;
143 }
144 
145 void
146 KLUSolver::Solve(void) const
147 {
148  if (bHasBeenReset) {
149  const_cast<KLUSolver *>(this)->Factor();
150  bHasBeenReset = false;
151  }
152 
153 
154  /*
155  * NOTE: Axp, Aip, App should have been set by * MakeCompactForm()
156  */
157 
158  int ok = klu_solve(Symbolic, Numeric, iSize, 1, LinearSolver::pdRhs,
159  &Control);
160  if (!ok || Control.status != KLU_OK) {
161  silent_cerr("KLUWRAP_solve failed" << std::endl);
162 
163  /* de-allocate memory */
164  if (Numeric) {
165  klu_free_numeric(&Numeric, &Control);
166  }
167  ASSERT(Numeric == 0);
168 
170  }
171 
172 }
173 
174 void
175 KLUSolver::Factor(void)
176 {
177 
178  /*
179  * NOTE: Axp, Aip, App should have been set by * MakeCompactForm()
180  */
181 
182  if (Symbolic == 0 && !bPrepareSymbolic()) {
184  }
185 
186  if (Numeric != 0) {
187  int ok;
188  ok = klu_refactor(App, Aip, Axp, Symbolic,
189  Numeric, &Control);
190  if (!ok) {
191  // FIXME: might be too late!
192  klu_free_numeric(&Numeric, &Control);
193  klu_free_symbolic(&Symbolic, &Control);
194  if (!bPrepareSymbolic()) {
196  }
197  }
198  }
199 
200  if (Numeric == 0) {
201  Numeric = klu_factor(App, Aip, Axp, Symbolic,
202  &Control);
203  }
204 
205  if (Control.status != KLU_OK) {
206  silent_cerr("KLUWRAP_numeric failed" << std::endl);
207 
208  /* de-allocate memory */
209  if (Symbolic) {
210  klu_free_symbolic(&Symbolic, &Control);
211  }
212  if (Numeric) {
213  klu_free_numeric(&Numeric, &Control);
214  }
215  ASSERT(Symbolic == 0);
216  ASSERT(Numeric == 0);
217 
219  }
220 }
221 
222 void
223 KLUSolver::MakeCompactForm(SparseMatrixHandler& mh,
224  std::vector<doublereal>& Ax,
225  std::vector<integer>& Ai,
226  std::vector<integer>& Ac,
227  std::vector<integer>& Ap) const
228 {
229  if (!bHasBeenReset) {
230  return;
231  }
232 
233  mh.MakeCompressedColumnForm(Ax, Ai, Ap, 0);
234 
235  Axp = &(Ax[0]);
236  Aip = &(Ai[0]);
237  App = &(Ap[0]);
238 }
239 
240 bool
241 KLUSolver::bPrepareSymbolic(void)
242 {
243  Symbolic = klu_analyze(iSize, App, Aip, &Control);
244  if (Control.status != KLU_OK) {
245  silent_cerr("KLUWRAP_symbolic failed" << std::endl);
246 
247  /* de-allocate memory */
248  if (Symbolic) {
249  klu_free_symbolic(&Symbolic, &Control);
250  }
251  ASSERT(Symbolic == 0);
252 
253  return false;
254  }
255 
256  return true;
257 }
258 
259 bool KLUSolver::bGetConditionNumber(doublereal& dCond)
260 {
261  if (Symbolic == 0 || Numeric == 0) { // FIXME: Is it really safe?
262  return false;
263  }
264 
265  const bool bOK = klu_condest(App, Axp, Symbolic, Numeric, &Control);
266 
267  dCond = Control.condest;
268 
269  return bOK;
270 }
271 
272 /* KLUSolver - end */
273 
274 /* KLUSparseSolutionManager - begin */
275 
276 KLUSparseSolutionManager::KLUSparseSolutionManager(integer Dim,
277  doublereal dPivot,
278  const ScaleOpt& s)
279 : A(Dim),
280 b(Dim),
281 bVH(Dim, &b[0]),
282 scale(s),
283 pMatScale(0)
284 {
285  KLUSolver::Scale kscale = KLUSolver::SCALE_UNDEF;
286 
287  switch (scale.algorithm) {
288  case SCALEA_UNDEF:
289  scale.when = SCALEW_NEVER;
290  break;
291 
292  case SCALEA_NONE:
293  kscale = KLUSolver::SCALE_NONE;
294  scale.when = SCALEW_NEVER;
295  break;
296 
297  case SCALEA_ROW_MAX:
298  kscale = KLUSolver::SCALE_MAX;
299  scale.when = SCALEW_NEVER; // Do not scale twice! Use built in scaling from KLU
300  break;
301 
302  case SCALEA_ROW_SUM:
303  kscale = KLUSolver::SCALE_SUM;
304  scale.when = SCALEW_NEVER; // Do not scale twice! Use built in scaling from KLU
305  break;
306 
307  default:
308  // Allocate MatrixScale<T> on demand
309  kscale = KLUSolver::SCALE_NONE; // Do not scale twice!
310  }
311 
312  SAFENEWWITHCONSTRUCTOR(pLS, KLUSolver,
313  KLUSolver(Dim, dPivot, kscale));
314 
315  (void)pLS->pdSetResVec(&b[0]);
316  (void)pLS->pdSetSolVec(&b[0]);
317  pLS->SetSolutionManager(this);
318 }
319 
320 
321 KLUSparseSolutionManager::~KLUSparseSolutionManager(void)
322 {
323  SAFEDELETE(pMatScale);
324 }
325 
326 void
327 KLUSparseSolutionManager::MatrReset(void)
328 {
329  pLS->Reset();
330 }
331 
332 void
333 KLUSparseSolutionManager::MakeCompressedColumnForm(void)
334 {
335  ScaleMatrixAndRightHandSide(A);
336 
337  pLS->MakeCompactForm(A, Ax, Ai, Adummy, Ap);
338 }
339 
340 template <typename MH>
341 void KLUSparseSolutionManager::ScaleMatrixAndRightHandSide(MH& mh)
342 {
343  if (scale.when != SCALEW_NEVER) {
344  MatrixScale<MH>& rMatScale = GetMatrixScale<MH>();
345 
346  if (pLS->bReset()) {
347  if (!rMatScale.bGetInitialized()
348  || scale.when == SolutionManager::SCALEW_ALWAYS) {
349  // (re)compute
350  rMatScale.ComputeScaleFactors(mh);
351  }
352  // in any case scale matrix and right-hand-side
353  rMatScale.ScaleMatrix(mh);
354 
355  if (silent_err) {
356  rMatScale.Report(std::cerr);
357  }
358  }
359 
360  rMatScale.ScaleRightHandSide(bVH);
361  }
362 }
363 
364 template <typename MH>
365 MatrixScale<MH>& KLUSparseSolutionManager::GetMatrixScale()
366 {
367  if (pMatScale == 0) {
368  pMatScale = MatrixScale<MH>::Allocate(scale);
369  }
370 
371  // Will throw std::bad_cast if the type does not match
372  return dynamic_cast<MatrixScale<MH>&>(*pMatScale);
373 }
374 
375 void KLUSparseSolutionManager::ScaleSolution(void)
376 {
377  if (scale.when != SCALEW_NEVER) {
378  ASSERT(pMatScale != 0);
379  // scale solution
380  pMatScale->ScaleSolution(bVH);
381  }
382 }
383 
384 /* Risolve il sistema Fattorizzazione + Backward Substitution */
385 void
386 KLUSparseSolutionManager::Solve(void)
387 {
388  MakeCompressedColumnForm();
389 
390  pLS->Solve();
391 
392  ScaleSolution();
393 }
394 
395 /* Rende disponibile l'handler per la matrice */
397 KLUSparseSolutionManager::pMatHdl(void) const
398 {
399  return &A;
400 }
401 
402 /* Rende disponibile l'handler per il termine noto */
404 KLUSparseSolutionManager::pResHdl(void) const
405 {
406  return &bVH;
407 }
408 
409 /* Rende disponibile l'handler per la soluzione */
411 KLUSparseSolutionManager::pSolHdl(void) const
412 {
413  return &bVH;
414 }
415 
416 /* KLUSparseSolutionManager - end */
417 
418 template <class CC>
419 KLUSparseCCSolutionManager<CC>::KLUSparseCCSolutionManager(integer Dim,
420  doublereal dPivot,
421  const ScaleOpt& scale)
422 : KLUSparseSolutionManager(Dim, dPivot, scale),
423 CCReady(false),
424 Ac(0)
425 {
426  NO_OP;
427 }
428 
429 template <class CC>
430 KLUSparseCCSolutionManager<CC>::~KLUSparseCCSolutionManager(void)
431 {
432  if (Ac) {
433  SAFEDELETE(Ac);
434  }
435 }
436 
437 template <class CC>
438 void
439 KLUSparseCCSolutionManager<CC>::MatrReset(void)
440 {
441  pLS->Reset();
442 }
443 
444 template <class CC>
445 void
446 KLUSparseCCSolutionManager<CC>::MakeCompressedColumnForm(void)
447 {
448  if (!CCReady) {
449  pLS->MakeCompactForm(A, Ax, Ai, Adummy, Ap);
450 
451  if (Ac == 0) {
452  SAFENEWWITHCONSTRUCTOR(Ac, CC, CC(Ax, Ai, Ap));
453  }
454 
455  CCReady = true;
456  }
457 
458  ScaleMatrixAndRightHandSide(*Ac);
459 
460 }
461 
462 /* Inizializzatore "speciale" */
463 template <class CC>
464 void
465 KLUSparseCCSolutionManager<CC>::MatrInitialize()
466 {
467  CCReady = false;
468 
469  if (Ac) {
470  // If a DirCColMatrixHandler is in use and matrix scaling is enabled
471  // an uncaught exception (MatrixHandler::ErrRebuildMatrix) will be thrown
472  // if zero entries in the matrix become nonzero.
473  // For that reason we have to reinitialize Ac!
474  SAFEDELETE(Ac);
475  Ac = 0;
476  }
477 
478  MatrReset();
479 }
480 
481 /* Rende disponibile l'handler per la matrice */
482 template <class CC>
484 KLUSparseCCSolutionManager<CC>::pMatHdl(void) const
485 {
486  if (!CCReady) {
487  return &A;
488  }
489 
490  ASSERT(Ac != 0);
491  return Ac;
492 }
493 
494 template class KLUSparseCCSolutionManager<CColMatrixHandler<0> >;
495 template class KLUSparseCCSolutionManager<DirCColMatrixHandler<0> >;
496 /* KLUSparseCCSolutionManager - end */
497 
498 #endif /* USE_KLU */
499 
virtual integer MakeCompressedColumnForm(doublereal *const Ax, integer *const Ai, integer *const Ap, int offset=0) const =0
#define MBDYN_EXCEPT_ARGS
Definition: except.h:63
bool bGetInitialized() const
Definition: dgeequ.h:59
bool ComputeScaleFactors(const T &mh)
Definition: dgeequ.h:305
#define NO_OP
Definition: myassert.h:74
std::ostream & Report(std::ostream &os) const
Definition: dgeequ.cc:52
void Reset(scalar_func_type &d)
Definition: gradient.h:2836
Definition: mbdyn.h:76
VectorHandler & ScaleRightHandSide(VectorHandler &bVH) const
Definition: dgeequ.h:216
#define ASSERT(expression)
Definition: colamd.c:977
#define SAFENEWWITHCONSTRUCTOR(pnt, item, constructor)
Definition: mynewmem.h:698
static MatrixScale< T > * Allocate(const SolutionManager::ScaleOpt &scale)
Definition: dgeequ.h:313
doublereal * pdRhs
Definition: ls.h:74
double doublereal
Definition: colamd.c:52
long int integer
Definition: colamd.c:51
#define SAFEDELETE(pnt)
Definition: mynewmem.h:710
T & ScaleMatrix(T &mh) const
Definition: dgeequ.h:274