51#include <Teuchos_ParameterList.hpp>
75 :
A_(FactoredMatrix.
A_),
122 bool cerr_warning_if_unused)
147 int * InI=0, * LI=0, * UI = 0;
148 double * InV=0, * LV=0, * UV = 0;
149 int NumIn, NumL, NumU;
151 int NumNonzeroDiags = 0;
155 if (
Graph_.LevelOverlap()>0 &&
Graph_.L_Graph().DomainMap().DistributedGlobal()) {
164 InI =
new int[MaxNumEntries];
165 UI =
new int[MaxNumEntries];
166 InV =
new double[MaxNumEntries];
167 UV =
new double[MaxNumEntries];
170 ierr =
D_->ExtractView(&DV);
185 for (j=0; j< NumIn; j++) {
193 else if (k < 0)
return(-1);
203 if (DiagFound) NumNonzeroDiags++;
204 if (NumU)
U_->ReplaceMyValues(i, NumU, UV, UI);
223 int TotalNonzeroDiags = 0;
246 double MaxDiagonalValue = 1.0/MinDiagonalValue;
255 int MaxNumEntries =
U_->MaxNumEntries() + 1;
257 int * InI =
new int[MaxNumEntries];
258 double * InV =
new double[MaxNumEntries];
262 ierr =
D_->ExtractView(&DV);
264#ifdef IFPACK_FLOPCOUNTERS
265 int current_madds = 0;
274 for (j=0; j<
NumMyCols(); j++) colflag[j] = - 1;
280 NumIn = MaxNumEntries;
281 IFPACK_CHK_ERR(L_->ExtractMyRowCopy(i, NumIn, NumL, InV, InI)==0);
288 IFPACK_CHK_ERR(
U_->ExtractMyRowCopy(i, NumIn-NumL-1, NumU, InV+NumL+1, InI+NumL+1));
294 for (j=0; j<NumIn; j++) colflag[InI[j]] = j;
296 double diagmod = 0.0;
298 for (
int jj=0; jj<NumL; jj++) {
300 double multiplier = InV[jj];
307 for (k=0; k<NumUU; k++) {
308 int kk = colflag[UUI[k]];
310 InV[kk] -= multiplier*UUV[k];
311#ifdef IFPACK_FLOPCOUNTERS
319 for (k=0; k<NumUU; k++) {
320 int kk = colflag[UUI[k]];
321 if (kk>-1) InV[kk] -= multiplier*UUV[k];
322 else diagmod -= multiplier*UUV[k];
323#ifdef IFPACK_FLOPCOUNTERS
339 if (fabs(DV[i]) > MaxDiagonalValue) {
340 if (DV[i] < 0) DV[i] = - MinDiagonalValue;
341 else DV[i] = MinDiagonalValue;
346 for (j=0; j<NumU; j++) UV[j] *= DV[i];
353 for (j=0; j<NumIn; j++) colflag[InI[j]] = -1;
357#ifdef IFPACK_FLOPCOUNTERS
360 double current_flops = 2 * current_madds;
361 double total_flops = 0;
363 Graph_.L_Graph().RowMap().Comm().SumAll(¤t_flops, &total_flops, 1);
366 total_flops += (double) L_->NumGlobalNonzeros();
367 total_flops += (double)
D_->GlobalLength();
368 if (
RelaxValue_!=0.0) total_flops += 2 * (double)
D_->GlobalLength();
392 bool UnitDiagonal =
true;
397 if (
Graph_.LevelOverlap()>0 &&
Graph_.L_Graph().DomainMap().DistributedGlobal()) {
407#ifdef IFPACK_FLOPCOUNTERS
410 L_->SetFlopCounter(*counter);
412 U_->SetFlopCounter(*counter);
418 L_->Solve(Lower, Trans, UnitDiagonal, *X1, *Y1);
420 U_->Solve(Upper, Trans, UnitDiagonal, *Y1, *Y1);
424 U_->Solve(Upper, Trans, UnitDiagonal, *X1, *Y1);
426 L_->Solve(Lower, Trans, UnitDiagonal, *Y1, *Y1);
431 if (
Graph_.LevelOverlap()>0 &&
Graph_.L_Graph().DomainMap().DistributedGlobal())
448 bool UnitDiagonal =
true;
453 if (
Graph_.LevelOverlap()>0 &&
Graph_.L_Graph().DomainMap().DistributedGlobal()) {
470#ifdef IFPACK_FLOPCOUNTERS
473 L_->SetFlopCounter(*counter);
475 U_->SetFlopCounter(*counter);
481 L_->Solve(Lower, Trans, UnitDiagonal, *X1, *Y1);
483 U_->Solve(Upper, Trans, UnitDiagonal, *Y1, *Y1);
487 U_->Solve(Upper, Trans, UnitDiagonal, *X1, *Y1);
489 L_->Solve(Lower, Trans, UnitDiagonal, *Y1, *Y1);
494 if (
Graph_.LevelOverlap()>0 &&
Graph_.L_Graph().DomainMap().DistributedGlobal())
509 bool UnitDiagonal =
true;
514 if (
Graph_.LevelOverlap()>0 &&
Graph_.L_Graph().DomainMap().DistributedGlobal()) {
531#ifdef IFPACK_FLOPCOUNTERS
534 L_->SetFlopCounter(*counter);
536 U_->SetFlopCounter(*counter);
542 L_->Multiply(Trans, *X1, *Y1);
543 Y1->
Update(1.0, *X1, 1.0);
546 U_->Multiply(Trans, Y1temp, *Y1);
547 Y1->
Update(1.0, Y1temp, 1.0);
550 U_->Multiply(Trans, *X1, *Y1);
551 Y1->
Update(1.0, *X1, 1.0);
554 L_->Multiply(Trans, Y1temp, *Y1);
555 Y1->
Update(1.0, Y1temp, 1.0);
559 if (
Graph_.LevelOverlap()>0 &&
Graph_.L_Graph().DomainMap().DistributedGlobal())
597 os <<
" Level of Fill = "; os << LevelFill;
600 os <<
" Level of Overlap = "; os << LevelOverlap;
604 os <<
" Lower Triangle = ";
610 os <<
" Inverse of Diagonal = ";
616 os <<
" Upper Triangle = ";
#define EPETRA_CHK_ERR(a)
const double Epetra_MinDouble
#define IFPACK_CHK_ERR(ifpack_err)
std::ostream & operator<<(std::ostream &os, const Ifpack_CrsRick &A)
<< operator will work for Ifpack_CrsRick.
bool DistributedGlobal() const
const Epetra_Comm & Comm() const
virtual int SumAll(double *PartialSums, double *GlobalSums, int Count) const=0
Epetra_Flops * GetFlopCounter() const
void SetFlopCounter(const Epetra_Flops &FlopCounter_in)
void UpdateFlops(int Flops_in) const
const Epetra_BlockMap & DomainMap() const
const Epetra_BlockMap & RowMap() const
int FillComplete(bool OptimizeDataStorage=true)
int ExtractMyRowCopy(int MyRow, int Length, int &NumEntries, double *Values, int *Indices) const
int MaxNumEntries() const
int Import(const Epetra_SrcDistObject &A, const Epetra_Import &Importer, Epetra_CombineMode CombineMode, const Epetra_OffsetIndex *Indexor=0)
int Export(const Epetra_SrcDistObject &A, const Epetra_Import &Importer, Epetra_CombineMode CombineMode, const Epetra_OffsetIndex *Indexor=0)
int Abs(const Epetra_MultiVector &A)
int Multiply(char TransA, char TransB, double ScalarAB, const Epetra_MultiVector &A, const Epetra_MultiVector &B, double ScalarThis)
int Update(double ScalarA, const Epetra_MultiVector &A, double ScalarThis)
int MaxValue(double *Result) const
int PutScalar(double ScalarConstant)
int ReciprocalMultiply(double ScalarAB, const Epetra_MultiVector &A, const Epetra_MultiVector &B, double ScalarThis)
int NumMyCols() const
Returns the number of local matrix columns.
const Epetra_CrsMatrix & A_
const Ifpack_IlukGraph & Graph_
Ifpack_CrsRick(const Epetra_CrsMatrix &A, const Ifpack_IlukGraph &Graph)
Ifpack_CrsRick constuctor with variable number of indices per row.
virtual int NumGlobalDiagonals() const
Returns the number of diagonal entries found in the global input graph.
Epetra_MultiVector * OverlapX_
Epetra_CombineMode OverlapMode_
virtual int NumMyDiagonals() const
Returns the number of diagonal entries found in the local input graph.
int SetParameters(const Teuchos::ParameterList ¶meterlist, bool cerr_warning_if_unused=false)
Set parameters using a Teuchos::ParameterList object.
int Factor()
Compute IC factor L using the specified graph, diagonal perturbation thresholds and relaxation parame...
const Epetra_Vector & D() const
Returns the address of the D factor associated with this factored matrix.
Epetra_MultiVector * OverlapY_
int NumMyRows() const
Returns the number of local matrix rows.
int InitValues()
Initialize L and U with values from user matrix A.
int NumGlobalRows() const
Returns the number of global matrix rows.
const Ifpack_IlukGraph & Graph() const
Returns the address of the Ifpack_IlukGraph associated with this factored matrix.
virtual ~Ifpack_CrsRick()
Ifpack_CrsRick Destructor.
bool Factored() const
If factor is completed, this query returns true, otherwise it returns false.
void SetFactored(bool Flag)
const Epetra_CrsMatrix & U() const
Returns the address of the U factor associated with this factored matrix.
int Solve(bool Trans, const Epetra_Vector &x, Epetra_Vector &y) const
Returns the result of a Ifpack_CrsRick forward/back solve on a Epetra_Vector x in y.
int SetAllocated(bool Flag)
int Condest(bool Trans, double &ConditionNumberEstimate) const
Returns the maximum over all the condition number estimate for each local IC set of factors.
int Multiply(bool Trans, const Epetra_MultiVector &X, Epetra_MultiVector &Y) const
Returns the result of multiplying U, D and U^T in that order on an Epetra_MultiVector X in Y.
void SetValuesInitialized(bool Flag)
bool ValuesInitialized() const
If values have been initialized, this query returns true, otherwise it returns false.
Ifpack_IlukGraph: A class for constructing level filled graphs for use with ILU(k) class precondition...
virtual int LevelFill() const
Returns the level of fill used to construct this graph.
virtual Epetra_CrsGraph & U_Graph()
Returns the graph of upper triangle of the ILU(k) graph as a Epetra_CrsGraph.
virtual int LevelOverlap() const
Returns the level of overlap used to construct this graph.
virtual Epetra_CrsGraph & L_Graph()
Returns the graph of lower triangle of the ILU(k) graph as a Epetra_CrsGraph.
void set_parameters(const Teuchos::ParameterList ¶meterlist, param_struct ¶ms, bool cerr_warning_if_unused)
Epetra_CombineMode overlap_mode
double double_params[FIRST_INT_PARAM]