44 #include "MsqHessian.hpp"
45 #include "MsqTimer.hpp"
48 #ifdef MSQ_USE_OLD_C_HEADERS
54 #ifdef MSQ_USE_OLD_IO_HEADERS
55 # include <iostream.h>
63 origin_pd(0), mEntries(0), mRowStart(0), mColIndex(0),
64 mAccumulation(0), mAccumElemStart(0), mSize(0),
65 mPreconditioner(0), precondArraySize(0),
66 mR(0), mZ(0), mP(0), mW(0), cgArraySizes(0), maxCGiter(50)
70 MsqHessian::~MsqHessian()
76 delete[] mAccumulation;
77 delete[] mAccumElemStart;
79 delete[] mPreconditioner;
91 void MsqHessian::initialize(PatchData &pd, MsqError &err)
97 delete[] mAccumulation;
98 delete[] mAccumElemStart;
100 size_t num_vertices = pd.num_vertices();
101 size_t num_elements = pd.num_elements();
102 size_t const * vtx_list;
103 size_t e, r,
rs, re, c, cs, ce, nz, nnz, nve,
i,
j;
104 patchElemArray = pd.get_element_array(err);
MSQ_CHKERR(err);
106 if (num_vertices == 0) {
107 MSQ_SETERR( err )(
"No vertices in PatchData", MsqError::INVALID_ARG);
111 mSize = num_vertices;
116 size_t* col_start =
new size_t[num_vertices + 1];
117 mAccumElemStart =
new size_t[num_elements+1];
118 mAccumElemStart[0] = 0;
120 for (i = 0; i < num_vertices; ++
i) {
124 for (e = 0; e < num_elements; ++e) {
125 nve = patchElemArray[e].vertex_count();
126 vtx_list = patchElemArray[e].get_vertex_index_array();
127 mAccumElemStart[e+1] = mAccumElemStart[e] + (nve+1)*nve/2;
129 for (i = 0; i < nve; ++
i) {
132 for (j = i; j < nve; ++
j) {
146 for (i = 0; i < num_vertices; ++
i) {
155 int* row_instr =
new int[5*nz];
156 size_t* row_index =
new size_t[nz];
159 for (e = 0; e < num_elements; ++e) {
160 nve = patchElemArray[e].vertex_count();
161 vtx_list = patchElemArray[e].get_vertex_index_array();
163 for (i = 0; i < nve; ++
i) {
166 for (j = i; j < nve; ++
j) {
170 row_index[col_start[c]] = r;
171 row_instr[col_start[c]] = nz;
175 row_index[col_start[r]] = c;
176 row_instr[col_start[r]] = -nz;
185 for (i = num_vertices-1; i > 0; --
i) {
186 col_start[i+1] = col_start[
i];
188 col_start[1] = col_start[0];
208 size_t* row_start =
new size_t[num_vertices + 1];
210 for (i = 0; i < num_vertices; ++
i) {
214 for (i = 0; i < nz; ++
i) {
215 ++row_start[row_index[
i]];
219 for (i = 0; i < num_vertices; ++
i) {
228 size_t* col_index =
new size_t[nz];
229 int* col_instr =
new int[nz];
231 for (i = 0; i < num_vertices; ++
i) {
238 col_index[row_start[r]] =
i;
239 col_instr[row_start[r]] = row_instr[cs];
246 for (i = num_vertices-1; i > 0; --
i) {
247 row_start[i+1] = row_start[
i];
249 row_start[1] = row_start[0];
258 mRowStart = col_start;
259 mAccumulation = row_instr;
261 for (i = 0; i <= num_vertices; ++
i) {
266 for (i = 0; i < num_vertices; ++
i) {
272 if (c != col_index[rs]) {
280 if (col_instr[rs] >= 0) {
281 mAccumulation[col_instr[
rs]] = nnz - 1;
284 mAccumulation[-col_instr[
rs]] = 1 - nnz;
292 for (i = 0; i < num_vertices; ++
i) {
303 mColIndex =
new size_t[nnz];
305 for (i = 0; i < num_vertices; ++
i) {
311 if (c != col_index[rs]) {
315 mColIndex[mRowStart[
i]] = c;
322 for (i = num_vertices-1; i > 0; --
i) {
323 mRowStart[i+1] = mRowStart[
i];
325 mRowStart[1] = mRowStart[0];
331 mEntries =
new Matrix3D[nnz];
332 for (i=0;i<nnz;++
i) mEntries[i] = 0.;
341 void MsqHessian::get_diagonal_blocks(msq_std::vector<Matrix3D> &diag,
345 if (diag.size() != size()) {
346 diag.reserve(size());
349 for (
size_t i=0; i<size(); ++
i) {
350 diag[
i] = mEntries[mRowStart[
i]];
358 void MsqHessian::compute_preconditioner(MsqError &)
361 if (mSize > precondArraySize || mSize < precondArraySize/10 ) {
362 delete[] mPreconditioner;
363 mPreconditioner =
new Matrix3D[mSize];
366 Matrix3D* diag_block;
371 for (m=0; m<mSize; ++m) {
372 diag_block = mEntries + mRowStart[m];
374 #if DIAGONAL_PRECONDITIONER
376 sum = (*diag_block)[0][0] + (*diag_block)[1][1] + (*diag_block)[2][2];
383 mPreconditioner[m][0][0] = inv_sum;
384 mPreconditioner[m][1][1] = inv_sum;
385 mPreconditioner[m][2][2] = inv_sum;
395 if ((*diag_block)[0][0] == 0.0) {
400 sum = (*diag_block)[0][0] + (*diag_block)[1][1] + (*diag_block)[2][2];
404 mPreconditioner[m][0][0] = sum;
405 mPreconditioner[m][0][1] = 0.0;
406 mPreconditioner[m][0][2] = 0.0;
407 mPreconditioner[m][1][1] = sum;
408 mPreconditioner[m][1][2] = 0.0;
409 mPreconditioner[m][2][2] = sum;
412 mPreconditioner[m][0][0] = 1.0 / (*diag_block)[0][0];
413 mPreconditioner[m][0][1] = (*diag_block)[0][1] * mPreconditioner[m][0][0];
414 mPreconditioner[m][0][2] = (*diag_block)[0][2] * mPreconditioner[m][0][0];
416 mPreconditioner[m][1][1] =
417 1.0 / ((*diag_block)[1][1] -
418 (*diag_block)[0][1] * mPreconditioner[m][0][1]);
420 tmp = (*diag_block)[1][2] -
421 (*diag_block)[0][2] * mPreconditioner[m][0][1];
423 mPreconditioner[m][1][2] = mPreconditioner[m][1][1] * tmp;
425 mPreconditioner[m][2][2] =
426 1.0 / ((*diag_block)[2][2] -
427 (*diag_block)[0][2]*mPreconditioner[m][0][2] -
428 mPreconditioner[m][1][2]*tmp);
445 if (mSize > cgArraySizes || mSize < cgArraySizes/10 ) {
454 cgArraySizes = mSize;
458 double alpha_, alpha, beta;
459 double cg_tol = 1e-2;
460 double norm_g =
length(b, mSize);
461 double norm_r = norm_g;
465 this->compute_preconditioner(err);
MSQ_CHKERR(err);
467 for (i=0; i<mSize; ++
i) x[i] = 0. ;
468 for (i=0; i<mSize; ++
i) mR[i] = -b[i] ;
471 this->apply_preconditioner(mZ, mR, err);
472 for (i=0; i<mSize; ++
i) mP[i] = mZ[i] ;
473 rzm1 =
inner(mZ,mR,mSize);
476 while ((norm_r > norm_g) && (cg_iter < maxCGiter)) {
479 axpy(mW, mSize, *
this, mP, mSize, 0,0,err);
481 alpha_ =
inner(mP,mW,mSize);
484 for (i=0; i<mSize; ++
i) x[i] += mP[i];
489 alpha = rzm1 / alpha_;
491 for (i=0; i<mSize; ++
i) x[i] += alpha*mP[i];
492 for (i=0; i<mSize; ++
i) mR[i] -= alpha*mW[i];
493 norm_r =
length(mR, mSize);
495 this->apply_preconditioner(mZ, mR, err);
498 rzm1 =
inner(mZ,mR,mSize);
500 for (i=0; i<mSize; ++
i) mP[i] = mZ[i] + beta*mP[i];
507 msq_stdio::ostream&
operator<<(msq_stdio::ostream &
s,
const MsqHessian &h)
510 s <<
"MsqHessian of size: " << h.mSize <<
"x"<< h.mSize <<
"\n";
511 for (i=0; i<h.mSize; ++
i) {
512 s <<
" ROW " << i <<
" ------------------------\n";
513 for (j=h.mRowStart[i]; j<h.mRowStart[i+1]; ++j) {
514 s <<
" column " << h.mColIndex[
j] <<
" ----\n";
subroutine rs(nm, n, a, w, matz, z, fv1, fv2, ierr)
void axpy(Vector3D res[], size_t size_r, const MsqHessian &H, const Vector3D x[], size_t size_x, const Vector3D y[], size_t size_y, MsqError &)
double length(Vector3D *const v, int n)
NVec< 3, double > Vector3D
#define MSQ_CHKERR(err)
Mesquite's Error Checking macro.
#define MSQ_SETERR(err)
Macro to set error - use err.clear() to clear.
#define MSQ_FUNCTION_TIMER(NAME)
std::ostream & operator<<(std::ostream &os, const COM_exception &ex)
Print out a given exception.
double inner(const Vector3D lhs[], const Vector3D rhs[], int n)