28 #include <oomph-lib-config.h> 
   50     if (dist_matrix_pt != 0)
 
   67             "Zero diagonal in matrix --> Cannot use diagonal preconditioner.",
 
   68             OOMPH_CURRENT_FUNCTION,
 
   69             OOMPH_EXCEPTION_LOCATION);
 
   90       for (
unsigned i = 0; 
i < n_row; 
i++)
 
   96             "Zero diagonal in matrix --> Cannot use diagonal preconditioner.",
 
   97             OOMPH_CURRENT_FUNCTION,
 
   98             OOMPH_EXCEPTION_LOCATION);
 
  123       std::ostringstream error_message_stream;
 
  125         << 
"The r vector must have the same distribution as the " 
  127         << 
"(this is the same as the matrix passed to setup())";
 
  129                           OOMPH_CURRENT_FUNCTION,
 
  130                           OOMPH_EXCEPTION_LOCATION);
 
  136         std::ostringstream error_message_stream;
 
  138           << 
"The z vector distribution has been setup; it must have the " 
  139           << 
"same distribution as the r vector (and preconditioner).";
 
  141                             OOMPH_CURRENT_FUNCTION,
 
  142                             OOMPH_EXCEPTION_LOCATION);
 
  170     Nrow = matrix_pt()->nrow();
 
  173     if (Inv_lumped_diag_pt != 0)
 
  175       delete[] this->Inv_lumped_diag_pt;
 
  177     Inv_lumped_diag_pt = 
new double[this->Nrow];
 
  180     for (
unsigned i = 0; 
i < Nrow; 
i++)
 
  182       Inv_lumped_diag_pt[
i] = 0.0;
 
  189     int* m_row_index = cc_matrix_pt->
row_index();
 
  190     double* m_value = cc_matrix_pt->
value();
 
  191     unsigned m_nnz = cc_matrix_pt->
nnz();
 
  194     Positive_matrix = 
true;
 
  197     for (
unsigned i = 0; 
i < m_nnz; 
i++)
 
  200       if (m_value[
i] < 0.0)
 
  202         Positive_matrix = 
false;
 
  206       Inv_lumped_diag_pt[m_row_index[
i]] += m_value[
i];
 
  210     for (
unsigned i = 0; 
i < Nrow; 
i++)
 
  212       Inv_lumped_diag_pt[
i] = 1.0 / Inv_lumped_diag_pt[
i];
 
  217     this->build_distribution(dist);
 
  233     if (Inv_lumped_diag_pt != 0)
 
  235       delete[] this->Inv_lumped_diag_pt;
 
  237     Inv_lumped_diag_pt = 
new double[this->Nrow];
 
  240     for (
unsigned i = 0; 
i < Nrow; 
i++)
 
  242       Inv_lumped_diag_pt[
i] = 0.0;
 
  246     int* m_row_start = cr_matrix_pt->
row_start();
 
  247     double* m_value = cr_matrix_pt->
value();
 
  250     Positive_matrix = 
true;
 
  253     for (
unsigned i = 0; 
i < Nrow; 
i++)
 
  255       Inv_lumped_diag_pt[
i] = 0.0;
 
  256       for (
int j = m_row_start[
i]; j < m_row_start[
i + 1]; j++)
 
  259         if (m_value[j] < 0.0)
 
  261           Positive_matrix = 
false;
 
  265         Inv_lumped_diag_pt[
i] += m_value[j];
 
  268       Inv_lumped_diag_pt[
i] = 1.0 / Inv_lumped_diag_pt[
i];
 
  279   template<
typename MATRIX>
 
  286       std::ostringstream error_message_stream;
 
  288         << 
"The r vector must have teh same distribution as the " 
  290         << 
"(this is the same as the matrix passed to setup())";
 
  292                           OOMPH_CURRENT_FUNCTION,
 
  293                           OOMPH_EXCEPTION_LOCATION);
 
  299         std::ostringstream error_message_stream;
 
  301           << 
"The z vector distribution has been setup; it must have the " 
  302           << 
"same distribution as the r vector (and preconditioner).";
 
  304                             OOMPH_CURRENT_FUNCTION,
 
  305                             OOMPH_EXCEPTION_LOCATION);
 
  311     for (
unsigned i = 0; 
i < Nrow; 
i++)
 
  313       z[
i] = Inv_lumped_diag_pt[
i] * r[
i];
 
  332     if (cc_matrix_pt == 0)
 
  334       std::ostringstream error_msg;
 
  335       error_msg << 
"Failed to conver matrix_pt to CCDoubleMatrix*.";
 
  337         error_msg.str(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION);
 
  342     int n_row = cc_matrix_pt->
nrow();
 
  360     m_value = cc_matrix_pt->
value();
 
  363     for (
int i = 0; 
i < n_row; 
i++)
 
  365       for (
int j = m_column_start[
i]; j < m_column_start[
i + 1]; j++)
 
  367         if (m_row_index[j] > 
i)
 
  380     L_column_start.resize(n_row + 1);
 
  381     L_row_entry.resize(l_nz);
 
  384     U_column_start.resize(n_row + 1);
 
  385     U_row_entry.resize(u_nz);
 
  388     L_column_start[0] = 0;
 
  389     U_column_start[0] = 0;
 
  392     for (
int i = 0; 
i < n_row; 
i++)
 
  394       L_column_start[
i + 1] = L_column_start[
i];
 
  395       U_column_start[
i + 1] = U_column_start[
i];
 
  396       for (
int j = m_column_start[
i]; j < m_column_start[
i + 1]; j++)
 
  398         if (m_row_index[j] > 
i)
 
  400           int k = L_column_start[
i + 1]++;
 
  401           L_row_entry[k].index() = m_row_index[j];
 
  402           L_row_entry[k].value() = m_value[j];
 
  406           int k = U_column_start[
i + 1]++;
 
  407           U_row_entry[k].index() = m_row_index[j];
 
  408           U_row_entry[k].value() = m_value[j];
 
  416     for (
unsigned i = 0; 
i < unsigned(n_row); 
i++)
 
  419       std::sort(L_row_entry.begin() + L_column_start[
i],
 
  420                 L_row_entry.begin() + L_column_start[
i + 1]);
 
  423       std::sort(U_row_entry.begin() + U_column_start[
i],
 
  424                 U_row_entry.begin() + U_column_start[
i + 1]);
 
  430     unsigned j, pn, qn, rn;
 
  435     for (
i = 0; 
i < n_row - 1; 
i++)
 
  437       multiplier = U_row_entry[U_column_start[
i + 1] - 1].value();
 
  438       for (j = L_column_start[
i]; j < L_column_start[
i + 1]; j++)
 
  439         L_row_entry[j].value() /= multiplier;
 
  440       for (j = U_column_start[
i + 1]; j < U_column_start[
i + 2] - 1; j++)
 
  442         multiplier = U_row_entry[j].value();
 
  444         rn = L_column_start[
i + 1];
 
  445         for (pn = L_column_start[U_row_entry[j].index()];
 
  446              pn < L_column_start[U_row_entry[j].index() + 1] &&
 
  447              static_cast<int>(L_row_entry[pn].index()) <= 
i + 1;
 
  450           while (qn < U_column_start[
i + 2] &&
 
  451                  U_row_entry[qn].index() < L_row_entry[pn].index())
 
  453           if (qn < U_column_start[
i + 2] &&
 
  454               L_row_entry[pn].index() == U_row_entry[qn].index())
 
  455             U_row_entry[qn].value() -= multiplier * L_row_entry[pn].value();
 
  457         for (; pn < L_column_start[U_row_entry[j].index() + 1]; pn++)
 
  459           while (rn < L_column_start[
i + 2] &&
 
  460                  L_row_entry[rn].index() < L_row_entry[pn].index())
 
  462           if (rn < L_column_start[
i + 2] &&
 
  463               L_row_entry[pn].index() == L_row_entry[rn].index())
 
  464             L_row_entry[rn].value() -= multiplier * L_row_entry[pn].value();
 
  480     if (cr_matrix_pt == 0)
 
  482       std::ostringstream error_msg;
 
  483       error_msg << 
"Failed to conver matrix_pt to CRDoubleMatrix*.";
 
  485         error_msg.str(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION);
 
  490     bool built_global = 
false;
 
  497       cr_matrix_pt = global_matrix_pt;
 
  507     int n_row = cr_matrix_pt->
nrow();
 
  521     m_value = cr_matrix_pt->
value();
 
  524     for (
int i = 0; 
i < n_row; 
i++)
 
  526       for (
int j = m_row_start[
i]; j < m_row_start[
i + 1]; j++)
 
  528         if (m_column_index[j] < 
i)
 
  541     L_row_start.resize(n_row + 1);
 
  542     L_row_entry.resize(l_nz);
 
  545     U_row_start.resize(n_row + 1);
 
  546     U_row_entry.resize(u_nz);
 
  553     for (
int i = 0; 
i < n_row; 
i++)
 
  555       L_row_start[
i + 1] = L_row_start[
i];
 
  556       U_row_start[
i + 1] = U_row_start[
i];
 
  557       for (
int j = m_row_start[
i]; j < m_row_start[
i + 1]; j++)
 
  559         if (m_column_index[j] < 
i)
 
  561           int k = L_row_start[
i + 1]++;
 
  562           L_row_entry[k].value() = m_value[j];
 
  563           L_row_entry[k].index() = m_column_index[j];
 
  567           int k = U_row_start[
i + 1]++;
 
  568           U_row_entry[k].value() = m_value[j];
 
  569           U_row_entry[k].index() = m_column_index[j];
 
  576     unsigned i, j, pn, qn, rn;
 
  581     for (
i = 1; i < static_cast<unsigned>(n_row); 
i++)
 
  583       for (j = L_row_start[
i]; j < L_row_start[
i + 1]; j++)
 
  585         pn = U_row_start[L_row_entry[j].index()];
 
  586         multiplier = (L_row_entry[j].value() /= U_row_entry[pn].value());
 
  589         for (pn++; pn < U_row_start[L_row_entry[j].index() + 1] &&
 
  590                    U_row_entry[pn].index() < 
i;
 
  593           while (qn < L_row_start[
i + 1] &&
 
  594                  L_row_entry[qn].index() < U_row_entry[pn].index())
 
  596           if (qn < L_row_start[
i + 1] &&
 
  597               U_row_entry[pn].index() == L_row_entry[qn].index())
 
  598             L_row_entry[qn].value() -= multiplier * U_row_entry[pn].value();
 
  600         for (; pn < U_row_start[L_row_entry[j].index() + 1]; pn++)
 
  602           while (rn < U_row_start[
i + 1] &&
 
  603                  U_row_entry[rn].index() < U_row_entry[pn].index())
 
  605           if (rn < U_row_start[
i + 1] &&
 
  606               U_row_entry[pn].index() == U_row_entry[rn].index())
 
  607             U_row_entry[rn].value() -= multiplier * U_row_entry[pn].value();
 
  628     int n_row = r.
nrow();
 
  647     for (
unsigned i = 0; i < static_cast<unsigned>(n_row); 
i++)
 
  649       for (
unsigned j = L_column_start[
i]; j < L_column_start[
i + 1]; j++)
 
  651         z[L_row_entry[j].index()] =
 
  652           z[L_row_entry[j].index()] - z[
i] * L_row_entry[j].value();
 
  658     for (
int i = n_row - 1; 
i >= 0; 
i--)
 
  660       x = z[
i] / U_row_entry[U_column_start[
i + 1] - 1].value();
 
  662       for (
unsigned j = U_column_start[
i]; j < U_column_start[
i + 1] - 1; j++)
 
  664         z[U_row_entry[j].index()] =
 
  665           z[U_row_entry[j].index()] - x * U_row_entry[j].value();
 
  686     int n_row = r.
nrow();
 
  706     for (
int i = 0; 
i < n_row; 
i++)
 
  709       for (
unsigned j = L_row_start[
i]; j < L_row_start[
i + 1]; j++)
 
  711         t = 
t + L_row_entry[j].value() * z[L_row_entry[j].index()];
 
  717     for (
int i = n_row - 1; 
i >= 0; 
i--)
 
  720       for (
unsigned j = U_row_start[
i] + 1; j < U_row_start[
i + 1]; j++)
 
  722         t = 
t + U_row_entry[j].value() * z[U_row_entry[j].index()];
 
  725       z[
i] = z[
i] / U_row_entry[U_row_start[
i]].value();
 
//////////////////////////////////////////////////////////////// ////////////////////////////////////...
unsigned long nrow() const
Return the number of rows of the matrix.
int * column_start()
Access to C-style column_start array.
int * row_index()
Access to C-style row index array.
A class for compressed row matrices. This is a distributable object.
int * column_index()
Access to C-style column index array.
int * row_start()
Access to C-style row_start array.
CRDoubleMatrix * global_matrix() const
if this matrix is distributed then a the equivalent global matrix is built using new and returned....
double * value()
Access to C-style value array.
unsigned long nrow() const
Return the number of rows of the matrix.
Base class for any linear algebra object that is distributable. Just contains storage for the LinearA...
bool distributed() const
distribution is serial or distributed
LinearAlgebraDistribution * distribution_pt() const
access to the LinearAlgebraDistribution
unsigned nrow() const
access function to the number of global rows.
unsigned nrow_local() const
access function for the num of local rows on this processor.
unsigned first_row() const
access function for the first row on this processor
void build_distribution(const LinearAlgebraDistribution *const dist_pt)
setup the distribution of this distributable linear algebra object
virtual unsigned long nrow() const =0
Return the number of rows of the matrix.
A vector in the mathematical sense, initially developed for linear algebra type applications....
void build(const DoubleVector &old_vector)
Just copys the argument DoubleVector.
double * values_pt()
access function to the underlying values
void redistribute(const LinearAlgebraDistribution *const &dist_pt)
The contents of the vector are redistributed to match the new distribution. In a non-MPI rebuild this...
Describes the distribution of a distributable linear algebra type object. Typically this is a contain...
void preconditioner_solve(const DoubleVector &r, DoubleVector &z)
Apply preconditioner to z, i.e. z=D^-1.
Vector< double > Inv_diag
Vector of inverse diagonal entries.
void setup()
Setup the preconditioner (store diagonal) from the fully assembled matrix.
Matrix-based lumped preconditioner.
void preconditioner_solve(const DoubleVector &r, DoubleVector &z)
Apply preconditioner to z, i.e. z=D^-1.
void setup()
Setup the preconditioner (store diagonal) from the fully assembled matrix. Problem pointer is ignored...
An OomphLibError object which should be thrown when an run-time error is encountered....
virtual DoubleMatrixBase * matrix_pt() const
Get function for matrix pointer.
virtual void preconditioner_solve(const DoubleVector &r, DoubleVector &z)=0
Apply the preconditioner. Pure virtual generic interface function. This method should apply the preco...
virtual const OomphCommunicator * comm_pt() const
Get function for comm pointer.
virtual void setup()=0
Setup the preconditioner. Pure virtual generic interface function.
T * value()
Access to C-style value array.
unsigned long nnz() const
Return the number of nonzero entries.
//////////////////////////////////////////////////////////////////// ////////////////////////////////...