std::pair<unsigned int, Real>
AztecLinearSolver<T>::solve (SparseMatrix<T>& matrix_in,
                             SparseMatrix<T>& precond_in,
                             NumericVector<T>& solution_in,
                             NumericVector<T>& rhs_in,
                             const double tol,
                             const unsigned int m_its)
{
  START_LOG("solve()", "AztecLinearSolver");

  // Make sure the data passed in are really of Epetra types
  EpetraMatrix<T>* matrix   = cast_ptr<EpetraMatrix<T>*>(&matrix_in);
  EpetraMatrix<T>* precond  = cast_ptr<EpetraMatrix<T>*>(&precond_in);
  EpetraVector<T>* solution = cast_ptr<EpetraVector<T>*>(&solution_in);
  EpetraVector<T>* rhs      = cast_ptr<EpetraVector<T>*>(&rhs_in);

  this->init();

  // Close the matrices and vectors in case this wasn't already done.
  matrix->close ();
  precond->close ();
  solution->close ();
  rhs->close ();

  _linear_solver->SetAztecOption(AZ_max_iter,m_its);
  _linear_solver->SetAztecParam(AZ_tol,tol);

  Epetra_FECrsMatrix * emat = matrix->mat();
  Epetra_Vector * esol = solution->vec();
  Epetra_Vector * erhs = rhs->vec();

  _linear_solver->Iterate(emat, esol, erhs, m_its, tol);

  STOP_LOG("solve()", "AztecLinearSolver");

  // return the # of its. and the final residual norm.
  return std::make_pair(_linear_solver->NumIters(), _linear_solver->TrueResidual());
}
Ejemplo n.º 2
0
//-----------------------------------------------------------------------------
void EpetraVector::gather(std::vector<double>& x,
                          const std::vector<dolfin::la_index>& indices) const
{
  const std::size_t _size = indices.size();
  x.resize(_size);
  dolfin_assert(x.size() == _size);

  // Gather values into a vector
  EpetraVector y;
  this->gather(y, indices);

  dolfin_assert(y.size() == _size);
  const Epetra_FEVector& _y = *(y.vec());

  // Copy values into x
  for (std::size_t i = 0; i < _size; ++i)
    x[i] = (_y)[0][i];
}
Ejemplo n.º 3
0
//-----------------------------------------------------------------------------
std::size_t EpetraKrylovSolver::solve(EpetraVector& x, const EpetraVector& b)
{
  dolfin_assert(solver);
  dolfin_assert(_A);
  dolfin_assert(_P);

  // Check dimensions
  const std::size_t M = _A->size(0);
  const std::size_t N = _A->size(1);
  if (N != b.size())
  {
    dolfin_error("EpetraKrylovSolver.cpp",
                 "solve linear system using Epetra Krylov solver",
                 "Non-matching dimensions for linear system");
  }

  // Write a message
  const bool report = parameters["report"];
  if (report && dolfin::MPI::process_number() == 0)
    info("Solving linear system of size %d x %d (Epetra Krylov solver).", M, N);

  // Reinitialize solution vector if necessary
  if (x.size() != M)
  {
    _A->resize(x, 1);
    x.zero();
  }
  else if (!parameters["nonzero_initial_guess"])
    x.zero();

  // Create linear problem
  Epetra_LinearProblem linear_problem(_A->mat().get(), x.vec().get(),
                                      b.vec().get());
  // Set-up linear solver
  solver->SetProblem(linear_problem);

  // Set output level
  if (parameters["monitor_convergence"])
  {
    const std::size_t interval = parameters["monitor_interval"];
    solver->SetAztecOption(AZ_output, interval);
  }
  else
    solver->SetAztecOption(AZ_output, AZ_none);

  dolfin_assert(_P);
  _preconditioner->set(*this, *_P);

  // Set covergence check method
  solver->SetAztecOption(AZ_conv, AZ_rhs);

  // Start solve
  solver->Iterate(static_cast<int>(parameters["maximum_iterations"]),
                  parameters["relative_tolerance"]);

  // Check solve status
  const double* status = solver->GetAztecStatus();
  if ((int) status[AZ_why] != AZ_normal)
  {
    std::string errorDescription;
    if( status[AZ_why] == AZ_maxits )
      errorDescription = "maximum iters reached";
    else if( status[AZ_why] == AZ_loss )
      errorDescription = "loss of accuracy";
    else if( status[AZ_why] == AZ_ill_cond  )
      errorDescription = "ill-conditioned";
    else if( status[AZ_why] == AZ_breakdown )
      errorDescription = "breakdown";
    else
      errorDescription = "unknown error";

    std::stringstream message;
    message << "Epetra (AztecOO) Krylov solver (" << _method << ", " << _preconditioner->name() << ") "
            << "failed to converge after " << (int)status[AZ_its] << " iterations "
            << "(" << errorDescription << ", error code " << (int)status[AZ_why] << ")";

    if (parameters["error_on_nonconvergence"])
    {
      dolfin_error("EpetraKrylovSolver.cpp",
                   "solve linear system using Epetra Krylov solver",
                   message.str());
    }
    else
      warning(message.str());
  }
  else
  {
    info("Epetra (AztecOO) Krylov solver (%s, %s) converged in %d iterations.",
          _method.c_str(), _preconditioner->name().c_str(), solver->NumIters());
  }

  // Update residuals
  absolute_residual = solver->TrueResidual();
  relative_residual = solver->ScaledResidual();

  // Return number of iterations
  return solver->NumIters();
}