PaStiX Handbook  6.3.2
Refinement internal function documentation

TODO. More...

Typedefs

typedef pastix_int_t(* refine_fct_t) (pastix_data_t *, pastix_rhs_t, pastix_rhs_t)
 Select the refinement function to call depending on the matrix type and the precision.
 

Functions

void d_refine_init (struct d_solver *solver, pastix_data_t *pastix_data)
 Initiate functions pointers to define basic operations. More...
 
pastix_int_t d_gmres_smp (pastix_data_t *pastix_data, pastix_rhs_t xp, pastix_rhs_t bp)
 
pastix_int_t d_grad_smp (pastix_data_t *pastix_data, pastix_rhs_t xp, pastix_rhs_t bp)
 
pastix_int_t d_pivot_smp (pastix_data_t *pastix_data, pastix_rhs_t xp, pastix_rhs_t bp)
 
pastix_int_t d_bicgstab_smp (pastix_data_t *pastix_data, pastix_rhs_t xp, pastix_rhs_t bp)
 
void c_refine_init (struct c_solver *solver, pastix_data_t *pastix_data)
 Initiate functions pointers to define basic operations. More...
 
pastix_int_t c_gmres_smp (pastix_data_t *pastix_data, pastix_rhs_t xp, pastix_rhs_t bp)
 
pastix_int_t c_grad_smp (pastix_data_t *pastix_data, pastix_rhs_t xp, pastix_rhs_t bp)
 
pastix_int_t c_pivot_smp (pastix_data_t *pastix_data, pastix_rhs_t xp, pastix_rhs_t bp)
 
pastix_int_t c_bicgstab_smp (pastix_data_t *pastix_data, pastix_rhs_t xp, pastix_rhs_t bp)
 
void s_refine_init (struct s_solver *solver, pastix_data_t *pastix_data)
 Initiate functions pointers to define basic operations. More...
 
pastix_int_t s_gmres_smp (pastix_data_t *pastix_data, pastix_rhs_t xp, pastix_rhs_t bp)
 
pastix_int_t s_grad_smp (pastix_data_t *pastix_data, pastix_rhs_t xp, pastix_rhs_t bp)
 
pastix_int_t s_pivot_smp (pastix_data_t *pastix_data, pastix_rhs_t xp, pastix_rhs_t bp)
 
pastix_int_t s_bicgstab_smp (pastix_data_t *pastix_data, pastix_rhs_t xp, pastix_rhs_t bp)
 
void z_refine_init (struct z_solver *solver, pastix_data_t *pastix_data)
 Initiate functions pointers to define basic operations. More...
 
pastix_int_t z_gmres_smp (pastix_data_t *pastix_data, pastix_rhs_t xp, pastix_rhs_t bp)
 
pastix_int_t z_grad_smp (pastix_data_t *pastix_data, pastix_rhs_t xp, pastix_rhs_t bp)
 
pastix_int_t z_pivot_smp (pastix_data_t *pastix_data, pastix_rhs_t xp, pastix_rhs_t bp)
 
pastix_int_t z_bicgstab_smp (pastix_data_t *pastix_data, pastix_rhs_t xp, pastix_rhs_t bp)
 
void z_refine_output_oneiter (pastix_fixdbl_t t0, pastix_fixdbl_t tf, double err, pastix_int_t nb_iters)
 Print statistics about one iteration. More...
 
void z_refine_output_final (pastix_data_t *pastix_data, pastix_complex64_t err, pastix_int_t nb_iters, pastix_fixdbl_t tf, void *x, pastix_complex64_t *gmresx)
 Final output. More...
 
void c_refine_output_oneiter (pastix_fixdbl_t t0, pastix_fixdbl_t tf, float err, pastix_int_t nb_iters)
 Print statistics about one iteration. More...
 
void c_refine_output_final (pastix_data_t *pastix_data, pastix_complex32_t err, pastix_int_t nb_iters, pastix_fixdbl_t tf, void *x, pastix_complex32_t *gmresx)
 Final output. More...
 
void d_refine_output_oneiter (pastix_fixdbl_t t0, pastix_fixdbl_t tf, double err, pastix_int_t nb_iters)
 Print statistics about one iteration. More...
 
void d_refine_output_final (pastix_data_t *pastix_data, double err, pastix_int_t nb_iters, pastix_fixdbl_t tf, void *x, double *gmresx)
 Final output. More...
 
void s_refine_output_oneiter (pastix_fixdbl_t t0, pastix_fixdbl_t tf, float err, pastix_int_t nb_iters)
 Print statistics about one iteration. More...
 
void s_refine_output_final (pastix_data_t *pastix_data, float err, pastix_int_t nb_iters, pastix_fixdbl_t tf, void *x, float *gmresx)
 Final output. More...
 

Detailed Description

TODO.

Function Documentation

◆ d_refine_init()

void d_refine_init ( struct d_solver *  solver,
pastix_data_t pastix_data 
)

◆ d_gmres_smp()

pastix_int_t d_gmres_smp ( pastix_data_t pastix_data,
pastix_rhs_t  xp,
pastix_rhs_t  bp 
)

d_gmres_smp - Function computing GMRES iterative refinement.

Parameters
[in]pastix_dataThe PaStiX data structure that describes the solver instance.
[out]xpThe solution vector.
[in]bpThe right hand side member (only one).
Returns
Number of iterations

H stores the h_{i,j} elements of the upper hessenberg matrix H (See Alg. 9.5 p 270) V stores the v_{i} vectors W stores the M^{-1} v_{i} vectors to avoid the application of the preconditioner on the output result (See line 11 of Alg 9.5)

If no preconditioner is applied, or the user wants to save memory, W stores only temporarily one vector for the Ax product (ldw is set to 0 to reuse the same vector at each iteration)

Algorithm from Iterative Methods for Sparse Linear systems, Y. Saad, Second Ed. p267-273

The version implemented is the Right preconditioned algorithm.

Compute x_m = x_0 + M^{-1} V_m y_m = x_0 + W_m y_m

Since we saved memory, we do not have (M^{-1} V_m) stored, thus we compute: w = V_m y_m w = M^{-1} (V_m y_m) x = x0 + (M^{-1} (V_m y_m))

Since we did not saved memory, we do have (M^{-1} V_m) stored in W_m if precond is true, thus we compute: x = x0 + W_m y_m, if precond x = x0 + V_m y_m, if not precond

Exit only if maximum number of iteration is reached. Exit on residual if checked at the beginning of the outer loop to be sure that the final residual of Ax-b is equal to the estimator computed within the inner loop.

Definition at line 49 of file d_refine_gmres.c.

References pastix_rhs_s::b, pastix_data_s::bcsc, d_refine_init(), pastix_data_s::dparm, DPARM_EPSILON_REFINEMENT, pastix_data_s::iparm, IPARM_GMRES_IM, IPARM_ITERMAX, IPARM_MIXED, IPARM_VERBOSE, pastix_int_t, PastixNoTrans, PastixVerboseNot, pastix_data_s::procnum, and pastix_data_s::steps.

◆ d_grad_smp()

pastix_int_t d_grad_smp ( pastix_data_t pastix_data,
pastix_rhs_t  xp,
pastix_rhs_t  bp 
)

d_grad_smp - Refine the solution using conjugate gradian method.

Parameters
[in]pastix_dataThe PaStiX data structure that describes the solver instance.
[out]xpThe solution vector.
[in]bpThe right hand side member (only one).
Returns
Number of iterations

Definition at line 49 of file d_refine_grad.c.

References pastix_rhs_s::b, pastix_data_s::bcsc, d_refine_init(), pastix_data_s::dparm, DPARM_EPSILON_REFINEMENT, pastix_data_s::iparm, IPARM_ITERMAX, IPARM_MIXED, IPARM_VERBOSE, pastix_int_t, PastixNoTrans, PastixVerboseNot, pastix_data_s::procnum, and pastix_data_s::steps.

◆ d_pivot_smp()

pastix_int_t d_pivot_smp ( pastix_data_t pastix_data,
pastix_rhs_t  xp,
pastix_rhs_t  bp 
)

d_pivot_smp - Refine the solution using static pivoting method.

Parameters
[in]pastix_dataThe PaStiX data structure that describes the solver instance.
[out]xpThe solution vector.
[in]bpThe right hand side member (only one).
Returns
Number of iterations

Definition at line 48 of file d_refine_pivot.c.

References pastix_rhs_s::b, pastix_data_s::bcsc, d_refine_init(), pastix_data_s::dparm, DPARM_EPSILON_REFINEMENT, pastix_data_s::iparm, IPARM_ITERMAX, IPARM_MIXED, IPARM_VERBOSE, pastix_int_t, PastixNoTrans, PastixVerboseNot, pastix_data_s::procnum, and pastix_data_s::steps.

◆ d_bicgstab_smp()

pastix_int_t d_bicgstab_smp ( pastix_data_t pastix_data,
pastix_rhs_t  xp,
pastix_rhs_t  bp 
)

d_bicgstab_smp - Function computing bicgstab iterative refinement.

Parameters
[in]pastix_dataThe PaStiX data structure that describes the solver instance.
[out]xpThe solution vector.
[in]bpThe right hand side member (only one).
Returns
Number of iterations

Definition at line 49 of file d_refine_bicgstab.c.

References pastix_rhs_s::b, pastix_data_s::bcsc, d_refine_init(), pastix_data_s::dparm, DPARM_EPSILON_REFINEMENT, pastix_data_s::iparm, IPARM_ITERMAX, IPARM_MIXED, IPARM_VERBOSE, pastix_int_t, PastixNoTrans, PastixVerboseNot, pastix_data_s::procnum, and pastix_data_s::steps.

◆ c_refine_init()

void c_refine_init ( struct c_solver *  solver,
pastix_data_t pastix_data 
)

◆ c_gmres_smp()

pastix_int_t c_gmres_smp ( pastix_data_t pastix_data,
pastix_rhs_t  xp,
pastix_rhs_t  bp 
)

c_gmres_smp - Function computing GMRES iterative refinement.

Parameters
[in]pastix_dataThe PaStiX data structure that describes the solver instance.
[out]xpThe solution vector.
[in]bpThe right hand side member (only one).
Returns
Number of iterations

H stores the h_{i,j} elements of the upper hessenberg matrix H (See Alg. 9.5 p 270) V stores the v_{i} vectors W stores the M^{-1} v_{i} vectors to avoid the application of the preconditioner on the output result (See line 11 of Alg 9.5)

If no preconditioner is applied, or the user wants to save memory, W stores only temporarily one vector for the Ax product (ldw is set to 0 to reuse the same vector at each iteration)

Algorithm from Iterative Methods for Sparse Linear systems, Y. Saad, Second Ed. p267-273

The version implemented is the Right preconditioned algorithm.

Compute x_m = x_0 + M^{-1} V_m y_m = x_0 + W_m y_m

Since we saved memory, we do not have (M^{-1} V_m) stored, thus we compute: w = V_m y_m w = M^{-1} (V_m y_m) x = x0 + (M^{-1} (V_m y_m))

Since we did not saved memory, we do have (M^{-1} V_m) stored in W_m if precond is true, thus we compute: x = x0 + W_m y_m, if precond x = x0 + V_m y_m, if not precond

Exit only if maximum number of iteration is reached. Exit on residual if checked at the beginning of the outer loop to be sure that the final residual of Ax-b is equal to the estimator computed within the inner loop.

Definition at line 49 of file c_refine_gmres.c.

References pastix_rhs_s::b, pastix_data_s::bcsc, c_refine_init(), pastix_data_s::dparm, DPARM_EPSILON_REFINEMENT, pastix_data_s::iparm, IPARM_GMRES_IM, IPARM_ITERMAX, IPARM_MIXED, pastix_int_t, PastixNoTrans, and pastix_data_s::steps.

◆ c_grad_smp()

pastix_int_t c_grad_smp ( pastix_data_t pastix_data,
pastix_rhs_t  xp,
pastix_rhs_t  bp 
)

c_grad_smp - Refine the solution using conjugate gradian method.

Parameters
[in]pastix_dataThe PaStiX data structure that describes the solver instance.
[out]xpThe solution vector.
[in]bpThe right hand side member (only one).
Returns
Number of iterations

Definition at line 49 of file c_refine_grad.c.

References pastix_rhs_s::b, pastix_data_s::bcsc, c_refine_init(), pastix_data_s::dparm, DPARM_EPSILON_REFINEMENT, pastix_data_s::iparm, IPARM_ITERMAX, IPARM_MIXED, IPARM_VERBOSE, pastix_int_t, PastixNoTrans, PastixVerboseNot, pastix_data_s::procnum, and pastix_data_s::steps.

◆ c_pivot_smp()

pastix_int_t c_pivot_smp ( pastix_data_t pastix_data,
pastix_rhs_t  xp,
pastix_rhs_t  bp 
)

c_pivot_smp - Refine the solution using static pivoting method.

Parameters
[in]pastix_dataThe PaStiX data structure that describes the solver instance.
[out]xpThe solution vector.
[in]bpThe right hand side member (only one).
Returns
Number of iterations

Definition at line 48 of file c_refine_pivot.c.

References pastix_rhs_s::b, pastix_data_s::bcsc, c_refine_init(), pastix_data_s::dparm, DPARM_EPSILON_REFINEMENT, pastix_data_s::iparm, IPARM_ITERMAX, IPARM_MIXED, IPARM_VERBOSE, pastix_int_t, PastixNoTrans, PastixVerboseNot, pastix_data_s::procnum, and pastix_data_s::steps.

◆ c_bicgstab_smp()

pastix_int_t c_bicgstab_smp ( pastix_data_t pastix_data,
pastix_rhs_t  xp,
pastix_rhs_t  bp 
)

c_bicgstab_smp - Function computing bicgstab iterative refinement.

Parameters
[in]pastix_dataThe PaStiX data structure that describes the solver instance.
[out]xpThe solution vector.
[in]bpThe right hand side member (only one).
Returns
Number of iterations

Definition at line 49 of file c_refine_bicgstab.c.

References pastix_rhs_s::b, pastix_data_s::bcsc, c_refine_init(), pastix_data_s::dparm, DPARM_EPSILON_REFINEMENT, pastix_data_s::iparm, IPARM_ITERMAX, IPARM_MIXED, IPARM_VERBOSE, pastix_int_t, PastixNoTrans, PastixVerboseNot, pastix_data_s::procnum, and pastix_data_s::steps.

◆ s_refine_init()

void s_refine_init ( struct s_solver *  solver,
pastix_data_t pastix_data 
)

◆ s_gmres_smp()

pastix_int_t s_gmres_smp ( pastix_data_t pastix_data,
pastix_rhs_t  xp,
pastix_rhs_t  bp 
)

s_gmres_smp - Function computing GMRES iterative refinement.

Parameters
[in]pastix_dataThe PaStiX data structure that describes the solver instance.
[out]xpThe solution vector.
[in]bpThe right hand side member (only one).
Returns
Number of iterations

H stores the h_{i,j} elements of the upper hessenberg matrix H (See Alg. 9.5 p 270) V stores the v_{i} vectors W stores the M^{-1} v_{i} vectors to avoid the application of the preconditioner on the output result (See line 11 of Alg 9.5)

If no preconditioner is applied, or the user wants to save memory, W stores only temporarily one vector for the Ax product (ldw is set to 0 to reuse the same vector at each iteration)

Algorithm from Iterative Methods for Sparse Linear systems, Y. Saad, Second Ed. p267-273

The version implemented is the Right preconditioned algorithm.

Compute x_m = x_0 + M^{-1} V_m y_m = x_0 + W_m y_m

Since we saved memory, we do not have (M^{-1} V_m) stored, thus we compute: w = V_m y_m w = M^{-1} (V_m y_m) x = x0 + (M^{-1} (V_m y_m))

Since we did not saved memory, we do have (M^{-1} V_m) stored in W_m if precond is true, thus we compute: x = x0 + W_m y_m, if precond x = x0 + V_m y_m, if not precond

Exit only if maximum number of iteration is reached. Exit on residual if checked at the beginning of the outer loop to be sure that the final residual of Ax-b is equal to the estimator computed within the inner loop.

Definition at line 49 of file s_refine_gmres.c.

References pastix_rhs_s::b, pastix_data_s::bcsc, pastix_data_s::dparm, DPARM_EPSILON_REFINEMENT, pastix_data_s::iparm, IPARM_GMRES_IM, IPARM_ITERMAX, IPARM_MIXED, IPARM_VERBOSE, pastix_int_t, PastixNoTrans, PastixVerboseNot, pastix_data_s::procnum, s_refine_init(), and pastix_data_s::steps.

◆ s_grad_smp()

pastix_int_t s_grad_smp ( pastix_data_t pastix_data,
pastix_rhs_t  xp,
pastix_rhs_t  bp 
)

s_grad_smp - Refine the solution using conjugate gradian method.

Parameters
[in]pastix_dataThe PaStiX data structure that describes the solver instance.
[out]xpThe solution vector.
[in]bpThe right hand side member (only one).
Returns
Number of iterations

Definition at line 49 of file s_refine_grad.c.

References pastix_rhs_s::b, pastix_data_s::bcsc, pastix_data_s::dparm, DPARM_EPSILON_REFINEMENT, pastix_data_s::iparm, IPARM_ITERMAX, IPARM_MIXED, IPARM_VERBOSE, pastix_int_t, PastixNoTrans, PastixVerboseNot, pastix_data_s::procnum, s_refine_init(), and pastix_data_s::steps.

◆ s_pivot_smp()

pastix_int_t s_pivot_smp ( pastix_data_t pastix_data,
pastix_rhs_t  xp,
pastix_rhs_t  bp 
)

s_pivot_smp - Refine the solution using static pivoting method.

Parameters
[in]pastix_dataThe PaStiX data structure that describes the solver instance.
[out]xpThe solution vector.
[in]bpThe right hand side member (only one).
Returns
Number of iterations

Definition at line 48 of file s_refine_pivot.c.

References pastix_rhs_s::b, pastix_data_s::bcsc, pastix_data_s::dparm, DPARM_EPSILON_REFINEMENT, pastix_data_s::iparm, IPARM_ITERMAX, IPARM_MIXED, IPARM_VERBOSE, pastix_int_t, PastixNoTrans, PastixVerboseNot, pastix_data_s::procnum, s_refine_init(), and pastix_data_s::steps.

◆ s_bicgstab_smp()

pastix_int_t s_bicgstab_smp ( pastix_data_t pastix_data,
pastix_rhs_t  xp,
pastix_rhs_t  bp 
)

s_bicgstab_smp - Function computing bicgstab iterative refinement.

Parameters
[in]pastix_dataThe PaStiX data structure that describes the solver instance.
[out]xpThe solution vector.
[in]bpThe right hand side member (only one).
Returns
Number of iterations

Definition at line 49 of file s_refine_bicgstab.c.

References pastix_rhs_s::b, pastix_data_s::bcsc, pastix_data_s::dparm, DPARM_EPSILON_REFINEMENT, pastix_data_s::iparm, IPARM_ITERMAX, IPARM_MIXED, IPARM_VERBOSE, pastix_int_t, PastixNoTrans, PastixVerboseNot, pastix_data_s::procnum, s_refine_init(), and pastix_data_s::steps.

◆ z_refine_init()

void z_refine_init ( struct z_solver *  solver,
pastix_data_t pastix_data 
)

◆ z_gmres_smp()

pastix_int_t z_gmres_smp ( pastix_data_t pastix_data,
pastix_rhs_t  xp,
pastix_rhs_t  bp 
)

z_gmres_smp - Function computing GMRES iterative refinement.

Parameters
[in]pastix_dataThe PaStiX data structure that describes the solver instance.
[out]xpThe solution vector.
[in]bpThe right hand side member (only one).
Returns
Number of iterations

H stores the h_{i,j} elements of the upper hessenberg matrix H (See Alg. 9.5 p 270) V stores the v_{i} vectors W stores the M^{-1} v_{i} vectors to avoid the application of the preconditioner on the output result (See line 11 of Alg 9.5)

If no preconditioner is applied, or the user wants to save memory, W stores only temporarily one vector for the Ax product (ldw is set to 0 to reuse the same vector at each iteration)

Algorithm from Iterative Methods for Sparse Linear systems, Y. Saad, Second Ed. p267-273

The version implemented is the Right preconditioned algorithm.

Compute x_m = x_0 + M^{-1} V_m y_m = x_0 + W_m y_m

Since we saved memory, we do not have (M^{-1} V_m) stored, thus we compute: w = V_m y_m w = M^{-1} (V_m y_m) x = x0 + (M^{-1} (V_m y_m))

Since we did not saved memory, we do have (M^{-1} V_m) stored in W_m if precond is true, thus we compute: x = x0 + W_m y_m, if precond x = x0 + V_m y_m, if not precond

Exit only if maximum number of iteration is reached. Exit on residual if checked at the beginning of the outer loop to be sure that the final residual of Ax-b is equal to the estimator computed within the inner loop.

Definition at line 49 of file z_refine_gmres.c.

References pastix_rhs_s::b, pastix_data_s::bcsc, pastix_data_s::dparm, DPARM_EPSILON_REFINEMENT, pastix_data_s::iparm, IPARM_GMRES_IM, IPARM_ITERMAX, IPARM_MIXED, pastix_int_t, PastixNoTrans, pastix_data_s::steps, and z_refine_init().

◆ z_grad_smp()

pastix_int_t z_grad_smp ( pastix_data_t pastix_data,
pastix_rhs_t  xp,
pastix_rhs_t  bp 
)

z_grad_smp - Refine the solution using conjugate gradian method.

Parameters
[in]pastix_dataThe PaStiX data structure that describes the solver instance.
[out]xpThe solution vector.
[in]bpThe right hand side member (only one).
Returns
Number of iterations

Definition at line 49 of file z_refine_grad.c.

References pastix_rhs_s::b, pastix_data_s::bcsc, pastix_data_s::dparm, DPARM_EPSILON_REFINEMENT, pastix_data_s::iparm, IPARM_ITERMAX, IPARM_MIXED, IPARM_VERBOSE, pastix_int_t, PastixNoTrans, PastixVerboseNot, pastix_data_s::procnum, pastix_data_s::steps, and z_refine_init().

◆ z_pivot_smp()

pastix_int_t z_pivot_smp ( pastix_data_t pastix_data,
pastix_rhs_t  xp,
pastix_rhs_t  bp 
)

z_pivot_smp - Refine the solution using static pivoting method.

Parameters
[in]pastix_dataThe PaStiX data structure that describes the solver instance.
[out]xpThe solution vector.
[in]bpThe right hand side member (only one).
Returns
Number of iterations

Definition at line 48 of file z_refine_pivot.c.

References pastix_rhs_s::b, pastix_data_s::bcsc, pastix_data_s::dparm, DPARM_EPSILON_REFINEMENT, pastix_data_s::iparm, IPARM_ITERMAX, IPARM_MIXED, IPARM_VERBOSE, pastix_int_t, PastixNoTrans, PastixVerboseNot, pastix_data_s::procnum, pastix_data_s::steps, and z_refine_init().

◆ z_bicgstab_smp()

pastix_int_t z_bicgstab_smp ( pastix_data_t pastix_data,
pastix_rhs_t  xp,
pastix_rhs_t  bp 
)

z_bicgstab_smp - Function computing bicgstab iterative refinement.

Parameters
[in]pastix_dataThe PaStiX data structure that describes the solver instance.
[out]xpThe solution vector.
[in]bpThe right hand side member (only one).
Returns
Number of iterations

Definition at line 49 of file z_refine_bicgstab.c.

References pastix_rhs_s::b, pastix_data_s::bcsc, pastix_data_s::dparm, DPARM_EPSILON_REFINEMENT, pastix_data_s::iparm, IPARM_ITERMAX, IPARM_MIXED, IPARM_VERBOSE, pastix_int_t, PastixNoTrans, PastixVerboseNot, pastix_data_s::procnum, pastix_data_s::steps, and z_refine_init().

◆ z_refine_output_oneiter()

void z_refine_output_oneiter ( pastix_fixdbl_t  t0,
pastix_fixdbl_t  tf,
double  err,
pastix_int_t  nb_iters 
)

Print statistics about one iteration.

Parameters
[in]t0The clock value at the beginning of the iteration
[in]tfThe clock value at the end of the iteration
[in]errThe backward error after the iteration
[in]nb_itersCurrent number of refinement iterations

Definition at line 52 of file z_refine_functions.c.

Referenced by z_refine_init().

◆ z_refine_output_final()

void z_refine_output_final ( pastix_data_t pastix_data,
pastix_complex64_t  err,
pastix_int_t  nb_iters,
pastix_fixdbl_t  tf,
void *  x,
pastix_complex64_t *  gmresx 
)

Final output.

Parameters
[in]pastix_dataThe PaStiX data structure that describes the solver instance.
[in]errThe final backward error
[in]nb_itersThe final number of iterations
[in]tfThe final clock value
[in,out]xThe vector that is to be overwritten by gmresx
[in]gmresxThe final solution

Definition at line 90 of file z_refine_functions.c.

Referenced by z_refine_init().

◆ c_refine_output_oneiter()

void c_refine_output_oneiter ( pastix_fixdbl_t  t0,
pastix_fixdbl_t  tf,
float  err,
pastix_int_t  nb_iters 
)

Print statistics about one iteration.

Parameters
[in]t0The clock value at the beginning of the iteration
[in]tfThe clock value at the end of the iteration
[in]errThe backward error after the iteration
[in]nb_itersCurrent number of refinement iterations

Definition at line 52 of file c_refine_functions.c.

Referenced by c_refine_init().

◆ c_refine_output_final()

void c_refine_output_final ( pastix_data_t pastix_data,
pastix_complex32_t  err,
pastix_int_t  nb_iters,
pastix_fixdbl_t  tf,
void *  x,
pastix_complex32_t gmresx 
)

Final output.

Parameters
[in]pastix_dataThe PaStiX data structure that describes the solver instance.
[in]errThe final backward error
[in]nb_itersThe final number of iterations
[in]tfThe final clock value
[in,out]xThe vector that is to be overwritten by gmresx
[in]gmresxThe final solution

Definition at line 90 of file c_refine_functions.c.

Referenced by c_refine_init().

◆ d_refine_output_oneiter()

void d_refine_output_oneiter ( pastix_fixdbl_t  t0,
pastix_fixdbl_t  tf,
double  err,
pastix_int_t  nb_iters 
)

Print statistics about one iteration.

Parameters
[in]t0The clock value at the beginning of the iteration
[in]tfThe clock value at the end of the iteration
[in]errThe backward error after the iteration
[in]nb_itersCurrent number of refinement iterations

Definition at line 52 of file d_refine_functions.c.

Referenced by d_refine_init().

◆ d_refine_output_final()

void d_refine_output_final ( pastix_data_t pastix_data,
double  err,
pastix_int_t  nb_iters,
pastix_fixdbl_t  tf,
void *  x,
double *  gmresx 
)

Final output.

Parameters
[in]pastix_dataThe PaStiX data structure that describes the solver instance.
[in]errThe final backward error
[in]nb_itersThe final number of iterations
[in]tfThe final clock value
[in,out]xThe vector that is to be overwritten by gmresx
[in]gmresxThe final solution

Definition at line 90 of file d_refine_functions.c.

Referenced by d_refine_init().

◆ s_refine_output_oneiter()

void s_refine_output_oneiter ( pastix_fixdbl_t  t0,
pastix_fixdbl_t  tf,
float  err,
pastix_int_t  nb_iters 
)

Print statistics about one iteration.

Parameters
[in]t0The clock value at the beginning of the iteration
[in]tfThe clock value at the end of the iteration
[in]errThe backward error after the iteration
[in]nb_itersCurrent number of refinement iterations

Definition at line 52 of file s_refine_functions.c.

Referenced by s_refine_init().

◆ s_refine_output_final()

void s_refine_output_final ( pastix_data_t pastix_data,
float  err,
pastix_int_t  nb_iters,
pastix_fixdbl_t  tf,
void *  x,
float *  gmresx 
)

Final output.

Parameters
[in]pastix_dataThe PaStiX data structure that describes the solver instance.
[in]errThe final backward error
[in]nb_itersThe final number of iterations
[in]tfThe final clock value
[in,out]xThe vector that is to be overwritten by gmresx
[in]gmresxThe final solution

Definition at line 90 of file s_refine_functions.c.

Referenced by s_refine_init().