Functions dedicated to to the linear algebra settings and operations in case of CDO face-based schemes with a monolithic velocity-pressure coupling. More...
#include "cs_defs.h"
#include <stdio.h>
#include <stdlib.h>
#include <math.h>
#include <float.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_fp_exception.h"
#include "cs_matrix_default.h"
#include "cs_parall.h"
#include "cs_param_sles_setup.h"
#include "cs_saddle_solver_setup.h"
#include "cs_saddle_system.h"
#include "cs_timer.h"
#include "cs_cdofb_monolithic_sles.h"
Functions | |
void | cs_cdofb_monolithic_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_cdofb_monolithic_sles_init_system_helper (const cs_navsto_param_t *nsp, const cs_param_saddle_t *saddlep, cs_cdofb_monolithic_t *sc) |
Define the system helper for a CDO-Fb scheme solving the Navier-Stokes equation using a monolithic approach for the velocity-pressure coupling. More... | |
void | cs_cdofb_monolithic_sles_init_solver (const cs_navsto_param_t *nsp, const cs_param_saddle_t *saddlep, cs_cdofb_monolithic_t *sc) |
Define the saddle solver and its context for a CDO-Fb scheme solving the Navier-Stokes equation using a monolithic approach for the velocity-pressure coupling. More... | |
int | cs_cdofb_monolithic_sles_alu (const cs_navsto_param_t *nsp, cs_saddle_solver_t *solver, cs_real_t *u_f, cs_real_t *p_c) |
Solve a linear system arising from the discretization of the Navier-Stokes equation using a monolithic velocity-pressure coupling with a CDO face-based approach. Solve this system using the Augmented Lagrangian-Uzawa algorithm. More... | |
int | cs_cdofb_monolithic_sles_block_krylov (const cs_navsto_param_t *nsp, cs_saddle_solver_t *solver, cs_real_t *u_f, cs_real_t *p_c) |
Solve a linear system arising from the discretization of the Navier-Stokes equation with a CDO face-based approach. The system is split into a velocity 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_cdofb_monolithic_sles_full_system (const cs_navsto_param_t *nsp, cs_saddle_solver_t *solver, cs_real_t *u_f, cs_real_t *p_c) |
Solve a linear system arising from the discretization of the Navier-Stokes equation with a CDO face-based approach. 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_cdofb_monolithic_sles_gkb_inhouse (const cs_navsto_param_t *nsp, cs_saddle_solver_t *solver, cs_real_t *u_f, cs_real_t *p_c) |
Solve a linear system arising from the discretization of the Navier-Stokes equation using a monolithic velocity-pressure coupling with a CDO face-based approach. 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_cdofb_monolithic_sles_notay (const cs_navsto_param_t *nsp, cs_saddle_solver_t *solver, cs_real_t *u_f, cs_real_t *p_c) |
Solve a linear system arising from the discretization of the Navier-Stokes equation using a monolithic velocity-pressure coupling with a CDO face-based approach. 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_cdofb_monolithic_sles_uzawa_cg (const cs_navsto_param_t *nsp, cs_saddle_solver_t *solver, cs_real_t *u_f, cs_real_t *p_c) |
Solve a linear system arising from the discretization of the Navier-Stokes equation using a monolithic velocity-pressure coupling with a CDO face-based approach. Solve this system using the Uzawa-CG algorithm. More... | |
int | cs_cdofb_monolithic_sles_simple (const cs_navsto_param_t *nsp, cs_saddle_solver_t *solver, cs_real_t *u_f, cs_real_t *p_c) |
Solve a linear system arising from the discretization of the Navier-Stokes equation using a monolithic velocity-pressure coupling with a CDO face-based approach. Solve this system using the SIMPLE algorithm. More... | |
Functions dedicated to to the linear algebra settings and operations in case of CDO face-based schemes with a monolithic velocity-pressure coupling.
int cs_cdofb_monolithic_sles_alu | ( | const cs_navsto_param_t * | nsp, |
cs_saddle_solver_t * | solver, | ||
cs_real_t * | u_f, | ||
cs_real_t * | p_c | ||
) |
Solve a linear system arising from the discretization of the Navier-Stokes equation using a monolithic velocity-pressure coupling with a CDO face-based approach. Solve this system using the Augmented Lagrangian-Uzawa algorithm.
[in] | nsp | set of parameters related to the Navier-Stokes eqs. |
[in,out] | solver | pointer to a cs_saddle_solver_t structure |
[in,out] | u_f | values of the velocity at faces (3 components) |
[in,out] | p_c | values of the pressure in cells |
int cs_cdofb_monolithic_sles_block_krylov | ( | const cs_navsto_param_t * | nsp, |
cs_saddle_solver_t * | solver, | ||
cs_real_t * | u_f, | ||
cs_real_t * | p_c | ||
) |
Solve a linear system arising from the discretization of the Navier-Stokes equation with a CDO face-based approach. The system is split into a velocity 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] | nsp | set of parameters related to the Navier-Stokes eqs. |
[in,out] | solver | pointer to a saddle-point solver |
[in,out] | u_f | values of the velocity at faces (3 components) |
[in,out] | p_c | values of the pressure in cells |
int cs_cdofb_monolithic_sles_full_system | ( | const cs_navsto_param_t * | nsp, |
cs_saddle_solver_t * | solver, | ||
cs_real_t * | u_f, | ||
cs_real_t * | p_c | ||
) |
Solve a linear system arising from the discretization of the Navier-Stokes equation with a CDO face-based approach. The full system is treated as one block and solved as it is. In this situation, PETSc or MUMPS are usually considered.
[in] | nsp | set of parameters related to the Navier-Stokes eqs. |
[in,out] | solver | pointer to a saddle-point solver |
[in,out] | u_f | values of the velocity at faces (3 components) |
[in,out] | p_c | values of the pressure in cells |
int cs_cdofb_monolithic_sles_gkb_inhouse | ( | const cs_navsto_param_t * | nsp, |
cs_saddle_solver_t * | solver, | ||
cs_real_t * | u_f, | ||
cs_real_t * | p_c | ||
) |
Solve a linear system arising from the discretization of the Navier-Stokes equation using a monolithic velocity-pressure coupling with a CDO face-based approach. 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] | nsp | set of parameters related to the Navier-Stokes eqs. |
[in,out] | solver | pointer to a cs_saddle_solver_t structure |
[in,out] | u_f | values of the velocity at faces (3 components) |
[in,out] | p_c | values of the pressure in cells |
void cs_cdofb_monolithic_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_cdofb_monolithic_sles_init_solver | ( | const cs_navsto_param_t * | nsp, |
const cs_param_saddle_t * | saddlep, | ||
cs_cdofb_monolithic_t * | sc | ||
) |
Define the saddle solver and its context for a CDO-Fb scheme solving the Navier-Stokes equation using a monolithic approach for the velocity-pressure coupling.
[in] | nsp | set of parameters for the Navier-Stokes system |
[in] | saddlep | parameters for solving a saddle-point problem |
[in,out] | sc | pointer to a context structure cast on-the-fly |
void cs_cdofb_monolithic_sles_init_system_helper | ( | const cs_navsto_param_t * | nsp, |
const cs_param_saddle_t * | saddlep, | ||
cs_cdofb_monolithic_t * | sc | ||
) |
Define the system helper for a CDO-Fb scheme solving the Navier-Stokes equation using a monolithic approach for the velocity-pressure coupling.
[in] | nsp | Navier-Stokes paremeters |
[in] | saddlep | parameters for solving a saddle-point problem |
[in,out] | sc | pointer to a context structure cast on-the-fly |
int cs_cdofb_monolithic_sles_notay | ( | const cs_navsto_param_t * | nsp, |
cs_saddle_solver_t * | solver, | ||
cs_real_t * | u_f, | ||
cs_real_t * | p_c | ||
) |
Solve a linear system arising from the discretization of the Navier-Stokes equation using a monolithic velocity-pressure coupling with a CDO face-based approach. 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] | nsp | set of parameters related to the Navier-Stokes eqs. |
[in,out] | solver | pointer to a cs_saddle_solver_t structure |
[in,out] | u_f | values of the velocity at faces (3 components) |
[in,out] | p_c | values of the pressure in cells |
int cs_cdofb_monolithic_sles_simple | ( | const cs_navsto_param_t * | nsp, |
cs_saddle_solver_t * | solver, | ||
cs_real_t * | u_f, | ||
cs_real_t * | p_c | ||
) |
Solve a linear system arising from the discretization of the Navier-Stokes equation using a monolithic velocity-pressure coupling with a CDO face-based approach. Solve this system using the SIMPLE algorithm.
[in] | nsp | set of parameters related to the Navier-Stokes eqs. |
[in,out] | solver | pointer to a cs_saddle_solver_t structure |
[in,out] | u_f | values of the velocity at faces (3 components) |
[in,out] | p_c | values of the pressure in cells |
int cs_cdofb_monolithic_sles_uzawa_cg | ( | const cs_navsto_param_t * | nsp, |
cs_saddle_solver_t * | solver, | ||
cs_real_t * | u_f, | ||
cs_real_t * | p_c | ||
) |
Solve a linear system arising from the discretization of the Navier-Stokes equation using a monolithic velocity-pressure coupling with a CDO face-based approach. Solve this system using the Uzawa-CG algorithm.
[in] | nsp | set of parameters related to the Navier-Stokes eqs. |
[in,out] | solver | pointer to a cs_saddle_solver_t structure |
[in,out] | u_f | values of the velocity at faces (3 components) |
[in,out] | p_c | values of the pressure in cells |