// ====================================================================== void Eig(const Operator& Op, MultiVector& ER, MultiVector& EI) { int ierr; if (Op.GetDomainSpace() != Op.GetRangeSpace()) ML_THROW("Matrix is not square", -1); ER.Reshape(Op.GetDomainSpace()); EI.Reshape(Op.GetDomainSpace()); Epetra_LinearProblem Problem; Problem.SetOperator(const_cast<Epetra_RowMatrix*>(Op.GetRowMatrix())); Amesos_Lapack Lapack(Problem); Epetra_Vector ER_Epetra(Op.GetRowMatrix()->RowMatrixRowMap()); Epetra_Vector EI_Epetra(Op.GetRowMatrix()->RowMatrixRowMap()); ierr = Lapack.GEEV(ER_Epetra, EI_Epetra); if (ierr) ML_THROW("GEEV returned error code = " + GetString(ierr), -1); for (int i = 0 ; i < ER.GetMyLength() ; ++i) { ER(i) = ER_Epetra[i]; EI(i) = EI_Epetra[i]; } }
// ====================================================================== // FIXME: Add List void Eigs(const Operator& A, int NumEigenvalues, MultiVector& ER, MultiVector& EI) { if (A.GetDomainSpace() != A.GetRangeSpace()) ML_THROW("Input Operator is not square", -1); double time; time = GetClock(); int length = NumEigenvalues; double tol = 1e-3; int restarts = 1; int output = 10; bool PrintStatus = true; // 1.- set parameters for Anasazi Teuchos::ParameterList AnasaziList; // MatVec should be either "A" or "ML^{-1}A" AnasaziList.set("eigen-analysis: matrix operation", "A"); AnasaziList.set("eigen-analysis: use diagonal scaling", false); AnasaziList.set("eigen-analysis: symmetric problem", false); AnasaziList.set("eigen-analysis: length", length); AnasaziList.set("eigen-analysis: block-size",1); AnasaziList.set("eigen-analysis: tolerance", tol); AnasaziList.set("eigen-analysis: restart", restarts); AnasaziList.set("eigen-analysis: output", output); AnasaziList.get("eigen-analysis: print current status",PrintStatus); // data to hold real and imag for eigenvalues and eigenvectors Space ESpace(-1, NumEigenvalues); ER.Reshape(ESpace); EI.Reshape(ESpace); // this is the starting value -- random Epetra_MultiVector EigenVectors(A.GetRowMatrix()->OperatorDomainMap(), NumEigenvalues); EigenVectors.Random(); #ifdef HAVE_ML_ANASAxI //int NumRealEigenvectors, NumImagEigenvectors; #endif AnasaziList.set("eigen-analysis: action", "LM"); #ifdef HAVE_ML_ANASAxI ML_THROW("fixme...", -1); /* FIXME ML_Anasazi::Interface(A.GetRowMatrix(),EigenVectors,ER.GetValues(), EI.GetValues(), AnasaziList, 0, 0, &NumRealEigenvectors, &NumImagEigenvectors, 0); */ #else ML_THROW("Anasazi is no longer supported", -1); #endif return; }
// ====================================================================== void Krylov(const Operator& A, const MultiVector& LHS, const MultiVector& RHS, const BaseOperator& Prec, Teuchos::ParameterList& List) { #ifndef HAVE_ML_AZTECOO std::cerr << "Please configure ML with --enable-aztecoo to use" << std::endl; std::cerr << "MLAPI Krylov solvers" << std::endl; exit(EXIT_FAILURE); #else if (LHS.GetNumVectors() != 1) ML_THROW("FIXME: only one vector is currently supported", -1); Epetra_LinearProblem Problem; const Epetra_RowMatrix& A_Epetra = *(A.GetRowMatrix()); Epetra_Vector LHS_Epetra(View,A_Epetra.OperatorDomainMap(), (double*)&(LHS(0))); Epetra_Vector RHS_Epetra(View,A_Epetra.OperatorRangeMap(), (double*)&(RHS(0))); // FIXME: this works only for Epetra-based operators Problem.SetOperator((const_cast<Epetra_RowMatrix*>(&A_Epetra))); Problem.SetLHS(&LHS_Epetra); Problem.SetRHS(&RHS_Epetra); AztecOO solver(Problem); EpetraBaseOperator Prec_Epetra(A_Epetra.OperatorDomainMap(),Prec); solver.SetPrecOperator(&Prec_Epetra); // get options from List int NumIters = List.get("krylov: max iterations", 1550); double Tol = List.get("krylov: tolerance", 1e-9); std::string type = List.get("krylov: type", "gmres"); int output = List.get("krylov: output level", GetPrintLevel()); // set options in `solver' if (type == "cg") solver.SetAztecOption(AZ_solver, AZ_cg); else if (type == "cg_condnum") solver.SetAztecOption(AZ_solver, AZ_cg_condnum); else if (type == "gmres") solver.SetAztecOption(AZ_solver, AZ_gmres); else if (type == "gmres_condnum") solver.SetAztecOption(AZ_solver, AZ_gmres_condnum); else if (type == "fixed point") solver.SetAztecOption(AZ_solver, AZ_fixed_pt); else ML_THROW("krylov: type has incorrect value (" + type + ")", -1); solver.SetAztecOption(AZ_output, output); solver.Iterate(NumIters, Tol); #endif }
// ====================================================================== void PrintSparsity(const Operator& A, int NumPDEEquations) { std::string FileName = A.GetLabel() + ".ps"; Ifpack_PrintSparsity(*(A.GetRowMatrix()), FileName.c_str(), NumPDEEquations); }
// ====================================================================== void AnalyzeCheap(const Operator& A) { Ifpack_Analyze(*(A.GetRowMatrix())); }