Functions dedicated to to the linear algebra settings and operations in case of CDO cell-based schemes with a scalar-valued equations. More...
#include "cs_defs.h"
#include <assert.h>
#include <string.h>
#include <bft_mem.h>
#include "cs_array.h"
#include "cs_blas.h"
#include "cs_cdo_blas.h"
#include "cs_cdo_solve.h"
#include "cs_equation.h"
#include "cs_cdocb_scaleq.h"
#include "cs_fp_exception.h"
#include "cs_matrix_default.h"
#include "cs_parall.h"
#include "cs_saddle_solver.h"
#include "cs_saddle_system.h"
#include "cs_timer.h"
#include "cs_cdocb_scaleq_sles.h"
Functions | |
void | cs_cdocb_scaleq_sles_init_sharing (const cs_mesh_t *mesh, const cs_cdo_connect_t *connect, const cs_cdo_quantities_t *quant) |
Set pointers to shared structures. More... | |
void | cs_cdocb_scaleq_sles_init_system_helper (const cs_param_saddle_t *saddlep, cs_cdocb_scaleq_t *eqc) |
Define the system helper for a CDO-Cb scheme solving a scalar-valued equation (saddle-point system) More... | |
void | cs_cdocb_scaleq_sles_init_solver (const cs_equation_param_t *eqp, const cs_param_saddle_t *saddlep, cs_cdocb_scaleq_t *eqc) |
Define the saddle solver and its context for a CDO-Cb scheme solving a scalar-valued equation (saddle-point problem) More... | |
int | cs_cdocb_scaleq_sles_alu (cs_saddle_solver_t *solver, cs_real_t *flux, cs_real_t *pot) |
Solve a linear system arising from the discretization of a scalar-valed equation stemming from CDO cell-based schemes. Solve the saddle-point system using the Augmented Lagrangian-Uzawa algorithm. More... | |
int | cs_cdocb_scaleq_sles_block_krylov (cs_saddle_solver_t *solver, cs_real_t *flux, cs_real_t *pot) |
Solve a linear system arising from the discretization of a scalar-valed equation stemming from CDO cell-based schemes. The system is split into a flux block and the (unassembled) divergence operator. Block preconditioning using a Schur approximation on a Krylov solver such as the GCR or MINRES is available. More... | |
int | cs_cdocb_scaleq_sles_full_system (cs_saddle_solver_t *solver, cs_real_t *flux, cs_real_t *pot) |
Solve the saddle-point linear system arising from the discretization of the scalar-valued CDO cell-based scheme. The full system is treated as one block and solved as it is. In this situation, PETSc or MUMPS are usually considered. More... | |
int | cs_cdocb_scaleq_sles_gkb_inhouse (cs_saddle_solver_t *solver, cs_real_t *flux, cs_real_t *pot) |
Solve a linear system arising from the discretization of a scalar-valed equation stemming from CDO cell-based schemes. Solve this system using the Golub-Kahan Bidiagonalization algorithm. In-house implementation. The PETSc implementation is also available but appears less efficient in our tests. More... | |
int | cs_cdocb_scaleq_sles_notay (cs_saddle_solver_t *solver, cs_real_t *flux, cs_real_t *pot) |
Solve a linear system arising from the discretization of a scalar-valed equation stemming from CDO cell-based schemes. Solve this system using the Notay's algebraic transformation. The full system is treated as one block and solved as it is. In this situation, PETSc or MUMPS are usually considered. More... | |
int | cs_cdocb_scaleq_sles_uzawa_cg (cs_saddle_solver_t *solver, cs_real_t *flux, cs_real_t *pot) |
Solve a linear system arising from the discretization of a scalar-valed equation stemming from CDO cell-based schemes. Solve this system using the Uzawa-CG algorithm. More... | |
Functions dedicated to to the linear algebra settings and operations in case of CDO cell-based schemes with a scalar-valued equations.
int cs_cdocb_scaleq_sles_alu | ( | cs_saddle_solver_t * | solver, |
cs_real_t * | flux, | ||
cs_real_t * | pot | ||
) |
Solve a linear system arising from the discretization of a scalar-valed equation stemming from CDO cell-based schemes. Solve the saddle-point system using the Augmented Lagrangian-Uzawa algorithm.
[in,out] | solver | pointer to a cs_saddle_solver_t structure |
[in,out] | flux | values of the flux at faces (scalar) |
[in,out] | pot | values of the potential in cells |
int cs_cdocb_scaleq_sles_block_krylov | ( | cs_saddle_solver_t * | solver, |
cs_real_t * | flux, | ||
cs_real_t * | pot | ||
) |
Solve a linear system arising from the discretization of a scalar-valed equation stemming from CDO cell-based schemes. The system is split into a flux block and the (unassembled) divergence operator. Block preconditioning using a Schur approximation on a Krylov solver such as the GCR or MINRES is available.
[in,out] | solver | pointer to a cs_saddle_solver_t structure |
[in,out] | flux | values of the flux at faces (scalar) |
[in,out] | pot | values of the potential in cells |
int cs_cdocb_scaleq_sles_full_system | ( | cs_saddle_solver_t * | solver, |
cs_real_t * | flux, | ||
cs_real_t * | pot | ||
) |
Solve the saddle-point linear system arising from the discretization of the scalar-valued CDO cell-based scheme. The full system is treated as one block and solved as it is. In this situation, PETSc or MUMPS are usually considered.
[in,out] | solver | pointer to a saddle-point solver |
[in,out] | flux | values of the flux at faces (scalar) |
[in,out] | pot | values of the potential in cells |
int cs_cdocb_scaleq_sles_gkb_inhouse | ( | cs_saddle_solver_t * | solver, |
cs_real_t * | flux, | ||
cs_real_t * | pot | ||
) |
Solve a linear system arising from the discretization of a scalar-valed equation stemming from CDO cell-based schemes. Solve this system using the Golub-Kahan Bidiagonalization algorithm. In-house implementation. The PETSc implementation is also available but appears less efficient in our tests.
[in,out] | solver | pointer to a cs_saddle_solver_t structure |
[in,out] | flux | values of the flux at faces (scalar) |
[in,out] | pot | values of the potential in cells |
void cs_cdocb_scaleq_sles_init_sharing | ( | const cs_mesh_t * | mesh, |
const cs_cdo_connect_t * | connect, | ||
const cs_cdo_quantities_t * | quant | ||
) |
Set pointers to shared structures.
[in] | mesh | pointer to the mesh structure |
[in] | connect | pointer to additional CDO connectivities |
[in] | quant | pointer to additional CDO mesh quantities |
void cs_cdocb_scaleq_sles_init_solver | ( | const cs_equation_param_t * | eqp, |
const cs_param_saddle_t * | saddlep, | ||
cs_cdocb_scaleq_t * | eqc | ||
) |
Define the saddle solver and its context for a CDO-Cb scheme solving a scalar-valued equation (saddle-point problem)
[in] | eqp | set of equation parameters |
[in] | saddlep | parameters for solving a saddle-point problem |
[in,out] | eqc | pointer to a context structure cast on-the-fly |
void cs_cdocb_scaleq_sles_init_system_helper | ( | const cs_param_saddle_t * | saddlep, |
cs_cdocb_scaleq_t * | eqc | ||
) |
Define the system helper for a CDO-Cb scheme solving a scalar-valued equation (saddle-point system)
[in] | saddlep | parameters for solving a saddle-point problem |
[in,out] | eqc | pointer to a context structure cast on-the-fly |
int cs_cdocb_scaleq_sles_notay | ( | cs_saddle_solver_t * | solver, |
cs_real_t * | flux, | ||
cs_real_t * | pot | ||
) |
Solve a linear system arising from the discretization of a scalar-valed equation stemming from CDO cell-based schemes. Solve this system using the Notay's algebraic transformation. The full system is treated as one block and solved as it is. In this situation, PETSc or MUMPS are usually considered.
[in,out] | solver | pointer to a cs_saddle_solver_t structure |
[in,out] | flux | values of the flux at faces (scalar) |
[in,out] | pot | values of the potential in cells |
int cs_cdocb_scaleq_sles_uzawa_cg | ( | cs_saddle_solver_t * | solver, |
cs_real_t * | flux, | ||
cs_real_t * | pot | ||
) |
Solve a linear system arising from the discretization of a scalar-valed equation stemming from CDO cell-based schemes. Solve this system using the Uzawa-CG algorithm.
[in,out] | solver | pointer to a cs_saddle_solver_t structure |
[in,out] | flux | values of the flux at faces (scalar) |
[in,out] | pot | values of the potential in cells |