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; }