https://mooseframework.inl.gov
SegregatedSolverUtils.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 // Libmesh includes
13 #include "MooseTypes.h"
14 #include "MooseMesh.h"
15 #include "libmesh/petsc_vector.h"
16 #include "libmesh/petsc_matrix.h"
17 
18 namespace NS
19 {
20 namespace FV
21 {
34 void relaxMatrix(SparseMatrix<Number> & matrix_in,
35  const Real relaxation_parameter,
36  NumericVector<Number> & diff_diagonal);
37 
52  const NumericVector<Number> & solution_in,
53  const NumericVector<Number> & diff_diagonal);
54 
64  const NumericVector<Number> & vec_old,
65  const Real relaxation_factor);
66 
76  const Real min_limit = std::numeric_limits<Real>::epsilon(),
77  const Real max_limit = 1e10);
78 
93  const SparseMatrix<Number> & mat,
94  const NumericVector<Number> & rhs);
95 
106  NumericVector<Number> & rhs,
107  const Real desired_value,
108  const dof_id_type dof_id);
109 
118  const MooseMesh & mesh,
119  const Point & point);
120 
126 bool converged(const std::vector<std::pair<unsigned int, Real>> & residuals,
127  const std::vector<Real> & abs_tolerances);
128 } // End FV namespace
129 } // End Moose namespace
void relaxMatrix(SparseMatrix< Number > &matrix_in, const Real relaxation_parameter, NumericVector< Number > &diff_diagonal)
Relax the matrix to ensure diagonal dominance, we hold onto the difference in diagonals for later use...
Real computeNormalizationFactor(const NumericVector< Number > &solution, const SparseMatrix< Number > &mat, const NumericVector< Number > &rhs)
Compute a normalization factor which is applied to the linear residual to determine convergence...
MeshBase & mesh
void relaxRightHandSide(NumericVector< Number > &rhs_in, const NumericVector< Number > &solution_in, const NumericVector< Number > &diff_diagonal)
Relax the right hand side of an equation, this needs to be called once and the system matrix has been...
bool converged(const std::vector< std::pair< unsigned int, Real >> &residuals, const std::vector< Real > &abs_tolerances)
Based on the residuals, determine if the iterative process converged or not.
dof_id_type findPointDoFID(const MooseVariableFieldBase &variable, const MooseMesh &mesh, const Point &point)
Find the ID of the degree of freedom which corresponds to the variable and a given point on the mesh...
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
void constrainSystem(SparseMatrix< Number > &mx, NumericVector< Number > &rhs, const Real desired_value, const dof_id_type dof_id)
Implicitly constrain the system by adding a factor*(u-u_desired) to it at a desired dof value...
void limitSolutionUpdate(NumericVector< Number > &solution, const Real min_limit=std::numeric_limits< Real >::epsilon(), const Real max_limit=1e10)
Limit a solution to its minimum and maximum bounds: $u = min(max(u, min_limit), max_limit)$.
void relaxSolutionUpdate(NumericVector< Number > &vec_new, const NumericVector< Number > &vec_old, const Real relaxation_factor)
Relax the update on a solution field using the following approach: $u = u_{old}+ (u - u_{old})$...
uint8_t dof_id_type