void ReuseTpetraPreconditioner(const Teuchos::RCP<Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> >& inA,
                                 MueLu::TpetraOperator<Scalar,LocalOrdinal,GlobalOrdinal,Node>& Op) {
    typedef Scalar          SC;
    typedef LocalOrdinal    LO;
    typedef GlobalOrdinal   GO;
    typedef Node            NO;

    typedef Xpetra::Matrix<SC,LO,GO,NO>     Matrix;
    typedef Xpetra::Operator<SC,LO,GO,NO>   Operator;
    typedef MueLu ::Hierarchy<SC,LO,GO,NO>  Hierarchy;

    RCP<Hierarchy> H = Op.GetHierarchy();

    TEUCHOS_TEST_FOR_EXCEPTION(!H->GetNumLevels(), Exceptions::RuntimeError,
                               "ReuseTpetraPreconditioner: Hierarchy has no levels in it");
    TEUCHOS_TEST_FOR_EXCEPTION(!H->GetLevel(0)->IsAvailable("A"), Exceptions::RuntimeError,
                               "ReuseTpetraPreconditioner: Hierarchy has no fine level operator");
    RCP<Level> level0 = H->GetLevel(0);

    RCP<Operator> O0 = level0->Get<RCP<Operator> >("A");
    RCP<Matrix>   A0 = Teuchos::rcp_dynamic_cast<Matrix>(O0);

    RCP<Matrix> A = TpetraCrs_To_XpetraMatrix<SC,LO,GO,NO>(inA);
    if (!A0.is_null()) {
      // If a user provided a "number of equations" argument in a parameter list
      // during the initial setup, we must honor that settings and reuse it for
      // all consequent setups.
      A->SetFixedBlockSize(A0->GetFixedBlockSize());
    }
    level0->Set("A", A);

    H->SetupRe();
  }
示例#2
0
int main(int argc, char *argv[]) {
#include "MueLu_UseShortNames.hpp"

  using Teuchos::RCP;
  using namespace MueLuExamples;

  Teuchos::oblackholestream blackhole;
  Teuchos::GlobalMPISession mpiSession(&argc,&argv,&blackhole);

  RCP<const Teuchos::Comm<int> > comm = Teuchos::DefaultComm<int>::getComm();
  RCP<Teuchos::FancyOStream> out = Teuchos::fancyOStream(Teuchos::rcpFromRef(std::cout));
  out->setOutputToRootOnly(0);
  *out << MueLu::MemUtils::PrintMemoryUsage() << std::endl;

  // Timing
  Teuchos::Time myTime("global");
  Teuchos::TimeMonitor Monitor(myTime);


#ifndef HAVE_TEUCHOS_LONG_LONG_INT
  *out << "Warning: scaling test was not compiled with long long int support" << std::endl;
#endif

  // custom parameters
  LO maxLevels = 4;

  int maxCoarseSize=1; //FIXME clp doesn't like long long int
  std::string aggOrdering = "natural";
  int minPerAgg=3;
  int maxNbrAlreadySelected=0;

  // read in problem
  Epetra_Map emap(10201,0,*Xpetra::toEpetra(comm));
  Epetra_CrsMatrix * ptrA = 0;
  Epetra_Vector * ptrf = 0;

  std::cout << "Reading matrix market file" << std::endl;
  EpetraExt::MatrixMarketFileToCrsMatrix("Condif2Mat.mat",emap,emap,emap,ptrA);
  EpetraExt::MatrixMarketFileToVector("Condif2Rhs.mat",emap,ptrf);
  RCP<Epetra_CrsMatrix> epA = Teuchos::rcp(ptrA);
  RCP<Epetra_Vector> epv = Teuchos::rcp(ptrf);

  // Epetra_CrsMatrix -> Xpetra::Matrix
  RCP<Xpetra::CrsMatrix<double, int, int> > exA = Teuchos::rcp(new Xpetra::EpetraCrsMatrix(epA));
  RCP<Xpetra::CrsMatrixWrap<double, int, int> > crsOp = Teuchos::rcp(new Xpetra::CrsMatrixWrap<double, int, int>(exA));
  RCP<Xpetra::Matrix<double, int, int> > Op = Teuchos::rcp_dynamic_cast<Xpetra::Matrix<double, int, int> >(crsOp);

  // Epetra_Vector -> Xpetra::Vector
  RCP<Xpetra::Vector<double,int,int> > xRhs = Teuchos::rcp(new Xpetra::EpetraVector(epv));

  // Epetra_Map -> Xpetra::Map
  const RCP< const Xpetra::Map<int, int> > map = Xpetra::toXpetra(emap);

  // build nullspace
  RCP<MultiVector> nullSpace = MultiVectorFactory::Build(map,1);
  nullSpace->putScalar( (SC) 1.0);

  RCP<MueLu::Hierarchy<SC,LO,GO,NO,LMO> > H = rcp ( new Hierarchy() );
  H->setDefaultVerbLevel(Teuchos::VERB_HIGH);
  H->SetMaxCoarseSize((GO) maxCoarseSize);;

  // build finest Level
  RCP<MueLu::Level> Finest = H->GetLevel();
  Finest->setDefaultVerbLevel(Teuchos::VERB_HIGH);
  Finest->Set("A",Op);
  Finest->Set("Nullspace",nullSpace);

  RCP<CoupledAggregationFactory> CoupledAggFact = rcp(new CoupledAggregationFactory());
  *out << "========================= Aggregate option summary  =========================" << std::endl;
  *out << "min DOFs per aggregate :                " << minPerAgg << std::endl;
  *out << "min # of root nbrs already aggregated : " << maxNbrAlreadySelected << std::endl;
  CoupledAggFact->SetMinNodesPerAggregate(minPerAgg); //TODO should increase if run anything other than 1D
  CoupledAggFact->SetMaxNeighAlreadySelected(maxNbrAlreadySelected);
  std::transform(aggOrdering.begin(), aggOrdering.end(), aggOrdering.begin(), ::tolower);
  if (aggOrdering == "natural" || aggOrdering == "graph" || aggOrdering == "random") {
    *out << "aggregate ordering :                    " << aggOrdering << std::endl;
    CoupledAggFact->SetOrdering(aggOrdering);
  } else {
    std::string msg = "main: bad aggregation option """ + aggOrdering + """.";
    throw(MueLu::Exceptions::RuntimeError(msg));
  }
  CoupledAggFact->SetPhase3AggCreation(0.5);
  Finest->Keep("Aggregates",CoupledAggFact.get());
  *out << "=============================================================================" << std::endl;

  // build transfer operators
//   RCP<NullspaceFactory> nspFact = rcp(new NullspaceFactory()); // make sure that we can keep nullspace!!!
//   RCP<TentativePFactory> TentPFact = rcp(new TentativePFactory(CoupledAggFact,nspFact));
//   //RCP<PgPFactory> Pfact = rcp( new PgPFactory(TentPFact) );
//   //RCP<FactoryBase2> Rfact  = rcp( new GenericRFactory(Pfact));
//   RCP<SaPFactory> Pfact  = rcp( new SaPFactory(TentPFact) );
//   RCP<FactoryBase2>   Rfact  = rcp( new TransPFactory(Pfact) );
//   RCP<RAPFactory> Acfact = rcp( new RAPFactory(Pfact, Rfact) );
//   Acfact->setVerbLevel(Teuchos::VERB_HIGH);

//   Finest->Keep("Aggregates",CoupledAggFact.get());
//   Finest->Keep("Nullspace",nspFact.get());

  // build level smoothers
  RCP<SmootherPrototype> smooProto;
  std::string ifpackType;
  Teuchos::ParameterList ifpackList;
  ifpackList.set("relaxation: sweeps", (LO) 1);
  ifpackList.set("relaxation: damping factor", (SC) 0.9); // 0.7
  ifpackType = "RELAXATION";
  ifpackList.set("relaxation: type", "Gauss-Seidel");

  smooProto = Teuchos::rcp( new TrilinosSmoother(ifpackType, ifpackList) );
  RCP<SmootherFactory> SmooFact;
  if (maxLevels > 1)
    SmooFact = rcp( new SmootherFactory(smooProto) );

  RCP<FactoryManager> M = rcp(new FactoryManager());
  M->SetFactory("Aggregates", CoupledAggFact);
  M->SetFactory("Smoother", SmooFact);

  H->Setup(*M,0,maxLevels);

  // print out aggregation information
  for(LocalOrdinal l=0; l<H->GetNumLevels()-1;l++) {
    RCP<Level> level = H->GetLevel((int)l);
    ExportAggregates(level, CoupledAggFact.get(),comm);
  }

  return EXIT_SUCCESS;
}
  void MueLuPreconditionerFactory<Scalar,LocalOrdinal,GlobalOrdinal,Node>::
  initializePrec(const RCP<const LinearOpSourceBase<Scalar> >& fwdOpSrc, PreconditionerBase<Scalar>* prec, const ESupportSolveUse supportSolveUse) const {
    using Teuchos::rcp_dynamic_cast;

    // we are using typedefs here, since we are using objects from different packages (Xpetra, Thyra,...)
    typedef Xpetra::Map<LocalOrdinal,GlobalOrdinal,Node>                     XpMap;
    typedef Xpetra::Operator<Scalar, LocalOrdinal, GlobalOrdinal, Node>      XpOp;
    typedef Xpetra::ThyraUtils<Scalar,LocalOrdinal,GlobalOrdinal,Node>       XpThyUtils;
    typedef Xpetra::CrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node>        XpCrsMat;
    typedef Xpetra::BlockedCrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node> XpBlockedCrsMat;
    typedef Xpetra::Matrix<Scalar,LocalOrdinal,GlobalOrdinal,Node>           XpMat;
    typedef Xpetra::MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node>      XpMultVec;
    typedef Xpetra::MultiVector<double,LocalOrdinal,GlobalOrdinal,Node>      XpMultVecDouble;
    typedef Thyra::LinearOpBase<Scalar>                                      ThyLinOpBase;
#ifdef HAVE_MUELU_TPETRA
    typedef MueLu::TpetraOperator<Scalar,LocalOrdinal,GlobalOrdinal,Node> MueTpOp;
    typedef Tpetra::Operator<Scalar,LocalOrdinal,GlobalOrdinal,Node>      TpOp;
    typedef Thyra::TpetraLinearOp<Scalar,LocalOrdinal,GlobalOrdinal,Node> ThyTpLinOp;
#endif

    // Check precondition
    TEUCHOS_ASSERT(Teuchos::nonnull(fwdOpSrc));
    TEUCHOS_ASSERT(this->isCompatible(*fwdOpSrc));
    TEUCHOS_ASSERT(prec);

    // Create a copy, as we may remove some things from the list
    ParameterList paramList = *paramList_;

    // Retrieve wrapped concrete Xpetra matrix from FwdOp
    const RCP<const ThyLinOpBase> fwdOp = fwdOpSrc->getOp();
    TEUCHOS_TEST_FOR_EXCEPT(Teuchos::is_null(fwdOp));

    // Check whether it is Epetra/Tpetra
    bool bIsEpetra  = XpThyUtils::isEpetra(fwdOp);
    bool bIsTpetra  = XpThyUtils::isTpetra(fwdOp);
    bool bIsBlocked = XpThyUtils::isBlockedOperator(fwdOp);
    TEUCHOS_TEST_FOR_EXCEPT((bIsEpetra == true  && bIsTpetra == true));
    TEUCHOS_TEST_FOR_EXCEPT((bIsEpetra == bIsTpetra) && bIsBlocked == false);
    TEUCHOS_TEST_FOR_EXCEPT((bIsEpetra != bIsTpetra) && bIsBlocked == true);

    RCP<XpMat> A = Teuchos::null;
    if(bIsBlocked) {
      Teuchos::RCP<const Thyra::BlockedLinearOpBase<Scalar> > ThyBlockedOp =
          Teuchos::rcp_dynamic_cast<const Thyra::BlockedLinearOpBase<Scalar> >(fwdOp);
      TEUCHOS_TEST_FOR_EXCEPT(Teuchos::is_null(ThyBlockedOp));

      TEUCHOS_TEST_FOR_EXCEPT(ThyBlockedOp->blockExists(0,0)==false);

      Teuchos::RCP<const LinearOpBase<Scalar> > b00 = ThyBlockedOp->getBlock(0,0);
      TEUCHOS_TEST_FOR_EXCEPT(Teuchos::is_null(b00));

      RCP<const XpCrsMat > xpetraFwdCrsMat00 = XpThyUtils::toXpetra(b00);
      TEUCHOS_TEST_FOR_EXCEPT(Teuchos::is_null(xpetraFwdCrsMat00));

      // MueLu needs a non-const object as input
      RCP<XpCrsMat> xpetraFwdCrsMatNonConst00 = Teuchos::rcp_const_cast<XpCrsMat>(xpetraFwdCrsMat00);
      TEUCHOS_TEST_FOR_EXCEPT(Teuchos::is_null(xpetraFwdCrsMatNonConst00));

      // wrap the forward operator as an Xpetra::Matrix that MueLu can work with
      RCP<XpMat> A00 = rcp(new Xpetra::CrsMatrixWrap<Scalar,LocalOrdinal,GlobalOrdinal,Node>(xpetraFwdCrsMatNonConst00));
      TEUCHOS_TEST_FOR_EXCEPT(Teuchos::is_null(A00));

      RCP<const XpMap> rowmap00 = A00->getRowMap();
      RCP< const Teuchos::Comm< int > > comm = rowmap00->getComm();

      // create a Xpetra::BlockedCrsMatrix which derives from Xpetra::Matrix that MueLu can work with
      RCP<XpBlockedCrsMat> bMat = Teuchos::rcp(new XpBlockedCrsMat(ThyBlockedOp, comm));
      TEUCHOS_TEST_FOR_EXCEPT(Teuchos::is_null(bMat));

      // save blocked matrix
      A = bMat;
    } else {
      RCP<const XpCrsMat > xpetraFwdCrsMat = XpThyUtils::toXpetra(fwdOp);
      TEUCHOS_TEST_FOR_EXCEPT(Teuchos::is_null(xpetraFwdCrsMat));

      // MueLu needs a non-const object as input
      RCP<XpCrsMat> xpetraFwdCrsMatNonConst = Teuchos::rcp_const_cast<XpCrsMat>(xpetraFwdCrsMat);
      TEUCHOS_TEST_FOR_EXCEPT(Teuchos::is_null(xpetraFwdCrsMatNonConst));

      // wrap the forward operator as an Xpetra::Matrix that MueLu can work with
      A = rcp(new Xpetra::CrsMatrixWrap<Scalar,LocalOrdinal,GlobalOrdinal,Node>(xpetraFwdCrsMatNonConst));
    }
    TEUCHOS_TEST_FOR_EXCEPT(Teuchos::is_null(A));

    // Retrieve concrete preconditioner object
    const Teuchos::Ptr<DefaultPreconditioner<Scalar> > defaultPrec = Teuchos::ptr(dynamic_cast<DefaultPreconditioner<Scalar> *>(prec));
    TEUCHOS_TEST_FOR_EXCEPT(Teuchos::is_null(defaultPrec));

    // extract preconditioner operator
    RCP<ThyLinOpBase> thyra_precOp = Teuchos::null;
    thyra_precOp = rcp_dynamic_cast<Thyra::LinearOpBase<Scalar> >(defaultPrec->getNonconstUnspecifiedPrecOp(), true);

    // Variable for multigrid hierarchy: either build a new one or reuse the existing hierarchy
    RCP<MueLu::Hierarchy<Scalar,LocalOrdinal,GlobalOrdinal,Node> > H = Teuchos::null;

    // make a decision whether to (re)build the multigrid preconditioner or reuse the old one
    // rebuild preconditioner if startingOver == true
    // reuse preconditioner if startingOver == false
    const bool startingOver = (thyra_precOp.is_null() || !paramList.isParameter("reuse: type") || paramList.get<std::string>("reuse: type") == "none");

    if (startingOver == true) {
      // extract coordinates from parameter list
      Teuchos::RCP<XpMultVecDouble> coordinates = Teuchos::null;
      coordinates = MueLu::Utilities<Scalar,LocalOrdinal,GlobalOrdinal,Node>::ExtractCoordinatesFromParameterList(paramList);

      // TODO check for Xpetra or Thyra vectors?
      RCP<XpMultVec> nullspace = Teuchos::null;
#ifdef HAVE_MUELU_TPETRA
      if (bIsTpetra) {
        typedef Tpetra::MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node> tMV;
        RCP<tMV> tpetra_nullspace = Teuchos::null;
        if (paramList.isType<Teuchos::RCP<tMV> >("Nullspace")) {
          tpetra_nullspace = paramList.get<RCP<tMV> >("Nullspace");
          paramList.remove("Nullspace");
          nullspace = MueLu::TpetraMultiVector_To_XpetraMultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node>(tpetra_nullspace);
          TEUCHOS_TEST_FOR_EXCEPT(Teuchos::is_null(nullspace));
        }
      }
#endif
      // build a new MueLu hierarchy
      H = MueLu::CreateXpetraPreconditioner(A, paramList, coordinates, nullspace);

    } else {
      // reuse old MueLu hierarchy stored in MueLu Tpetra/Epetra operator and put in new matrix

      // get old MueLu hierarchy
#if defined(HAVE_MUELU_TPETRA)
      if (bIsTpetra) {

        RCP<ThyTpLinOp> tpetr_precOp = rcp_dynamic_cast<ThyTpLinOp>(thyra_precOp);
        RCP<MueTpOp>    muelu_precOp = rcp_dynamic_cast<MueTpOp>(tpetr_precOp->getTpetraOperator(),true);

        H = muelu_precOp->GetHierarchy();
      }
#endif
      // TODO add the blocked matrix case here...

      TEUCHOS_TEST_FOR_EXCEPTION(!H->GetNumLevels(), MueLu::Exceptions::RuntimeError,
                                 "Thyra::MueLuPreconditionerFactory: Hierarchy has no levels in it");
      TEUCHOS_TEST_FOR_EXCEPTION(!H->GetLevel(0)->IsAvailable("A"), MueLu::Exceptions::RuntimeError,
                                 "Thyra::MueLuPreconditionerFactory: Hierarchy has no fine level operator");
      RCP<MueLu::Level> level0 = H->GetLevel(0);
      RCP<XpOp>    O0 = level0->Get<RCP<XpOp> >("A");
      RCP<XpMat>   A0 = rcp_dynamic_cast<XpMat>(O0);

      if (!A0.is_null()) {
        // If a user provided a "number of equations" argument in a parameter list
        // during the initial setup, we must honor that settings and reuse it for
        // all consequent setups.
        A->SetFixedBlockSize(A0->GetFixedBlockSize());
      }

      // set new matrix
      level0->Set("A", A);

      H->SetupRe();
    }

    // wrap hierarchy H in thyraPrecOp
    RCP<ThyLinOpBase > thyraPrecOp = Teuchos::null;
#if defined(HAVE_MUELU_TPETRA)
    if (bIsTpetra) {
      RCP<MueTpOp> muelu_tpetraOp = rcp(new MueTpOp(H));
      TEUCHOS_TEST_FOR_EXCEPT(Teuchos::is_null(muelu_tpetraOp));
      RCP<TpOp> tpOp = Teuchos::rcp_dynamic_cast<TpOp>(muelu_tpetraOp);
      thyraPrecOp = Thyra::createLinearOp<Scalar, LocalOrdinal, GlobalOrdinal, Node>(tpOp);
    }
#endif

    if(bIsBlocked) {
      TEUCHOS_TEST_FOR_EXCEPT(Teuchos::nonnull(thyraPrecOp));

      typedef MueLu::XpetraOperator<Scalar,LocalOrdinal,GlobalOrdinal,Node>    MueXpOp;
      //typedef Thyra::XpetraLinearOp<Scalar,LocalOrdinal,GlobalOrdinal,Node>    ThyXpLinOp; // unused
      const RCP<MueXpOp> muelu_xpetraOp = rcp(new MueXpOp(H));

      RCP<const VectorSpaceBase<Scalar> > thyraRangeSpace  = Xpetra::ThyraUtils<Scalar,LocalOrdinal,GlobalOrdinal,Node>::toThyra(muelu_xpetraOp->getRangeMap());
      RCP<const VectorSpaceBase<Scalar> > thyraDomainSpace = Xpetra::ThyraUtils<Scalar,LocalOrdinal,GlobalOrdinal,Node>::toThyra(muelu_xpetraOp->getDomainMap());

      RCP <Xpetra::Operator<Scalar, LocalOrdinal, GlobalOrdinal, Node> > xpOp = Teuchos::rcp_dynamic_cast<Xpetra::Operator<Scalar,LocalOrdinal,GlobalOrdinal,Node> >(muelu_xpetraOp);
      thyraPrecOp = Thyra::xpetraLinearOp(thyraRangeSpace, thyraDomainSpace,xpOp);
    }

    TEUCHOS_TEST_FOR_EXCEPT(Teuchos::is_null(thyraPrecOp));

    defaultPrec->initializeUnspecified(thyraPrecOp);

  }
示例#4
0
int main(int argc, char *argv[]) {

  // RCPs
  using Teuchos::RCP;
  using Teuchos::rcp;

  // MPI initialization using Teuchos
  Teuchos::GlobalMPISession mpiSession(&argc, &argv, NULL);
  RCP< const Teuchos::Comm<int> > comm = Teuchos::DefaultComm<int>::getComm();

  // predetermined DOF and parameter info
  // for case 1
  int nAcousticDOFs = 126;
  int nPaddedDOFs   = 3*nAcousticDOFs;
  int nElasticDOFs  = 243;
  int nTotalDOFs    = nPaddedDOFs+nElasticDOFs;
  int nDOFsPerNode  = 3;
  int nTotalNodes   = nTotalDOFs/nDOFsPerNode;
  int nnzeros_stiff = 12869;
  int nnzeros_damp  = 98;
  int nnzeros_mass  = 5635;
  int maxDOFsPerRow = 100;
  double freq = 11.0;
  double h=0.5;
  int nx=3;
  int ny=3;
  int nz=23;
  // for case 2
  /*int nAcousticDOFs = 600;
  int nPaddedDOFs   = 3*nAcousticDOFs;
  int nElasticDOFs  = 1425;
  int nTotalDOFs    = nPaddedDOFs+nElasticDOFs;
  int nDOFsPerNode  = 3;
  int nTotalNodes   = nTotalDOFs/nDOFsPerNode;
  int nnzeros_stiff = 94557;
  int nnzeros_damp  = 338;
  int nnzeros_mass  = 39715;
  int maxDOFsPerRow = 100;
  double freq = 22.0;
  double h=0.25;
  int nx=5;
  int ny=5;
  int nz=43;*/
  // for case 3
  /*int nAcousticDOFs = 3564;
  int nPaddedDOFs   = 3*nAcousticDOFs;
  int nElasticDOFs  = 9477;
  int nTotalDOFs    = nPaddedDOFs+nElasticDOFs;
  int nDOFsPerNode  = 3;
  int nTotalNodes   = nTotalDOFs/nDOFsPerNode;
  int nnzeros_stiff = 719287;
  int nnzeros_damp  = 1250;
  int nnzeros_mass  = 296875;
  int maxDOFsPerRow = 100;
  double freq = 44.0;
  double h=0.125;
  int nx=9;
  int ny=9;
  int nz=83;*/

  // Some constants
  SC omega  = (SC) 2.0*M_PI*freq;
  SC omega2 = omega*omega;
  SC one(1.0,0.0);
  SC ii(0.0,1.0);

  // Construct a Map that puts approximately the same number of mesh nodes per processor
  RCP<const Tpetra::Map<LO, GO, NO> > map = Tpetra::createUniformContigMap<LO, GO>(nTotalNodes, comm);
  // Tpetra map into Xpetra map
  RCP<const Map> xmap = Xpetra::toXpetra(map);
  // Map takes constant number of DOFs per node
  xmap = MapFactory::Build(xmap,nDOFsPerNode);
  map = Xpetra::toTpetra(xmap);

  // Create a CrsMatrix using the map
  // K - only stiffness matrix
  // M - only mass matrix
  // A - original Helmholtz operator
  // S - shifted Helmholtz operator
  RCP<TCRS> K = rcp(new TCRS(map,maxDOFsPerRow));
  RCP<TCRS> M = rcp(new TCRS(map,maxDOFsPerRow));
  RCP<TCRS> A = rcp(new TCRS(map,maxDOFsPerRow));
  RCP<TCRS> S = rcp(new TCRS(map,maxDOFsPerRow));

  // Read data from .txt files and put into appropriate matrices
  // Note: must pad acoustic DOFs with identity matrix
  // stiffness matrix
  std::ifstream matfile_stiff;
  matfile_stiff.open("coupled_stiff.txt");
  for (int i = 0; i < nnzeros_stiff; i++) {
    int current_row, current_column;
    double current_value;
    matfile_stiff >> current_row >> current_column >> current_value ;
    SC cpx_current_value(current_value,0.0);
    if(current_row < nAcousticDOFs)    { current_row    = current_row*3;                             }
    else                               { current_row    = current_row-nAcousticDOFs+nPaddedDOFs;     }
    if(current_column < nAcousticDOFs) { current_column = current_column*3;                          }
    else                               { current_column = current_column-nAcousticDOFs+nPaddedDOFs;  }
    if(map->isNodeGlobalElement(current_row)==true) {
      K->insertGlobalValues(current_row,
			    Teuchos::tuple<GO> (current_column),
			    Teuchos::tuple<SC> (cpx_current_value));
    }
  }
  // pad with identity
  for(int i = 0; i < nAcousticDOFs; i++) {
    if(map->isNodeGlobalElement(3*i+1)==true) {
      K->insertGlobalValues(3*i+1,
			    Teuchos::tuple<GO> (3*i+1),
			    Teuchos::tuple<SC> (one));
    }
    if(map->isNodeGlobalElement(3*i+2)==true) {
      K->insertGlobalValues(3*i+2,
			    Teuchos::tuple<GO> (3*i+2),
			    Teuchos::tuple<SC> (one));
    }
  }
  // damping matrix - lump into stiffness matrix
  std::ifstream matfile_damp;
  matfile_damp.open("coupled_damp.txt");
  for (int i = 0; i < nnzeros_damp; i++) {
    int current_row, current_column;
    double current_value;
    matfile_damp >> current_row >> current_column >> current_value ;
    SC cpx_current_value(current_value,0.0);
    if(current_row < nAcousticDOFs)    { current_row    = current_row*3;                             }
    else                               { current_row    = current_row-nAcousticDOFs+nPaddedDOFs;     }
    if(current_column < nAcousticDOFs) { current_column = current_column*3;                          }
    else                               { current_column = current_column-nAcousticDOFs+nPaddedDOFs;  }
    if(map->isNodeGlobalElement(current_row)==true) {
      K->insertGlobalValues(current_row,
			    Teuchos::tuple<GO> (current_column),
			    Teuchos::tuple<SC> (ii*omega*cpx_current_value));
    }
  }
  // mass matrix
  std::ifstream matfile_mass;
  matfile_mass.open("coupled_mass.txt");
  for (int i = 0; i < nnzeros_mass; i++) {
    int current_row, current_column;
    double current_value;
    matfile_mass >> current_row >> current_column >> current_value ;
    SC cpx_current_value(current_value,0.0);
    if(current_row < nAcousticDOFs)    { current_row    = current_row*3;                             }
    else                               { current_row    = current_row-nAcousticDOFs+nPaddedDOFs;     }
    if(current_column < nAcousticDOFs) { current_column = current_column*3;                          }
    else                               { current_column = current_column-nAcousticDOFs+nPaddedDOFs;  }
    if(map->isNodeGlobalElement(current_row)==true) {
      M->insertGlobalValues(current_row,
			     Teuchos::tuple<GO> (current_column),
			     Teuchos::tuple<SC> (cpx_current_value));
    }
  }
  // pad with identity
  for(int i = 0; i < nAcousticDOFs; i++) {
    if(map->isNodeGlobalElement(3*i+1)==true) {
      M->insertGlobalValues(3*i+1,
			    Teuchos::tuple<GO> (3*i+1),
			    Teuchos::tuple<SC> (one));
    }
    if(map->isNodeGlobalElement(3*i+2)==true) {
      M->insertGlobalValues(3*i+2,
			    Teuchos::tuple<GO> (3*i+2),
			    Teuchos::tuple<SC> (one));
    }
  }
  // Complete fill
  K->fillComplete();
  M->fillComplete();

  // Turn Tpetra::CrsMatrix into MueLu::Matrix
  RCP<XCRS> mueluK_ = rcp(new XTCRS(K));
  RCP<XMAT> mueluK  = rcp(new XWRAP(mueluK_));
  RCP<XCRS> mueluM_ = rcp(new XTCRS(M));
  RCP<XMAT> mueluM  = rcp(new XWRAP(mueluM_));
  mueluK->SetFixedBlockSize(nDOFsPerNode);
  mueluM->SetFixedBlockSize(nDOFsPerNode);
  // combine to make Helmholtz and shifted Laplace operators
  RCP<XMAT> mueluA, mueluS;
  SC shift1(1.0,0.5);
  MueLu::Utils2<SC,LO,GO,NO,LMO>::TwoMatrixAdd(mueluK, false, (SC) 1.0, mueluM, false, -omega2, mueluA);
  MueLu::Utils2<SC,LO,GO,NO,LMO>::TwoMatrixAdd(mueluK, false, (SC) 1.0, mueluM, false, -shift1*omega2, mueluS);
  mueluA->fillComplete();
  mueluS->fillComplete();
  xmap=mueluA->getDomainMap();
  map=Xpetra::toTpetra(xmap);

  // MultiVector of coordinates
  // NOTE: must be of size equal to the number of DOFs!
  RCP<XMV> coordinates;
  coordinates = MultiVectorFactory::Build(xmap, 3);
  for(int k=0; k<nz; k++) {
    for(int j=0; j<ny; j++) {
      for(int i=0; i<nx; i++) {
	int curidx = i+nx*j+nx*ny*k;
        int curidx0 = curidx*3+0;
        int curidx1 = curidx*3+1;
        int curidx2 = curidx*3+2;
	SC ih = (SC) (i*h);
	SC jh = (SC) (j*h);
	SC kh = (SC) (k*h);
	if(xmap->isNodeGlobalElement(curidx0)==true) {
	  coordinates->replaceGlobalValue(curidx0,0,ih);
	  coordinates->replaceGlobalValue(curidx0,1,jh);
	  coordinates->replaceGlobalValue(curidx0,2,kh);
	}
	if(xmap->isNodeGlobalElement(curidx1)==true) {
	  coordinates->replaceGlobalValue(curidx1,0,ih);
	  coordinates->replaceGlobalValue(curidx1,1,jh);
	  coordinates->replaceGlobalValue(curidx1,2,kh);
	}
	if(xmap->isNodeGlobalElement(curidx2)==true) {
	  coordinates->replaceGlobalValue(curidx2,0,ih);
	  coordinates->replaceGlobalValue(curidx2,1,jh);
	  coordinates->replaceGlobalValue(curidx2,2,kh);
	}
      }
    }
  }

  // Multigrid Hierarchy
  RCP<Hierarchy> H = rcp(new Hierarchy(mueluK));
  FactoryManager Manager;

  // Prolongation/Restriction
  RCP<TPFactory>  TentPFact    = rcp(  new TPFactory       );
  RCP<SaPFactory> Pfact        = rcp(  new SaPFactory      );
  RCP<GRFactory>  Rfact        = rcp(  new GRFactory       );
  RCP<RAPShiftFactory> Acfact  = rcp(  new RAPShiftFactory );
  RCP<RBMFactory> RBMfact      = rcp(  new RBMFactory(3)   );
  RBMfact->setLastAcousticDOF(nPaddedDOFs);

  // Smoothers
  RCP<SmootherPrototype> smooProto;
  std::string ifpack2Type;
  Teuchos::ParameterList ifpack2List;
  // Krylov smoother
  /*ifpack2Type = "KRYLOV";
  ifpack2List.set("krylov: iteration type",1);
  ifpack2List.set("krylov: number of iterations",4);
  ifpack2List.set("krylov: residual tolerance",1e-6);
  ifpack2List.set("krylov: block size",1);
  ifpack2List.set("krylov: zero starting solution",true);
  ifpack2List.set("krylov: preconditioner type",1);*/
  // Additive Schwarz smoother
  //ifpack2Type = "SCHWARZ";
  //ifpack2List.set("schwarz: compute condest", false);
  //ifpack2List.set("schwarz: combine mode", "Add"); // use string mode for this
  //ifpack2List.set("schwarz: reordering type", "none");
  //ifpack2List.set("schwarz: filter singletons", false);
  //ifpack2List.set("schwarz: overlap level", 0);
  // ILUT smoother
  //ifpack2Type = "ILUT";
  //ifpack2List.set("fact: ilut level-of-fill", (double)1.0);
  //ifpack2List.set("fact: absolute threshold", (double)0.0);
  //ifpack2List.set("fact: relative threshold", (double)1.0);
  //ifpack2List.set("fact: relax value", (double)0.0);
  // Gauss-Seidel smoother
  ifpack2Type = "RELAXATION";
  ifpack2List.set("relaxation: sweeps", (LO) 4);
  ifpack2List.set("relaxation: damping factor", (SC) 1.0); // 0.7
  ifpack2List.set("relaxation: type", "Gauss-Seidel");

  smooProto = Teuchos::rcp( new Ifpack2Smoother(ifpack2Type,ifpack2List) );
  RCP<SmootherFactory> SmooFact;
  LO maxLevels = 6;
  if (maxLevels > 1)
    SmooFact = rcp( new SmootherFactory(smooProto) );

  // create coarsest smoother
  RCP<SmootherPrototype> coarsestSmooProto;
  // Direct Solver
  std::string type = "";
  Teuchos::ParameterList coarsestSmooList;
#if defined(HAVE_AMESOS_SUPERLU)
  coarsestSmooProto = rcp( new DirectSolver("Superlu",coarsestSmooList) );
#else
  coarsestSmooProto = rcp( new DirectSolver("Klu",coarsestSmooList) );
#endif
  RCP<SmootherFactory> coarsestSmooFact = rcp(new SmootherFactory(coarsestSmooProto, Teuchos::null));

  RCP<XMV> nullspace;
  RBMfact->BuildRBM(mueluA,coordinates,nullspace);

  // determine shifts for RAPShiftFactory
  std::vector<SC> shifts;
  for(int i=0; i<maxLevels; i++) {
    double alpha=1.0;
    double beta=0.5+((double) i)*0.2;
    SC shift(alpha,beta);
    shifts.push_back(-shift*omega2);
  }
  Acfact->SetShifts(shifts);

  // Set factory managers for prolongation/restriction
  Manager.SetFactory("P", Pfact);
  Manager.SetFactory("R", Rfact);
  Manager.SetFactory("Ptent", TentPFact);
  H->Keep("P", Pfact.get());
  H->Keep("R", Rfact.get());
  H->Keep("Ptent", TentPFact.get());
  H->GetLevel(0)->Set("Nullspace",nullspace);
  H->SetImplicitTranspose(true);
  H->Setup(Manager, 0, maxLevels);
  H->Delete("Smoother");
  H->Delete("CoarseSolver");
  H->print(*getFancyOStream(Teuchos::rcpFromRef(std::cout)), MueLu::High);

  // Set factories for smoothing and coarse grid
  Manager.SetFactory("Smoother", SmooFact);
  Manager.SetFactory("CoarseSolver", coarsestSmooFact);
  Manager.SetFactory("A", Acfact);
  Manager.SetFactory("K", Acfact);
  Manager.SetFactory("M", Acfact);
  H->GetLevel(0)->Set("A",mueluS);
  H->GetLevel(0)->Set("K",mueluK);
  H->GetLevel(0)->Set("M",mueluM);
  H->Setup(Manager, 0, H->GetNumLevels());

  // right hand side and left hand side vectors
  RCP<TVEC> X = Tpetra::createVector<SC,LO,GO,NO>(map);
  RCP<TVEC> B = Tpetra::createVector<SC,LO,GO,NO>(map);
  X->putScalar((SC) 0.0);
  if(comm->getRank()==0)
    B->replaceGlobalValue(0, 1.0);

  // Define Operator and Preconditioner
  RCP<OP> belosOp   = rcp(new Belos::XpetraOp<SC,LO,GO,NO,LMO>(mueluA));   // Turns a Xpetra::Matrix object into a Belos operator
  RCP<OP> belosPrec = rcp(new Belos::MueLuOp<SC,LO,GO,NO,LMO>(H));         // Turns a MueLu::Hierarchy object into a Belos operator

  // Construct a Belos LinearProblem object
  RCP<Problem> belosProblem = rcp(new Problem(belosOp,X,B));
  belosProblem->setRightPrec(belosPrec);
  bool set = belosProblem->setProblem();
  if (set == false) {
    if(comm->getRank()==0) {
      std::cout << std::endl << "ERROR:  Belos::LinearProblem failed to set up correctly!" << std::endl;
    }
    return EXIT_FAILURE;
  }

  // Belos parameter list
  int maxIts = 100;
  double tol = 1e-6;
  Teuchos::ParameterList belosList;
  belosList.set("Maximum Iterations",    maxIts); // Maximum number of iterations allowed
  belosList.set("Convergence Tolerance", tol);    // Relative convergence tolerance requested
  belosList.set("Flexible Gmres", false);         // set flexible GMRES on

  // Create a FGMRES solver manager
  RCP<BelosSolver> solver = rcp( new BelosGMRES(belosProblem, rcp(&belosList, false)) );

  // Perform solve
  Belos::ReturnType ret = solver->solve();

  // print solution entries
  //using Teuchos::VERB_EXTREME;
  //Teuchos::RCP<Teuchos::FancyOStream> out = Teuchos::getFancyOStream( Teuchos::rcpFromRef(std::cerr) );
  //X->describe(*out,VERB_EXTREME);

  // Get the number of iterations for this solve.
  if(comm->getRank()==0) {
    std::cout << "Number of iterations performed for this solve: " << solver->getNumIters() << std::endl;
  }
  // Compute actual residuals.
  int numrhs=1;
  bool badRes = false;
  std::vector<double> actual_resids(numrhs);
  std::vector<double> rhs_norm(numrhs);
  RCP<TMV> resid = Tpetra::createMultiVector<SC,LO,GO,NO>(map, numrhs);
  OPT::Apply(*belosOp, *X, *resid);
  MVT::MvAddMv(-1.0, *resid, 1.0, *B, *resid);
  MVT::MvNorm(*resid, actual_resids);
  MVT::MvNorm(*B, rhs_norm);
  if(comm->getRank()==0) {
    std::cout<< "---------- Actual Residuals (normalized) ----------"<<std::endl<<std::endl;
  }
  for (int i = 0; i < numrhs; i++) {
    double actRes = abs(actual_resids[i])/rhs_norm[i];
    if(comm->getRank()==0) {
      std::cout <<"Problem " << i << " : \t" << actRes <<std::endl;
    }
    if (actRes > tol) { badRes = true; }
  }

  // Check convergence
  if (ret != Belos::Converged || badRes) {
    if(comm->getRank()==0) {
      std::cout << std::endl << "ERROR:  Belos did not converge! " << std::endl;
    }
    return EXIT_FAILURE;
  }
  if(comm->getRank()==0) {
    std::cout << std::endl << "SUCCESS:  Belos converged!" << std::endl;
  }

  return EXIT_SUCCESS;
}