1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171 172 173 174 175 176 177 178 179 180 181 182 183 184 185 186 187 188 189 190 191 192 193
|
/******************************************************************************
* SOFA, Simulation Open-Framework Architecture, version 1.0 beta 4 *
* (c) 2006-2009 MGH, INRIA, USTL, UJF, CNRS *
* *
* This library is free software; you can redistribute it and/or modify it *
* under the terms of the GNU Lesser General Public License as published by *
* the Free Software Foundation; either version 2.1 of the License, or (at *
* your option) any later version. *
* *
* This library is distributed in the hope that it will be useful, but WITHOUT *
* ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or *
* FITNESS FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License *
* for more details. *
* *
* You should have received a copy of the GNU Lesser General Public License *
* along with this library; if not, write to the Free Software Foundation, *
* Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA. *
*******************************************************************************
* SOFA :: Modules *
* *
* Authors: The SOFA Team and external contributors (see Authors.txt) *
* *
* Contact information: contact@sofa-framework.org *
******************************************************************************/
// Author: Hadrien Courtecuisse
//
// Copyright: See COPYING file that comes with this distribution
#include <sofa/component/linearsolver/SparseCholeskySolver.h>
#include <sofa/core/ObjectFactory.h>
#include <sofa/helper/system/thread/CTime.h>
#include <iostream>
#include <math.h>
namespace sofa {
namespace component {
namespace linearsolver {
using namespace sofa::defaulttype;
using namespace sofa::core::componentmodel::behavior;
using namespace sofa::simulation;
using namespace sofa::core::objectmodel;
using sofa::helper::system::thread::CTime;
using sofa::helper::system::thread::ctime_t;
using std::cerr;
using std::endl;
template<class TMatrix, class TVector>
SparseCholeskySolver<TMatrix,TVector>::SparseCholeskySolver()
: f_verbose( initData(&f_verbose,false,"verbose","Dump system state at each iteration") )
, f_graph( initData(&f_graph,"graph","Graph of residuals at each iteration") )
, S(NULL), N(NULL), tmp(NULL)
{
f_graph.setWidget("graph");
f_graph.setReadOnly(true);
}
template<class TMatrix, class TVector>
SparseCholeskySolver<TMatrix,TVector>::~SparseCholeskySolver() {
if (S) cs_sfree (S);
if (N) cs_nfree (N);
if (tmp) cs_free (tmp);
}
template<class TMatrix, class TVector>
void SparseCholeskySolver<TMatrix,TVector>::solve (Matrix& /*M*/, Vector& z, Vector& r) {
int n = A.n;
cs_ipvec (n, S->Pinv, r.ptr(), tmp); //x = P*b
cs_lsolve (N->L, tmp); //x = L\x
cs_ltsolve (N->L, tmp); //x = L'\x/
cs_pvec (n, S->Pinv, tmp, z.ptr()); //b = P'*x
}
template<class TMatrix, class TVector>
void SparseCholeskySolver<TMatrix,TVector>::invert(Matrix& M) {
int order = -1; //?????
if (S) cs_sfree(S);
if (N) cs_nfree(N);
if (tmp) cs_free(tmp);
M.compress();
//remplir A avec M
A.nzmax = M.getColsValue().size(); // maximum number of entries
A.m = M.rowBSize(); // number of rows
A.n = M.colBSize(); // number of columns
A_p = M.getRowBegin();
A.p = (int *) &(A_p[0]); // column pointers (size n+1) or col indices (size nzmax)
A_i = M.getColsIndex();
A.i = (int *) &(A_i[0]); // row indices, size nzmax
A_x = M.getColsValue();
A.x = (double *) &(A_x[0]); // numerical values, size nzmax
A.nz = -1; // # of entries in triplet matrix, -1 for compressed-col
cs_dropzeros( &A );
//M.check_matrix();
//CompressedRowSparseMatrix<double>::check_matrix(-1 /*A.nzmax*/,A.m,A.n,A.p,A.i,A.x);
//sout << "diag =";
//for (int i=0;i<A.n;++i) sout << " " << M.element(i,i);
//sout << sendl;
//sout << "SparseCholeskySolver: start factorization, n = " << A.n << " nnz = " << A.p[A.n] << sendl;
tmp = (double *) cs_malloc (A.n, sizeof (double)) ;
S = cs_schol (&A, order) ; /* ordering and symbolic analysis */
N = cs_chol (&A, S) ; /* numeric Cholesky factorization */
//sout << "SparseCholeskySolver: factorization complete, nnz = " << N->L->p[N->L->n] << sendl;
}
template<class TMatrix, class TVector>
bool SparseCholeskySolver<TMatrix,TVector>::readFile(std::istream& in) {
std::cout << "Read SparseCholeskySolver" << std::endl;
/*
std::string s = "SparseCholeskySolver\n";
//in >> ss;
in >> s;
if (s.compare("SparseCholeskySolver\n")) {
std::cout << "File not contain a SparseLDLSolver" << std::endl;
return false;
}
in >> A.n;
in >> A_x;
in >> A_i;
in >> A_p;
in >> D;
in >> Parent;
in >> Lnz;
in >> Flag;
in >> Pattern;
in >> Lp;
in >> Lx;
in >> Li;
return true;
*/
return false;
}
template<class TMatrix, class TVector>
bool SparseCholeskySolver<TMatrix,TVector>::writeFile(std::ostream& out) {
std::string s = "SparseCholeskySolver\n";
out << s;
/*
out << A.n;
FullVector<double> v;
for (int i=0;i<n;i++) v[i]
out << A_i;
out << A_p;
out << D;
out << Y;
out << Parent;
out << Lnz;
out << Flag;
out << Pattern;
out << Lp;
out << Lx;
out << Li;
return true;
*/
return false;
}
SOFA_DECL_CLASS(SparseCholeskySolver)
int SparseCholeskySolverClass = core::RegisterObject("Linear system solver using the conjugate gradient iterative algorithm")
.add< SparseCholeskySolver< CompressedRowSparseMatrix<double>,FullVector<double> > >(true)
.addAlias("SparseCholeskySolverAlias")
;
template class SOFA_COMPONENT_LINEARSOLVER_API SparseCholeskySolver< CompressedRowSparseMatrix<double>,FullVector<double> >;
} // namespace linearsolver
} // namespace component
} // namespace sofa
|