Skip to content

Commit

Permalink
Address minor issues in NOX wrapper
Browse files Browse the repository at this point in the history
  • Loading branch information
jppelteret committed Jan 3, 2023
1 parent b5714a1 commit e0cb740
Show file tree
Hide file tree
Showing 6 changed files with 82 additions and 59 deletions.
2 changes: 2 additions & 0 deletions .clang-format
Original file line number Diff line number Diff line change
Expand Up @@ -117,6 +117,8 @@ IncludeCategories:
Priority: 300
- Regex: "deal.II/sundials/.*\\.h"
Priority: 310
- Regex: "deal.II/trilinos/.*\\.h"
Priority: 320
# put boost right after deal:
- Regex: "<boost.*>"
Priority: 500
Expand Down
1 change: 1 addition & 0 deletions doc/doxygen/options.dox.in
Original file line number Diff line number Diff line change
Expand Up @@ -240,6 +240,7 @@ PREDEFINED = DOXYGEN=1 \
DEAL_II_WITH_TRILINOS=1 \
DEAL_II_TRILINOS_WITH_EPETRAEXT=1 \
DEAL_II_TRILINOS_WITH_MUELU=1 \
DEAL_II_TRILINOS_WITH_NOX=1 \
DEAL_II_TRILINOS_WITH_ROL=1 \
DEAL_II_TRILINOS_WITH_SACADO=1 \
DEAL_II_TRILINOS_WITH_TPETRA=1 \
Expand Down
28 changes: 14 additions & 14 deletions include/deal.II/trilinos/nox.h
Original file line number Diff line number Diff line change
Expand Up @@ -38,7 +38,7 @@ namespace TrilinosWrappers

/**
* Wrapper around the nonlinear solver from the NOX
* packge (https://docs.trilinos.org/dev/packages/nox/doc/html/index.html),
* package (https://docs.trilinos.org/dev/packages/nox/doc/html/index.html),
* targeting deal.II data structures.
*
* The following code shows the steps how to use this class:
Expand All @@ -49,16 +49,15 @@ namespace TrilinosWrappers
* // create nonlinear solver
* TrilinosWrappers::NOXSolver<VectorType> solver(additional_data);
*
* // set user functions to compute residual,
* // set up Jacobian, and to apply the inverse of
* // Jacobian; note that there are more functions that
* // can be set
* // Set user functions to compute residual, to set up the Jacobian, and to
* // apply the inverse of the Jacobian.
* // Note that there are more functions that can be set.
* solver.residual = [](const auto &src, auto &dst) {...};
* solver.setup_jacobian = [](const auto &src) {...};
* solver.solve_with_jacobian =
* [](const auto &src, auto &dst, const auto) {...};
*
* // solver nonlinear system with solution containing the intial guess and
* // solver nonlinear system with solution containing the initial guess and
* // the final solution
* solver.solve(solution);
* @endcode
Expand Down Expand Up @@ -93,15 +92,15 @@ namespace TrilinosWrappers
/**
* Absolute l2 tolerance of the residual to be reached.
*
* @note Solver terminates successfully if either the absolut or
* @note Solver terminates successfully if either the absolute or
* the relative tolerance has been reached.
*/
double abs_tol;

/**
* Relative l2 tolerance of the residual to be reached.
*
* @note Solver terminates successfully if either the absolut or
* @note Solver terminates successfully if either the absolute or
* the relative tolerance has been reached.
*/
double rel_tol;
Expand All @@ -114,8 +113,9 @@ namespace TrilinosWrappers

/**
* Max number of linear iterations after which the preconditioner
* should be updated. This is only used if a lambda is attached to
* solve_with_jacobian_and_track_n_linear_iterations.
* should be updated. This is only used if
* solve_with_jacobian_and_track_n_linear_iterations has been given
* a target (i.e., it is not empty).
*/
unsigned int threshold_n_linear_iterations;

Expand Down Expand Up @@ -187,8 +187,8 @@ namespace TrilinosWrappers
*
* @note This function is optional and is used in the case of certain
* configurations. For instance, this function is required if the
* polynomial line search (@p NOX::LineSearch::Polynomial) is
* chosen, whereas for the full step case (@p NOX::LineSearch::FullStep)
* polynomial line search (`NOX::LineSearch::Polynomial`) is
* chosen, whereas for the full step case (`NOX::LineSearch::FullStep`)
* it won't be called.
*
* @note This function should return 0 in the case of success.
Expand All @@ -198,7 +198,7 @@ namespace TrilinosWrappers
/**
* A user function that applies the inverse of the Jacobian to
* @p x and writes the result in @p x. The parameter @p tolerance
* species the error reduction if a interative solver is used.
* specifies the error reduction if an iterative solver is used.
*
* @note This function is optional and is used in the case of certain
* configurations.
Expand Down Expand Up @@ -248,7 +248,7 @@ namespace TrilinosWrappers
* of linear iterations is exceeded.
*
* @note This function is optional. If no function is attached, this
* means implicitly a return value of false.
* means implicitly a return value of `false`.
*/
std::function<bool()> update_preconditioner_predicate;

Expand Down
86 changes: 53 additions & 33 deletions include/deal.II/trilinos/nox.templates.h
Original file line number Diff line number Diff line change
Expand Up @@ -18,10 +18,10 @@

#include <deal.II/base/config.h>

#include <deal.II/trilinos/nox.h>

#ifdef DEAL_II_TRILINOS_WITH_NOX

# include <deal.II/trilinos/nox.h>

# include <NOX_Abstract_Group.H>
# include <NOX_Abstract_Vector.H>
# include <NOX_Solver_Factory.H>
Expand Down Expand Up @@ -70,7 +70,7 @@ namespace TrilinosWrappers
}

/**
* Initialize every element of this vector with gamma.
* Initialize every element of this vector with @p gamma .
*/
NOX::Abstract::Vector &
init(double gamma) override
Expand All @@ -94,7 +94,7 @@ namespace TrilinosWrappers
}

/**
* Put element-wise absolute values of source vector y into this vector.
* Put element-wise absolute values of source vector @p y into this vector.
*/
NOX::Abstract::Vector &
abs(const NOX::Abstract::Vector &y) override
Expand All @@ -107,7 +107,7 @@ namespace TrilinosWrappers
}

/**
* Copy source vector y into this vector.
* Copy source vector @p y into this vector.
*/
NOX::Abstract::Vector &
operator=(const NOX::Abstract::Vector &y) override
Expand All @@ -127,7 +127,7 @@ namespace TrilinosWrappers
}

/**
* Put element-wise reciprocal of source vector y into this vector.
* Put element-wise reciprocal of source vector @p y into this vector.
*/
NOX::Abstract::Vector &
reciprocal(const NOX::Abstract::Vector &y) override
Expand All @@ -140,7 +140,7 @@ namespace TrilinosWrappers
}

/**
* Scale each element of this vector by gamma.
* Scale each element of this vector by @p gamma .
*/
NOX::Abstract::Vector &
scale(double gamma) override
Expand All @@ -151,7 +151,7 @@ namespace TrilinosWrappers
}

/**
* Scale this vector element-by-element by the vector a.
* Scale this vector element-by-element by the vector @p a.
*/
NOX::Abstract::Vector &
scale(const NOX::Abstract::Vector &a) override
Expand All @@ -166,7 +166,7 @@ namespace TrilinosWrappers
}

/**
* Compute x = (alpha * a) + (gamma * x) where x is this vector.
* Compute `x = (alpha * a) + (gamma * x)` where `x` is this vector.
*/
NOX::Abstract::Vector &
update(double alpha,
Expand All @@ -185,8 +185,8 @@ namespace TrilinosWrappers
}

/**
* Compute x = (alpha * a) + (beta * b) + (gamma * x) where x is this
* vector.
* Compute `x = (alpha * a) + (beta * b) + (gamma * x)` where `x` is
* this vector.
*/
NOX::Abstract::Vector &
update(double alpha,
Expand Down Expand Up @@ -250,7 +250,7 @@ namespace TrilinosWrappers
}

/**
* Return the vector's weighted 2-Norm.
* Return the vector's weighted 2-norm.
*/
double
norm(const NOX::Abstract::Vector &weights) const override
Expand Down Expand Up @@ -410,7 +410,7 @@ namespace TrilinosWrappers
}

/**
* Set the solution vector `x` to @y.
* Set the solution vector `x` to @p y.
*/
void
setX(const NOX::Abstract::Vector &y) override
Expand Down Expand Up @@ -458,7 +458,7 @@ namespace TrilinosWrappers
}

/**
* Return true if the residual vector `F` is valid.
* Return `true` if the residual vector `F` is valid.
*/
bool
isF() const override
Expand All @@ -484,7 +484,7 @@ namespace TrilinosWrappers
}

/**
* Return true if the Jacobian is valid.
* Return `true` if the Jacobian is valid.
*/
bool
isJacobian() const override
Expand Down Expand Up @@ -538,7 +538,7 @@ namespace TrilinosWrappers
}

/**
* Return an `RCP` to solution vector.
* Return a reference counting pointer to solution vector.
*/
Teuchos::RCP<const NOX::Abstract::Vector>
getXPtr() const override
Expand All @@ -548,7 +548,7 @@ namespace TrilinosWrappers
}

/**
* Return an RCP to the residual `F(x)`.
* Return a reference counting pointer to the residual `F(x)`.
*/
Teuchos::RCP<const NOX::Abstract::Vector>
getFPtr() const override
Expand All @@ -558,7 +558,7 @@ namespace TrilinosWrappers
}

/**
* Return an `RCP` to gradient.
* Return a reference counting pointer to gradient.
*/
Teuchos::RCP<const NOX::Abstract::Vector>
getGradientPtr() const override
Expand All @@ -568,7 +568,7 @@ namespace TrilinosWrappers
}

/**
* Return an `RCP` to the Newton descent direction.
* Return a reference counting pointer to the Newton descent direction.
*/
Teuchos::RCP<const NOX::Abstract::Vector>
getNewtonPtr() const override
Expand All @@ -579,8 +579,8 @@ namespace TrilinosWrappers

/**
* Create a new Group of the same derived type as this one by
* cloning this one, and return a reference counting
* pointer to the new group.
* cloning this one, and return a reference counting pointer to
* the new group.
*/
Teuchos::RCP<NOX::Abstract::Group>
clone(NOX::CopyType copy_type) const override
Expand Down Expand Up @@ -668,7 +668,7 @@ namespace TrilinosWrappers
}

/**
* Applies the Jacobian to the given input vector and assigns
* Applies the Jacobian to the given @p input vector and assigns
* the output to the @p result.
*/
NOX::Abstract::Group::ReturnType
Expand All @@ -693,7 +693,7 @@ namespace TrilinosWrappers

private:
/**
* Reset state.
* Reset the state of this object.
*/
void
reset()
Expand All @@ -702,36 +702,56 @@ namespace TrilinosWrappers
is_valid_j = false;
}

// internal vector for the current solution
/**
* An internal vector for the current solution.
*/
Vector<VectorType> x;

// internal vector for the residual
/**
* An internal vector for the residual.
*/
Vector<VectorType> f;

// internal vector for the solution gradient
/**
* An internal vector for the solution gradient.
*/
Vector<VectorType> gradient;

// internal vector for the newton step
/**
* An internal vector for the newton step.
*/
Vector<VectorType> newton;

// helper function to compute residual
/**
* A helper function to compute residual.
*/
std::function<int(const VectorType &x, VectorType &f)> residual;

// helper function to setup Jacobian
/**
* A helper function to setup Jacobian.
*/
std::function<int(const VectorType &x)> setup_jacobian;

// helper function to apply Jacobian
/**
* A helper function to apply Jacobian.
*/
std::function<int(const VectorType &x, VectorType &v)> apply_jacobian;

// helper function to solve jacobian
/**
* A helper function to solve Jacobian.
*/
std::function<
int(const VectorType &f, VectorType &x, const double tolerance)>
solve_with_jacobian;

// internal state (is residual computed?)
/**
* A flag that indicates if the has residual been computed.
*/
bool is_valid_f;

// internal state (is Jacobian computed?)
/**
* A flag that indicates if the Jacobian has been computed.
*/
bool is_valid_j;
};

Expand Down
24 changes: 12 additions & 12 deletions source/trilinos/nox.cc
Original file line number Diff line number Diff line change
Expand Up @@ -17,20 +17,20 @@

#include <deal.II/base/config.h>

#include <deal.II/lac/block_vector.h>
#include <deal.II/lac/la_parallel_block_vector.h>
#include <deal.II/lac/la_parallel_vector.h>
#include <deal.II/lac/la_vector.h>
#include <deal.II/lac/petsc_block_vector.h>
#include <deal.II/lac/petsc_vector.h>
#include <deal.II/lac/trilinos_parallel_block_vector.h>
#include <deal.II/lac/trilinos_vector.h>
#include <deal.II/lac/vector_memory.h>

#include <deal.II/trilinos/nox.templates.h>

#ifdef DEAL_II_TRILINOS_WITH_NOX

# include <deal.II/lac/block_vector.h>
# include <deal.II/lac/la_parallel_block_vector.h>
# include <deal.II/lac/la_parallel_vector.h>
# include <deal.II/lac/la_vector.h>
# include <deal.II/lac/petsc_block_vector.h>
# include <deal.II/lac/petsc_vector.h>
# include <deal.II/lac/trilinos_parallel_block_vector.h>
# include <deal.II/lac/trilinos_vector.h>
# include <deal.II/lac/vector_memory.h>

# include <deal.II/trilinos/nox.templates.h>

DEAL_II_NAMESPACE_OPEN

namespace TrilinosWrappers
Expand Down

0 comments on commit e0cb740

Please sign in to comment.