Anasazi Version of the Day
Loading...
Searching...
No Matches
MyIncompleteChol.cpp
1// @HEADER
2// *****************************************************************************
3// Anasazi: Block Eigensolvers Package
4//
5// Copyright 2004 NTESS and the Anasazi contributors.
6// SPDX-License-Identifier: BSD-3-Clause
7// *****************************************************************************
8// @HEADER
9
10// This software is a result of the research described in the report
11//
12// "A comparison of algorithms for modal analysis in the absence
13// of a sparse direct method", P. Arbenz, R. Lehoucq, and U. Hetmaniuk,
14// Sandia National Laboratories, Technical report SAND2003-1028J.
15//
16// It is based on the Epetra, AztecOO, and ML packages defined in the Trilinos
17// framework ( http://trilinos.org/ ).
18
19#include "MyIncompleteChol.h"
20
21
22MyIncompleteChol::MyIncompleteChol(const Epetra_Comm &_Comm, const Epetra_Operator *KK,
23 double tolDrop, int fillLevel, const Epetra_MultiVector *Z)
24 : MyComm(_Comm),
25 callBLAS(),
26 callLAPACK(),
27 K(KK),
28 Prec(0),
29 dropTol(tolDrop),
30 lFill(fillLevel),
31 Q(Z),
32 QtQ(0),
33 leftProjection(false),
34 rightProjection(false)
35 {
36
37 // Note the constness is cast away for ML
38 Epetra_CrsMatrix *Kcrs = dynamic_cast<Epetra_CrsMatrix*>(const_cast<Epetra_Operator*>(K));
39 if (Kcrs == 0) {
40 if (MyComm.MyPID() == 0) {
41 cerr << endl;
42 cerr << " !!! For incomplete Cholesky preconditioner, ";
43 cerr << "the matrix must be 'Epetra_CrsMatrix' object !!!\n";
44 cerr << endl;
45 }
46 return;
47 }
48
49 Prec = new Ifpack_CrsIct(*Kcrs, dropTol, Kcrs->GlobalMaxNumEntries() + lFill);
50 Prec->InitValues(*Kcrs);
51 Prec->Factor();
52
53 if (Q) {
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());
60 delete[] work;
61 int info = 0;
62 callLAPACK.POTRF('U', Q->NumVectors(), QtQ, Q->NumVectors(), &info);
63 if (info) {
64 cerr << endl;
65 cerr << " !!! In incomplete Cholesky, the null space vectors are linearly dependent !!!\n";
66 cerr << endl;
67 delete[] QtQ;
68 QtQ = 0;
69 }
70 }
71
72}
73
74
75MyIncompleteChol::~MyIncompleteChol() {
76
77 if (Prec) {
78 delete Prec;
79 Prec = 0;
80 }
81
82 if (QtQ) {
83 delete[] QtQ;
84 QtQ = 0;
85 }
86
87}
88
89
90int MyIncompleteChol::Apply(const Epetra_MultiVector& X, Epetra_MultiVector& Y) const {
91
92 Y.PutScalar(0.0);
93
94 return -1;
95
96}
97
98
99int MyIncompleteChol::ApplyInverse(const Epetra_MultiVector& X, Epetra_MultiVector& Y) const {
100
101 int xcol = X.NumVectors();
102
103 if (Y.NumVectors() < xcol)
104 return -1;
105
106 double *valX = (rightProjection == true) ? new double[xcol*X.MyLength()] : X.Values();
107
108 Epetra_MultiVector X2(View, X.Map(), valX, X.MyLength(), xcol);
109
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);
117 int info = 0;
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());
121 delete[] work;
122 }
123
124 Prec->ApplyInverse(X2, Y);
125
126 if (rightProjection)
127 delete[] valX;
128
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);
135 int info = 0;
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());
139 delete[] work;
140 }
141
142 return 0;
143
144}
145
146