Linear Solvers and Systems Interface

Linear Systems

class LinearSystem

Subclassed by sierra::nalu::HypreLinearSystem, sierra::nalu::TpetraLinearSystem, sierra::nalu::TpetraSegregatedLinearSystem

Public Functions

virtual void buildDirichletNodeGraph(const stk::mesh::PartVector&)

Process nodes that belong to Dirichlet-type BC.

virtual void buildDirichletNodeGraph(const std::vector<stk::mesh::Entity>&)

Process nodes as belonging to a Dirichlet-type row.

See the documentation/implementation of sierra::nalu::FixPressureAtNodeAlgorithm for an example of this use case.

virtual void resetRows(const std::vector<stk::mesh::Entity> &nodeList, const unsigned beginPos, const unsigned endPos, const double diag_value = 0.0, const double rhs_residual = 0.0) = 0

Reset LHS and RHS for the given set of nodes to 0.

Parameters
  • nodeList: A list of STK node entities whose rows are zeroed out

  • beginPos: Starting index (usually 0)

  • endPos: Terminating index (1 for scalar quantities; nDim for vectors)

class TpetraLinearSystem : public sierra::nalu::LinearSystem

Public Functions

virtual void resetRows(const std::vector<stk::mesh::Entity> &nodeList, const unsigned beginPos, const unsigned endPos, const double diag_value = 0.0, const double rhs_residual = 0.0)

Reset LHS and RHS for the given set of nodes to 0.

Parameters
  • nodeList: A list of STK node entities whose rows are zeroed out

  • beginPos: Starting index (usually 0)

  • endPos: Terminating index (1 for scalar quantities; nDim for vectors)

class HypreLinearSystem : public sierra::nalu::LinearSystem

Nalu interface to populate a Hypre Linear System.

This class provides an interface to the HYPRE IJMatrix and IJVector data structures. It is responsible for creating, resetting, and destroying the Hypre data structures and provides the HypreLinearSystem::sumInto interface used by Nalu Kernels and SupplementalAlgorithms to populate entries into the linear system. The HypreLinearSystem::solve method interfaces with sierra::nalu::HypreDirectSolver that is responsible for the actual solution of the system using the required solver and preconditioner combination.

Subclassed by sierra::nalu::HypreUVWLinearSystem

Public Functions

HypreLinearSystem(Realm &realm, const unsigned numDof, EquationSystem *eqSys, LinearSolver *linearSolver)

Parameters
  • [in] realm: The realm instance that holds the EquationSystem being solved

  • [in] numDof: The degrees of freedom for the equation system created (Default: 1)

  • [in] eqSys: The equation system instance

  • [in] linearSolver: Handle to the HypreDirectSolver instance

virtual void buildDirichletNodeGraph(const stk::mesh::PartVector&)

Tag rows that must be handled as a Dirichlet BC node.

Parameters
  • [in] partVec: List of parts that contain the Dirichlet nodes

virtual void buildDirichletNodeGraph(const std::vector<stk::mesh::Entity>&)

Tag rows that must be handled as a Dirichlet node.

See

sierra::nalu::FixPressureAtNodeAlgorithm

Parameters
  • [in] entities: List of nodes where Dirichlet conditions are applied

virtual void loadComplete()

Finalize construction of the linear system matrix and rhs vector.

This method calls the appropriate Hypre functions to assemble the matrix and rhs in a parallel run, as well as registers the matrix and rhs with the solver preconditioner.

virtual void zeroSystem()

Reset the matrix and rhs data structures for the next iteration/timestep.

virtual int solve(stk::mesh::FieldBase *linearSolutionField)

Solve the system Ax = b.

The solution vector is returned in linearSolutionField

Parameters
  • [out] linearSolutionField: STK field where the solution is populated

double copy_hypre_to_stk(stk::mesh::FieldBase *)

Helper method to transfer the solution from a HYPRE_IJVector instance to the STK field data instance.

virtual void applyDirichletBCs(stk::mesh::FieldBase *solutionField, stk::mesh::FieldBase *bcValuesField, const stk::mesh::PartVector &parts, const unsigned beginPos, const unsigned endPos)

Populate the LHS and RHS for the Dirichlet rows in linear system.

virtual void sumInto(const std::vector<stk::mesh::Entity> &sym_meshobj, std::vector<int> &scratchIds, std::vector<double> &scratchVals, const std::vector<double> &rhs, const std::vector<double> &lhs, const char *trace_tag)

Update coefficients of a particular row(s) in the linear system.

The core method of this class, it updates the matrix and RHS based on the inputs from the various algorithms. Note that, unlike TpetraLinearSystem, this method skips over the fringe points of Overset mesh and the Dirichlet nodes rather than resetting them afterward.

This overloaded method deals with classic SupplementalAlgorithms

Parameters
  • [in] sym_meshobj: A list of STK node entities

  • [in] scratchIds: Work array for row IDs

  • [in] scratchVals: Work array for row entries

  • [in] rhs: Array containing RHS entries to be summed into [numEntities * numDof]

  • [in] lhs: Array containing LHS entries to be summed into. [numEntities * numDof * numEntities * numDof]

  • [in] trace_tag: Debugging message

virtual void resetRows(const std::vector<stk::mesh::Entity> &nodeList, const unsigned beginPos, const unsigned endPos, const double diag_value, const double rhs_residual)

Reset LHS and RHS for the given set of nodes to 0.

Parameters
  • nodeList: A list of STK node entities whose rows are zeroed out

  • beginPos: Starting index (usually 0)

  • endPos: Terminating index (1 for scalar quantities; nDim for vectors)

class HypreLinSysCoeffApplier : public sierra::nalu::CoeffApplier

Public Members

stk::mesh::NgpMesh ngpMesh_

mesh

NGPHypreIDFieldType ngpHypreGlobalId_

stk mesh field for the Hypre Global Id

unsigned numDof_ = 0

number of degrees of freedom

unsigned nDim_ = 0

number of rhs vectors

HypreIntType iLower_ = 0

The lowest row owned by this MPI rank.

HypreIntType iUpper_ = 0

The highest row owned by this MPI rank.

HypreIntType num_rows_owned_

Data structures for the owned CSR Matrix and RHS Vector(s)

HypreIntType num_rows_shared_

Data structures for the shared CSR Matrix and RHS Vector(s)

UnsignedViewRA mat_row_start_owned_ra_

Random access views.

PeriodicNodeMap periodic_node_to_hypre_id_

Auxilliary Data structures.

map of the periodic nodes to hypre ids

HypreIntTypeViewScalar checkSkippedRows_

Flag indicating that sumInto should check to see if rows must be skipped.

HypreIntTypeUnorderedMap skippedRowsMap_

unordered map for skipped rows

HypreIntTypeUnorderedMap oversetRowsMap_

unordered map for overset rows

HypreLinSysCoeffApplier *devicePointer_

this is the pointer to the device function … that assembles the lists

HypreIntType num_mat_overset_pts_owned_

number of points in the overset data structures

Linear Solvers Interface

class LinearSolver

An abstract representation of a linear solver in Nalu.

Defines the basic API supported by the linear solvers for use within Nalu. See concrete implementations such as sierra::nalu::TpetraLinearSolver for more details.

Subclassed by sierra::nalu::HypreDirectSolver

Public Functions

virtual PetraType getType() = 0

Type of solver instance as defined in sierra::nalu::PetraType.

virtual void destroyLinearSolver() = 0

Utility method to cleanup solvers during simulation.

bool &recomputePreconditioner()

Flag indicating whether the preconditioner is recomputed on each invocation.

bool &reusePreconditioner()

Flag indicating whether the preconditioner is reused on each invocation.

void zero_timer_precond()

Reset the preconditioner timer to 0.0 for future accumulation.

double get_timer_precond()

Get the preconditioner timer for the last invocation.

bool &activeMueLu()

Flag indicating whether the user has activated MueLU.

LinearSolverConfig *getConfig()

Get the solver configuration specified in the input file.

Public Members

std::string name_

User-friendly identifier for this particular solver instance.

Warning

doxygenclass: Cannot find class “sierra::nalu::TpetraLinearSolver” in doxygen xml output for project “nalu” from directory: ./doxygen/xml

class HypreDirectSolver : public sierra::nalu::LinearSolver

Nalu interface to Hypre Solvers and Preconditioners.

This class is responsible creation, initialization, execution, and clean up of Hypre solver and preconditioner data structures during the simulation. It provides an abstraction layer so that the user can choose different Hypre solvers via input parameters. This class interacts with rest of Nalu solely via sierra::nalu::HypreLinearSystem. The configuration of Hypre solver is controlled via user input parameters processed in sierra::nalu::HypreLinearSolverConfig

Users are referred to the Hypre Reference Manual for detailed documentation on the Hypre functions and data structures used in this class.

Subclassed by sierra::nalu::HypreUVWSolver

Public Functions

virtual void destroyLinearSolver()

Clean up Hypre data structures during simulation.

int solve(int&, double&, bool)

Solves the linear system and updates the solution vector.

Parameters
  • iters: The number of linear iterations performed

  • norm: The norm of the final relative residual

virtual PetraType getType()

Return the type of solver instance.

virtual void set_initialize_solver_flag()

public API for resetting the flag for how often the preconditioner is recomputed

Public Members

HYPRE_ParCSRMatrix parMat_

Instance of the Hypre parallel matrix.

HYPRE_ParVector parRhs_

Instance of the Hypre parallel RHS vector.

HYPRE_ParVector parSln_

Instance of Hypre parallel solution vector.

class LinearSolvers

Collection of solvers and their associated configuration.

This class performs the following actions within a Nalu simulation:

  • Parse the linear_solvers section and create a mapping of user-defined configurations.

  • Create solvers for specific equation system and update the mapping

Public Functions

void load(const YAML::Node &node)

Parse the linear_solvers section from Nalu input file.

LinearSolver *create_solver(std::string solverBlockName, const std::string realmName, EquationType theEQ)

Create a solver for the EquationSystem.

Parameters
  • [in] solverBlockName: The name specified in the input file, e.g., solve_scalar

  • [in] theEQ: The type of equation

Public Members

SolverMap solvers_

Mapping of solver instances to the EquationType.

SolverTpetraConfigMap solverTpetraConfig_

A lookup table of solver configurations against the names provided in the input file when the type is tpetra

HypreSolverConfigMap solverHypreConfig_

A lookup table of solver configurations against the names provided in the input file when type is hypre or tpetra_hypre

Simulation &sim_

Reference to the sierra::nalu::Simulation instance.

Solver Configuration

class LinearSolverConfig

Subclassed by sierra::nalu::HypreLinearSolverConfig, sierra::nalu::TpetraLinearSolverConfig

Public Functions

bool reuseLinSysIfPossible() const

User flag indicating whether equation systems must attempt to reuse linear system data structures even for cases with mesh motion.

This option only affects decoupled overset system solves where the matrix graph doesn’t change, only the entries within the graph. This can be controlled on a per-solver basis.

class TpetraLinearSolverConfig : public sierra::nalu::LinearSolverConfig
class HypreLinearSolverConfig : public sierra::nalu::LinearSolverConfig

User configuration parmeters for Hypre solvers and preconditioners.

Public Functions

virtual void load(const YAML::Node&)

Process and validate the user inputs and register calls to appropriate Hypre functions to configure the solver and preconditioner.