Skip to content

Commit

Permalink
Tpetra: removed need for matrix to have isFillActive()=true in scale(…
Browse files Browse the repository at this point in the history
…) and setAllToScalar() (#9660)

* tpetra:  removed need for matrix to have isFillActive() in scale() and
setAllToScalar(), removing need for resumeFill/fillComplete around
calls to these functions #9640

tpetra:  removed caching of Frobenius norm value; there was no way to
guarantee that the cached value would be correct. #9658

muelu, panzer, teko, xpetra, zoltan2:  removed resumeFill/fillComplete
that are no longer needed in light of these changes.

* tpetra:  #9658 removed comments about caching the Frobenius norm
  • Loading branch information
kddevin authored Sep 27, 2021
1 parent 52be2c4 commit 9813022
Show file tree
Hide file tree
Showing 17 changed files with 37 additions and 197 deletions.
Original file line number Diff line number Diff line change
Expand Up @@ -106,9 +106,7 @@ namespace MueLu {
// R_0 = -A*X_0
R = Xpetra::MatrixFactory2<Scalar, LocalOrdinal, GlobalOrdinal, Node>::BuildCopy(T);

R->resumeFill();
R->scale(-one);
R->fillComplete(R->getDomainMap(), R->getRangeMap());

// Z_0 = M^{-1}R_0
Z = Xpetra::MatrixFactory2<Scalar, LocalOrdinal, GlobalOrdinal, Node>::BuildCopy(R);
Expand Down
Original file line number Diff line number Diff line change
Expand Up @@ -122,9 +122,7 @@ namespace MueLu {

rho = sqrt(Utilities::Frobenius(*V[0], *V[0]));

V[0]->resumeFill();
V[0]->scale(-one/rho);
V[0]->fillComplete(V[0]->getDomainMap(), V[0]->getRangeMap());
}

std::vector<SC> h((nIts_+1) * (nIts_+1));
Expand Down Expand Up @@ -182,9 +180,7 @@ namespace MueLu {
// Check for nonsymmetric case
if (h[I(i+1,i)] != zero) {
// Normalize V_i
V[i+1]->resumeFill();
V[i+1]->scale(one/h[I(i+1,i)]);
V[i+1]->fillComplete(V[i+1]->getDomainMap(), V[i+1]->getRangeMap());
}

if (i > 0)
Expand Down
Original file line number Diff line number Diff line change
Expand Up @@ -90,9 +90,7 @@ class TpetraLinearObjContainer : public LinearObjContainer
if(get_f()!=Teuchos::null) get_f()->putScalar(0.0);
if(get_A()!=Teuchos::null) {
Teuchos::RCP<CrsMatrixType> mat = get_A();
mat->resumeFill();
mat->setAllToScalar(0.0);
mat->fillComplete();
}
}

Expand All @@ -119,9 +117,7 @@ class TpetraLinearObjContainer : public LinearObjContainer

void initializeMatrix(ScalarT value)
{
A->resumeFill();
A->setAllToScalar(value);
A->fillComplete();
}

virtual void set_x_th(const Teuchos::RCP<Thyra::VectorBase<ScalarT> > & in)
Expand Down
Original file line number Diff line number Diff line change
Expand Up @@ -127,9 +127,7 @@ void fillComplete(CrsMatrixType & A)

void putScalar(ScalarT s,CrsMatrixType & A)
{
resumeFill(A);
A.setAllToScalar(s);
fillComplete(A);
}

template <typename Intrepid2Type>
Expand Down Expand Up @@ -470,9 +468,7 @@ TEUCHOS_UNIT_TEST(tBlockedTpetraLinearObjFactory, adjustDirichlet)
for(int i=0;i<numBlocks;i++) {
for(int j=0;j<numBlocks;j++) {
RCP<CrsMatrixType> M = getSubBlock_tp(i,j,*b_sys->get_A());
M->resumeFill();
M->setAllToScalar(-3.0);
M->fillComplete(M->getDomainMap(),M->getRangeMap());
}
}

Expand Down
Original file line number Diff line number Diff line change
Expand Up @@ -414,9 +414,7 @@ TEUCHOS_UNIT_TEST(tTpetraLinearObjFactory, adjustDirichlet)
TEST_ASSERT(!Teuchos::is_null(t_sys->get_A()));

t_sys->get_f()->putScalar(-3.0); // put some garbage in the systems
t_sys->get_A()->resumeFill();
t_sys->get_A()->setAllToScalar(-3.0);
t_sys->get_A()->fillComplete();

// there are 3 cases for adjustDirichlet
// 1. Local set only for GID
Expand Down
Original file line number Diff line number Diff line change
Expand Up @@ -135,9 +135,7 @@ Teko::LinearOp FullMaxwellPreconditionerFactory::buildPreconditionerOperator(Tek
RCP<const Thyra::TpetraLinearOp<Scalar,LocalOrdinal,GlobalOrdinal,Node> > tOp = rcp_dynamic_cast<const Thyra::TpetraLinearOp<Scalar,LocalOrdinal,GlobalOrdinal,Node> >(K2,true);
RCP<Thyra::TpetraLinearOp<Scalar,LocalOrdinal,GlobalOrdinal,Node> > tOp2 = Teuchos::rcp_const_cast<Thyra::TpetraLinearOp<Scalar,LocalOrdinal,GlobalOrdinal,Node>>(tOp);
RCP<Tpetra::CrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node> > crsOp = rcp_dynamic_cast<Tpetra::CrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node> >(tOp2->getTpetraOperator(),true);
crsOp->resumeFill();
crsOp->scale(dt);
crsOp->fillComplete(crsOp->getDomainMap(),crsOp->getRangeMap());
diffK = Teko::explicitAdd(K, Teko::scale(-1.0,K2));

writeOut("DiscreteCurl.mm",*C);
Expand Down
Original file line number Diff line number Diff line change
Expand Up @@ -130,9 +130,7 @@ int main(int argc,char * argv[])
RCP<TP_Crs> crsMat = Tpetra::MatrixMarket::Reader<TP_Crs>::readSparseFile("../data/nsjac.mm", Tpetra::getDefaultComm());

RCP<TP_Crs> zeroCrsMat = rcp(new TP_Crs(*crsMat, Teuchos::Copy));
zeroCrsMat->resumeFill();
zeroCrsMat->setAllToScalar(0.0);
zeroCrsMat->fillComplete();

RCP<TP_Op> Mat = crsMat;
RCP<TP_Op> zeroMat = zeroCrsMat;
Expand Down
2 changes: 0 additions & 2 deletions packages/teko/tests/src/Tpetra/tBlockedTpetraOperator.cpp
Original file line number Diff line number Diff line change
Expand Up @@ -236,9 +236,7 @@ bool tBlockedTpetraOperator::test_vector_constr(int verbosity,std::ostream & os)
<< "testing tolerance over many matrix vector multiplies ( " << max << " <= "
<< tolerance_ << " )");

A->resumeFill();
A->scale(2.0); //double everything
A->fillComplete();

ST afterNorm = A->getFrobeniusNorm();
TEST_ASSERT(beforeNorm!=afterNorm,
Expand Down
4 changes: 0 additions & 4 deletions packages/teko/tests/src/Tpetra/tStridedTpetraOperator.cpp
Original file line number Diff line number Diff line change
Expand Up @@ -209,9 +209,7 @@ bool tStridedTpetraOperator::test_numvars_constr(int verbosity,std::ostream & os
<< "testing tolerance over many matrix vector multiplies ( " << max << " <= "
<< tolerance_ << " )");

A->resumeFill();
A->scale(2.0); //double everything
A->fillComplete();

ST afterNorm = A->getFrobeniusNorm();
TEST_ASSERT(beforeNorm!=afterNorm,
Expand Down Expand Up @@ -320,9 +318,7 @@ bool tStridedTpetraOperator::test_vector_constr(int verbosity,std::ostream & os)
<< "testing tolerance over many matrix vector multiplies ( " << max << " <= "
<< tolerance_ << " )");

A->resumeFill();
A->scale(2.0); // double everything
A->fillComplete();

ST afterNorm = A->getFrobeniusNorm();
TEST_ASSERT(beforeNorm!=afterNorm,
Expand Down
8 changes: 0 additions & 8 deletions packages/teko/tests/src/tExplicitOps_tpetra.cpp
Original file line number Diff line number Diff line change
Expand Up @@ -274,14 +274,10 @@ bool tExplicitOps_tpetra::test_mult_modScaleMatProd(int verbosity,std::ostream &

RCP<Thyra::TpetraLinearOp<ST,LO,GO,NT> > tF = Teuchos::rcp_dynamic_cast<Thyra::TpetraLinearOp<ST,LO,GO,NT> >(F_);
RCP<Tpetra::CrsMatrix<ST,LO,GO,NT> > crsF = Teuchos::rcp_dynamic_cast<Tpetra::CrsMatrix<ST,LO,GO,NT> >(tF->getTpetraOperator(),true);
crsF->resumeFill();
crsF->scale(5.0);
crsF->fillComplete();
RCP<Thyra::TpetraLinearOp<ST,LO,GO,NT> > tG = Teuchos::rcp_dynamic_cast<Thyra::TpetraLinearOp<ST,LO,GO,NT> >(G_);
RCP<Tpetra::CrsMatrix<ST,LO,GO,NT> > crsG = Teuchos::rcp_dynamic_cast<Tpetra::CrsMatrix<ST,LO,GO,NT> >(tG->getTpetraOperator(),true);
crsG->resumeFill();
crsG->scale(2.0);
crsG->fillComplete();

// do some random violence (oh my brothers) to one row
size_t numEntries = crsF->getNumEntriesInLocalRow (3);
Expand Down Expand Up @@ -380,14 +376,10 @@ bool tExplicitOps_tpetra::test_add_mod(int verbosity,std::ostream & os)

RCP<Thyra::TpetraLinearOp<ST,LO,GO,NT> > tF = Teuchos::rcp_dynamic_cast<Thyra::TpetraLinearOp<ST,LO,GO,NT> >(F_);
RCP<Tpetra::CrsMatrix<ST,LO,GO,NT> > crsF = Teuchos::rcp_dynamic_cast<Tpetra::CrsMatrix<ST,LO,GO,NT> >(tF->getTpetraOperator(),true);
crsF->resumeFill();
crsF->scale(5.0);
crsF->fillComplete();
RCP<Thyra::TpetraLinearOp<ST,LO,GO,NT> > tG = Teuchos::rcp_dynamic_cast<Thyra::TpetraLinearOp<ST,LO,GO,NT> >(G_);
RCP<Tpetra::CrsMatrix<ST,LO,GO,NT> > crsG = Teuchos::rcp_dynamic_cast<Tpetra::CrsMatrix<ST,LO,GO,NT> >(tG->getTpetraOperator(),true);
crsG->resumeFill();
crsG->scale(2.0);
crsG->fillComplete();

// do some random violence (oh my brothers) to one row
size_t numEntries = crsF->getNumEntriesInLocalRow (3);
Expand Down
4 changes: 0 additions & 4 deletions packages/tpetra/core/src/Tpetra_CrsGraph_def.hpp
Original file line number Diff line number Diff line change
Expand Up @@ -2920,10 +2920,6 @@ namespace Tpetra {
allocateIndices (LocalIndices, verbose_);
}

// FIXME (mfh 13 Aug 2014) What if they haven't been cleared on
// all processes?
clearGlobalConstants ();

if (k_numRowEntries_.extent (0) != 0) {
this->k_numRowEntries_(lrow) = 0;
}
Expand Down
34 changes: 5 additions & 29 deletions packages/tpetra/core/src/Tpetra_CrsMatrix_decl.hpp
Original file line number Diff line number Diff line change
Expand Up @@ -2484,10 +2484,6 @@ namespace Tpetra {
/// \|A\|_F = \sqrt{\sum_{i,j} \|A(i,j)\|^2}.
/// \f\].
///
/// If the matrix is fill complete, then the computed value is
/// cached; the cache is cleared whenever resumeFill() is called.
/// Otherwise, the value is computed every time the method is
/// called.
mag_type getFrobeniusNorm () const override;

/// \brief Return \c true if getLocalRowView() and
Expand Down Expand Up @@ -3865,27 +3861,14 @@ namespace Tpetra {
sortAndMergeIndicesAndValues (const bool sorted,
const bool merged);

/// \brief Clear matrix properties that require collectives.
///
/// This clears whatever computeGlobalConstants() (which see)
/// computed, in preparation for changes to the matrix. The
/// current implementation of this method does nothing.
///
/// This method is called in resumeFill().
void clearGlobalConstants();

/// \brief Compute matrix properties that require collectives.
///
/// The corresponding Epetra_CrsGraph method computes things
/// like the global number of nonzero entries, that require
/// collectives over the matrix's communicator. The current
/// Tpetra implementation of this method does nothing.
///
public:
/// This method is called in fillComplete().
void computeGlobalConstants();
#ifdef TPETRA_ENABLE_DEPRECATED_CODE
TPETRA_DEPRECATED void computeGlobalConstants() {};
#endif // TPETRA_ENABLE_DEPRECATED_CODE

//! Returns true if globalConstants have been computed; false otherwise
bool haveGlobalConstants() const;

protected:
/// \brief Column Map MultiVector used in apply().
///
Expand Down Expand Up @@ -4116,13 +4099,6 @@ namespace Tpetra {
std::map<GlobalOrdinal, std::pair<Teuchos::Array<GlobalOrdinal>,
Teuchos::Array<Scalar> > > nonlocals_;

/// \brief Cached Frobenius norm of the (global) matrix.
///
/// The value -1 means that the norm has not yet been computed, or
/// that the values in the matrix may have changed and the norm
/// must be recomputed.
mutable mag_type frobNorm_ = -STM::one();

public:
// FIXME (mfh 24 Feb 2014) Is it _really_ necessary to make this a
// public inner class of CrsMatrix? It looks like it doesn't
Expand Down
Loading

0 comments on commit 9813022

Please sign in to comment.