19#include "MyIncompleteChol.h"
22MyIncompleteChol::MyIncompleteChol(
const Epetra_Comm &_Comm,
const Epetra_Operator *KK,
23 double tolDrop,
int fillLevel,
const Epetra_MultiVector *Z)
33 leftProjection(false),
34 rightProjection(false)
38 Epetra_CrsMatrix *Kcrs =
dynamic_cast<Epetra_CrsMatrix*
>(
const_cast<Epetra_Operator*
>(K));
40 if (MyComm.MyPID() == 0) {
42 cerr <<
" !!! For incomplete Cholesky preconditioner, ";
43 cerr <<
"the matrix must be 'Epetra_CrsMatrix' object !!!\n";
49 Prec =
new Ifpack_CrsIct(*Kcrs, dropTol, Kcrs->GlobalMaxNumEntries() + lFill);
50 Prec->InitValues(*Kcrs);
54 QtQ =
new double[Q->NumVectors()*Q->NumVectors()];
55 double *work =
new double[Q->NumVectors()*Q->NumVectors()];
56 callBLAS.GEMM(
'T',
'N', Q->NumVectors(), Q->NumVectors(), Q->MyLength(),
57 1.0, Q->Values(), Q->MyLength(), Q->Values(), Q->MyLength(),
58 0.0, work, Q->NumVectors());
59 MyComm.SumAll(work, QtQ, Q->NumVectors()*Q->NumVectors());
62 callLAPACK.POTRF(
'U', Q->NumVectors(), QtQ, Q->NumVectors(), &info);
65 cerr <<
" !!! In incomplete Cholesky, the null space vectors are linearly dependent !!!\n";
75MyIncompleteChol::~MyIncompleteChol() {
90int MyIncompleteChol::Apply(
const Epetra_MultiVector& X, Epetra_MultiVector& Y)
const {
99int MyIncompleteChol::ApplyInverse(
const Epetra_MultiVector& X, Epetra_MultiVector& Y)
const {
101 int xcol = X.NumVectors();
103 if (Y.NumVectors() < xcol)
106 double *valX = (rightProjection ==
true) ?
new double[xcol*X.MyLength()] : X.Values();
108 Epetra_MultiVector X2(View, X.Map(), valX, X.MyLength(), xcol);
110 if ((rightProjection ==
true) && (Q)) {
111 int qcol = Q->NumVectors();
112 double *work =
new double[2*qcol*xcol];
113 memcpy(X2.Values(), X.Values(), xcol*X.MyLength()*
sizeof(
double));
114 callBLAS.GEMM(
'T',
'N', qcol, xcol, Q->MyLength(), 1.0, Q->Values(), Q->MyLength(),
115 X2.Values(), X2.MyLength(), 0.0, work + qcol*xcol, qcol);
116 MyComm.SumAll(work + qcol*xcol, work, qcol*xcol);
118 callLAPACK.POTRS(
'U', qcol, xcol, QtQ, qcol, work, qcol, &info);
119 callBLAS.GEMM(
'N',
'N', X2.MyLength(), xcol, qcol, -1.0, Q->Values(), Q->MyLength(),
120 work, qcol, 1.0, X2.Values(), X2.MyLength());
124 Prec->ApplyInverse(X2, Y);
129 if ((leftProjection ==
true) && (Q)) {
130 int qcol = Q->NumVectors();
131 double *work =
new double[2*qcol*xcol];
132 callBLAS.GEMM(
'T',
'N', qcol, xcol, Q->MyLength(), 1.0, Q->Values(), Q->MyLength(),
133 Y.Values(), Y.MyLength(), 0.0, work + qcol*xcol, qcol);
134 MyComm.SumAll(work + qcol*xcol, work, qcol*xcol);
136 callLAPACK.POTRS(
'U', qcol, xcol, QtQ, qcol, work, qcol, &info);
137 callBLAS.GEMM(
'N',
'N', Y.MyLength(), xcol, qcol, -1.0, Q->Values(), Q->MyLength(),
138 work, qcol, 1.0, Y.Values(), Y.MyLength());