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" 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);
101 PetscMatrix<Number> & original_mat,
102 const std::vector<dof_id_type> & grows,
103 PetscMatrix<Number> & block_mat);
139 std::vector<PetscInt> & c);
171 std::unique_ptr<PetscMatrix<Number>>
_D,
_M,
_K;
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.
std::vector< dof_id_type > _global_lm_dofs
Vectors of DoFs: indices associated with lagrange multipliers, and its coupled primary variable globa...
std::vector< dof_id_type > _cols
std::unordered_map< dof_id_type, dof_id_type > _cols_to_idx
VariableCondensationPreconditioner(const InputParameters ¶ms)
std::unique_ptr< PetscMatrix< Number > > _M
std::vector< dof_id_type > _rows
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".
std::vector< dof_id_type > _global_primary_dofs
std::vector< unsigned int > _lm_var_ids
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.
virtual ~VariableCondensationPreconditioner()
bool _need_condense
Whether the DoFs associated the variable are to be condensed.
PerfID _init_timer
Timers.
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...
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
static InputParameters validParams()
const unsigned int _n_vars
Number of variables.
std::unordered_map< dof_id_type, dof_id_type > _map_global_primary_order
std::vector< dof_id_type > _global_cols
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.
std::vector< unsigned int > _primary_var_ids
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
std::vector< dof_id_type > _primary_dofs
std::vector< dof_id_type > _lm_dofs
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.