https://mooseframework.inl.gov
VariableCondensationPreconditioner.h
Go to the documentation of this file.
1 //* This file is part of the MOOSE framework
2 //* https://mooseframework.inl.gov
3 //*
4 //* All rights reserved, see COPYRIGHT for full restrictions
5 //* https://github.com/idaholab/moose/blob/master/COPYRIGHT
6 //*
7 //* Licensed under LGPL 2.1, please see LICENSE for details
8 //* https://www.gnu.org/licenses/lgpl-2.1.html
9 
10 #pragma once
11 
12 // MOOSE includes
13 #include "MoosePreconditioner.h"
14 
15 // libMesh includes
16 #include "libmesh/preconditioner.h"
17 #include "libmesh/linear_implicit_system.h"
18 #include "libmesh/enum_preconditioner_type.h"
19 #include "libmesh/mesh_base.h"
20 #include "libmesh/petsc_matrix.h"
21 #include "libmesh/parallel_object.h"
22 
23 // C++ includes
24 #include <vector>
25 
26 // Forward declarations
29 
34  public libMesh::Preconditioner<Number>
35 {
36 public:
38 
41 
45  virtual void init();
46 
52  virtual void setup();
53 
58  virtual void apply(const NumericVector<Number> & y, NumericVector<Number> & x);
59 
63  virtual void clear();
64 
65 private:
69  void getDofToCondense();
70 
74  void getDofColRow();
75 
79  void findZeroDiagonals(SparseMatrix<Number> & mat, std::vector<dof_id_type> & indices);
80 
84  void condenseSystem();
85 
89  void preallocateCondensedJacobian(PetscMatrix<Number> & condensed_mat,
90  PetscMatrix<Number> & original_mat,
91  const std::vector<dof_id_type> & rows,
92  const std::vector<dof_id_type> & cols,
93  const std::vector<dof_id_type> & grows,
94  const std::vector<dof_id_type> & gcols,
95  PetscMatrix<Number> & block_mat);
96 
100  void computeCondensedJacobian(PetscMatrix<Number> & condensed_mat,
101  PetscMatrix<Number> & original_mat,
102  const std::vector<dof_id_type> & grows,
103  PetscMatrix<Number> & block_mat);
107  void computeDInverse(Mat & mat);
108 
113  void computeDInverseDiag(Mat & mat);
114 
119 
124 
129 
135  void mergeArrays(const PetscInt * a,
136  const PetscInt * b,
137  const PetscInt & na,
138  const PetscInt & nb,
139  std::vector<PetscInt> & c);
140 
143 
146 
149 
152 
155 
157  const unsigned int _n_vars;
158 
160  const std::vector<std::string> _lm_var_names;
161  std::vector<unsigned int> _lm_var_ids;
162 
164  const std::vector<std::string> _primary_var_names;
165  std::vector<unsigned int> _primary_var_ids;
166 
171  std::unique_ptr<PetscMatrix<Number>> _D, _M, _K;
173  Mat _dinv;
174 
176  std::unique_ptr<PetscMatrix<Number>> _J_condensed;
177 
181  std::unique_ptr<NumericVector<Number>> _x_hat, _y_hat, _primary_rhs_vec, _lm_sol_vec;
182 
185  std::vector<dof_id_type> _zero_rows;
186 
189  mutable bool _need_condense;
190 
193 
195  std::unique_ptr<Preconditioner<Number>> _preconditioner;
196 
206 
208  std::vector<dof_id_type> _global_rows, _rows, _global_cols, _cols;
209 
211  std::unordered_map<dof_id_type, dof_id_type> _global_rows_to_idx, _rows_to_idx,
213 
217  std::unordered_map<dof_id_type, dof_id_type> _map_global_lm_primary, _map_global_primary_order;
218 
222 };
std::vector< dof_id_type > _global_rows
row and column indices for the condensed system
void getDofToCondense()
Get dofs for the variable to be condensed out.
NonlinearSystemBase & _nl
The nonlinear system this PC is associated with (convenience reference)
std::unordered_map< dof_id_type, dof_id_type > _rows_to_idx
virtual void clear()
Release all memory and clear data structures.
void getCondensedXY(const NumericVector< Number > &y, NumericVector< Number > &x)
Get condensed x and y.
void preallocateCondensedJacobian(PetscMatrix< Number > &condensed_mat, PetscMatrix< Number > &original_mat, const std::vector< dof_id_type > &rows, const std::vector< dof_id_type > &cols, const std::vector< dof_id_type > &grows, const std::vector< dof_id_type > &gcols, PetscMatrix< Number > &block_mat)
Preallocate memory for the condensed Jacobian matrix.
std::unique_ptr< PetscMatrix< Number > > _J_condensed
Condensed Jacobian.
std::unique_ptr< PetscMatrix< Number > > _K
void mergeArrays(const PetscInt *a, const PetscInt *b, const PetscInt &na, const PetscInt &nb, std::vector< PetscInt > &c)
Find the common part of arrays a and b and save it in c.
The main MOOSE class responsible for handling user-defined parameters in almost every MOOSE system...
std::vector< dof_id_type > _global_lm_dofs
Vectors of DoFs: indices associated with lagrange multipliers, and its coupled primary variable globa...
std::unordered_map< dof_id_type, dof_id_type > _cols_to_idx
VariableCondensationPreconditioner(const InputParameters &params)
std::unique_ptr< PetscMatrix< Number > > _M
void computeDInverse(Mat &mat)
Compute inverse of D using LU.
virtual void apply(const NumericVector< Number > &y, NumericVector< Number > &x)
Computes the preconditioned vector "x" based on input "y".
unsigned int PerfID
Definition: MooseTypes.h:212
Base class for MOOSE preconditioners.
std::unique_ptr< NumericVector< Number > > _y_hat
Nonlinear system to be solved.
std::unique_ptr< Preconditioner< Number > > _preconditioner
Holds one Preconditioner object for the condensed system to solve.
bool _need_condense
Whether the DoFs associated the variable are to be condensed.
void getFullSolution(const NumericVector< Number > &y, NumericVector< Number > &x)
Assemble the full solution vector.
MooseMesh wraps a libMesh::Mesh object and enhances its capabilities by caching additional data and s...
Definition: MooseMesh.h:88
virtual void init()
Initialize data structures if not done so already.
virtual void setup()
This is called every time the "operator might have changed".
std::unordered_map< dof_id_type, dof_id_type > _global_cols_to_idx
const unsigned int _n_vars
Number of variables.
std::unordered_map< dof_id_type, dof_id_type > _map_global_primary_order
PreconditionerType
void findZeroDiagonals(SparseMatrix< Number > &mat, std::vector< dof_id_type > &indices)
Check if the original jacobian has zero diagonal entries and save the row indices.
const bool _adaptive_condensation
Whether to condense all specified variable.
void computeCondensedJacobian(PetscMatrix< Number > &condensed_mat, PetscMatrix< Number > &original_mat, const std::vector< dof_id_type > &grows, PetscMatrix< Number > &block_mat)
The condensed Jacobian matrix is computed in this function.
std::unique_ptr< NumericVector< Number > > _lm_sol_vec
void getDofColRow()
Get row and col dofs for the condensed system.
const std::vector< std::string > _primary_var_names
Name and ID of the corresponding coupled variable.
void computeDInverseDiag(Mat &mat)
Compute (approximate) inverse of D by inverting its diagonal entries.
void condenseSystem()
Reconstruct the equation system.
const bool _is_lm_coupling_diagonal
Whether the coupling is diagonal.
std::unordered_map< dof_id_type, dof_id_type > _global_rows_to_idx
Maps to keep track of row and col indices from the original Jacobian matrix to the condensed Jacobian...
std::unique_ptr< NumericVector< Number > > _x_hat
_x_hat, _y_hat: condensed solution and RHS vectors _primary_rhs_vec: part of the RHS vector that corr...
std::unique_ptr< NumericVector< Number > > _primary_rhs_vec
const std::vector< std::string > _lm_var_names
Name and ID of the variables that are to be condensed out (usually the Lagrange multiplier variable) ...
const libMesh::DofMap & _dofmap
DofMap for easy reference.
MooseMesh & _mesh
Mesh object for easy reference.
Interface for condensing out LMs for the dual mortar approach.
std::unordered_map< dof_id_type, dof_id_type > _map_global_lm_primary
Maps to keep track of the dof orders for keeping nonzero diagonal entries of the condensed system _ma...
libMesh::PreconditionerType _pre_type
Which preconditioner to use for the solve.
std::unique_ptr< PetscMatrix< Number > > _D
Submatrices that are frequently needed while computing the condensed system _D: the submatrix that co...
std::vector< dof_id_type > _zero_rows
The row indices that correspond to the zero diagonal entries in the original Jacobian matrix This is ...
void computeCondensedVariables()
Compute condensed variables (Lagrange multipliers) values using updated solution vector.