|
IFPACK Development
|
00001 /*@HEADER 00002 // *********************************************************************** 00003 // 00004 // Ifpack: Object-Oriented Algebraic Preconditioner Package 00005 // Copyright (2002) Sandia Corporation 00006 // 00007 // Under terms of Contract DE-AC04-94AL85000, there is a non-exclusive 00008 // license for use of this work by or on behalf of the U.S. Government. 00009 // 00010 // This library is free software; you can redistribute it and/or modify 00011 // it under the terms of the GNU Lesser General Public License as 00012 // published by the Free Software Foundation; either version 2.1 of the 00013 // License, or (at your option) any later version. 00014 // 00015 // This library is distributed in the hope that it will be useful, but 00016 // WITHOUT ANY WARRANTY; without even the implied warranty of 00017 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU 00018 // Lesser General Public License for more details. 00019 // 00020 // You should have received a copy of the GNU Lesser General Public 00021 // License along with this library; if not, write to the Free Software 00022 // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 00023 // USA 00024 // Questions? Contact Michael A. Heroux (maherou@sandia.gov) 00025 // 00026 // *********************************************************************** 00027 //@HEADER 00028 */ 00029 00030 #include "Ifpack_ConfigDefs.h" 00031 #include "Ifpack_SparsityFilter.h" 00032 #include "Epetra_ConfigDefs.h" 00033 #include "Epetra_RowMatrix.h" 00034 #include "Epetra_Comm.h" 00035 #include "Epetra_Map.h" 00036 #include "Epetra_MultiVector.h" 00037 #include "Epetra_Vector.h" 00038 00039 //============================================================================== 00040 Ifpack_SparsityFilter::Ifpack_SparsityFilter(const Teuchos::RefCountPtr<Epetra_RowMatrix>& Matrix, 00041 int AllowedEntries, 00042 int AllowedBandwidth) : 00043 A_(Matrix), 00044 MaxNumEntries_(0), 00045 MaxNumEntriesA_(0), 00046 AllowedBandwidth_(AllowedBandwidth), 00047 AllowedEntries_(AllowedEntries), 00048 NumNonzeros_(0), 00049 NumRows_(0) 00050 { 00051 // use this filter only on serial matrices 00052 if (A_->Comm().NumProc() != 1) { 00053 cerr << "Ifpack_SparsityFilter can be used with Comm().NumProc() == 1" << endl; 00054 cerr << "only. This class is a tool for Ifpack_AdditiveSchwarz," << endl; 00055 cerr << "and it is not meant to be used otherwise." << endl; 00056 exit(EXIT_FAILURE); 00057 } 00058 00059 // only square serial matrices 00060 if ((A_->NumMyRows() != A_->NumMyCols()) || 00061 (A_->NumMyRows() != A_->NumGlobalRows())) 00062 IFPACK_CHK_ERRV(-1); 00063 00064 NumRows_ = A_->NumMyRows(); 00065 MaxNumEntriesA_ = A_->MaxNumEntries(); 00066 Indices_.resize(MaxNumEntriesA_); 00067 Values_.resize(MaxNumEntriesA_); 00068 00069 // default value is to not consider bandwidth 00070 if (AllowedBandwidth_ == -1) 00071 AllowedBandwidth_ = NumRows_; 00072 00073 // computes the number of nonzero elements per row in the 00074 // dropped matrix. Stores this number in NumEntries_. 00075 // Also, computes the global number of nonzeros. 00076 vector<int> Ind(MaxNumEntriesA_); 00077 vector<double> Val(MaxNumEntriesA_); 00078 00079 NumEntries_.resize(NumRows_); 00080 for (int i = 0 ; i < NumRows_ ; ++i) 00081 NumEntries_[i] = MaxNumEntriesA_; 00082 00083 for (int i = 0 ; i < A_->NumMyRows() ; ++i) { 00084 int Nnz; 00085 IFPACK_CHK_ERRV(ExtractMyRowCopy(i,MaxNumEntriesA_,Nnz, 00086 &Val[0], &Ind[0])); 00087 00088 NumEntries_[i] = Nnz; 00089 NumNonzeros_ += Nnz; 00090 if (Nnz > MaxNumEntries_) 00091 MaxNumEntries_ = Nnz; 00092 00093 } 00094 } 00095 00096 //============================================================================== 00097 int Ifpack_SparsityFilter:: 00098 ExtractMyRowCopy(int MyRow, int Length, int & NumEntries, 00099 double *Values, int * Indices) const 00100 { 00101 if (Length < NumEntries_[MyRow]) 00102 IFPACK_CHK_ERR(-1); 00103 00104 int Nnz; 00105 IFPACK_CHK_ERR(A_->ExtractMyRowCopy(MyRow,MaxNumEntriesA_,Nnz, 00106 &Values_[0],&Indices_[0])); 00107 00108 double Threshold = 0.0; 00109 00110 // this `if' is to define the cut-off value 00111 if (Nnz > AllowedEntries_) { 00112 00113 vector<double> Values2(Nnz); 00114 int count = 0; 00115 for (int i = 0 ; i < Nnz ; ++i) { 00116 // skip diagonal entry (which is always inserted) 00117 if (Indices_[i] == MyRow) 00118 continue; 00119 // put absolute value 00120 Values2[count] = IFPACK_ABS(Values_[i]); 00121 count++; 00122 } 00123 00124 for (int i = count ; i < Nnz ; ++i) 00125 Values2[i] = 0.0; 00126 00127 // sort in descending order 00128 sort(Values2.rbegin(),Values2.rend()); 00129 // get the cut-off value 00130 Threshold = Values2[AllowedEntries_ - 1]; 00131 00132 } 00133 00134 // loop over all nonzero elements of row MyRow, 00135 // and drop elements below specified threshold. 00136 // Also, add value to zero diagonal 00137 NumEntries = 0; 00138 00139 for (int i = 0 ; i < Nnz ; ++i) { 00140 00141 if (IFPACK_ABS(Indices_[i] - MyRow) > AllowedBandwidth_) 00142 continue; 00143 00144 if ((Indices_[i] != MyRow) && (IFPACK_ABS(Values_[i]) < Threshold)) 00145 continue; 00146 00147 Values[NumEntries] = Values_[i]; 00148 Indices[NumEntries] = Indices_[i]; 00149 00150 NumEntries++; 00151 if (NumEntries > AllowedEntries_) 00152 break; 00153 } 00154 00155 return(0); 00156 } 00157 00158 //============================================================================== 00159 int Ifpack_SparsityFilter:: 00160 ExtractDiagonalCopy(Epetra_Vector & Diagonal) const 00161 { 00162 IFPACK_RETURN(A_->ExtractDiagonalCopy(Diagonal)); 00163 } 00164 00165 //============================================================================== 00166 int Ifpack_SparsityFilter:: 00167 Multiply(bool TransA, const Epetra_MultiVector& X, 00168 Epetra_MultiVector& Y) const 00169 { 00170 00171 int NumVectors = X.NumVectors(); 00172 if (NumVectors != Y.NumVectors()) 00173 IFPACK_CHK_ERR(-1); 00174 00175 Y.PutScalar(0.0); 00176 00177 vector<int> Indices(MaxNumEntries_); 00178 vector<double> Values(MaxNumEntries_); 00179 00180 for (int i = 0 ; i < A_->NumMyRows() ; ++i) { 00181 00182 int Nnz; 00183 ExtractMyRowCopy(i,MaxNumEntries_,Nnz, 00184 &Values[0], &Indices[0]); 00185 if (!TransA) { 00186 // no transpose first 00187 for (int j = 0 ; j < NumVectors ; ++j) { 00188 for (int k = 0 ; k < Nnz ; ++k) { 00189 Y[j][i] += Values[k] * X[j][Indices[k]]; 00190 } 00191 } 00192 } 00193 else { 00194 // transpose here 00195 for (int j = 0 ; j < NumVectors ; ++j) { 00196 for (int k = 0 ; k < Nnz ; ++k) { 00197 Y[j][Indices[k]] += Values[k] * X[j][i]; 00198 } 00199 } 00200 } 00201 } 00202 00203 return(0); 00204 } 00205 00206 //============================================================================== 00207 int Ifpack_SparsityFilter:: 00208 Solve(bool Upper, bool Trans, bool UnitDiagonal, 00209 const Epetra_MultiVector& X, Epetra_MultiVector& Y) const 00210 { 00211 IFPACK_CHK_ERR(-98); 00212 } 00213 00214 //============================================================================== 00215 int Ifpack_SparsityFilter:: 00216 Apply(const Epetra_MultiVector& X, Epetra_MultiVector& Y) const 00217 { 00218 IFPACK_RETURN(Multiply(UseTranspose(),X,Y)); 00219 } 00220 00221 //============================================================================== 00222 int Ifpack_SparsityFilter:: 00223 ApplyInverse(const Epetra_MultiVector& X, Epetra_MultiVector& Y) const 00224 { 00225 IFPACK_CHK_ERR(-98); 00226 }
1.7.4