1#ifndef __CS_SLES_IT_PRIV_H__
2#define __CS_SLES_IT_PRIV_H__
87#define HUGE_VAL 1.E+12
128typedef struct _cs_sles_it_setup_t {
130 double initial_residual;
146typedef struct _cs_sles_it_add_t {
155struct _cs_sles_it_t {
164 bool ignore_convergence;
168 int restart_interval;
173 cs_sles_it_solve_t *solve;
185 unsigned n_iterations_last;
187 unsigned n_iterations_min;
189 unsigned n_iterations_max;
191 unsigned long long n_iterations_tot;
207# if defined(HAVE_MPI)
209 MPI_Comm caller_comm;
215 const struct _cs_sles_it_t *shared;
219 cs_sles_it_add_t *add_data;
221 cs_sles_it_setup_t *setup_data;
227 int fallback_n_max_iter;
238struct _cs_sles_it_convergence_t {
244 unsigned n_iterations;
245 unsigned n_iterations_max;
274 double s =
cs_dot(c->setup_data->n_rows, x, y);
278 if (c->comm != MPI_COMM_NULL) {
280 MPI_Allreduce(&s, &_sum, 1, MPI_DOUBLE, MPI_SUM, c->comm);
310 if (c->comm != MPI_COMM_NULL) {
312 MPI_Allreduce(&s, &_sum, 1, MPI_DOUBLE, MPI_SUM, c->comm);
345 if (c->comm != MPI_COMM_NULL) {
347 MPI_Allreduce(s, _sum, 2, MPI_DOUBLE, MPI_SUM, c->comm);
384 if (c->comm != MPI_COMM_NULL) {
386 MPI_Allreduce(s, _sum, 2, MPI_DOUBLE, MPI_SUM, c->comm);
425 if (c->comm != MPI_COMM_NULL) {
428 MPI_Allreduce(s, _sum, 3, MPI_DOUBLE, MPI_SUM, c->comm);
473 if (c->comm != MPI_COMM_NULL) {
475 MPI_Allreduce(s, _sum, 5, MPI_DOUBLE, MPI_SUM, c->comm);
476 memcpy(s, _sum, 5*
sizeof(
double));
507 aux[0] = (c[0] - b[0]);
508 aux[1] = (c[1] - b[1]) - aux[0]*mat[3];
509 aux[2] = (c[2] - b[2]) - aux[0]*mat[6] - aux[1]*mat[7];
511 x[2] = aux[2]/mat[8];
512 x[1] = (aux[1] - mat[5]*x[2])/mat[4];
513 x[0] = (aux[0] - mat[1]*x[1] - mat[2]*x[2])/mat[0];
535 assert(db_size <= DB_SIZE_MAX);
539 for (
int ii = 0; ii < db_size; ii++) {
540 aux[ii] = (c[ii] - b[ii]);
541 for (
int jj = 0; jj < ii; jj++) {
542 aux[ii] -= aux[jj]*mat[ii*db_size + jj];
547 for (
int ii = db_size - 1; ii >= 0; ii-=1) {
549 for (
int jj = db_size - 1; jj > ii; jj-=1) {
550 x[ii] -= x[jj]*mat[ii*db_size + jj];
552 x[ii] /= mat[ii*(db_size + 1)];
573 assert(db_size <= DB_SIZE_MAX);
576 for (
int ii = 0; ii < db_size; ii++) {
578 for (
int jj = 0; jj < ii; jj++)
579 x[ii] -= x[jj]*mat[ii*db_size + jj];
583 for (
int ii = db_size - 1; ii >= 0; ii--) {
584 for (
int jj = db_size - 1; jj > ii; jj--)
585 x[ii] -= x[jj]*mat[ii*db_size + jj];
586 x[ii] /= mat[ii*(db_size + 1)];
612 const char *solver_name,
640 bool block_nn_inverse);
double cs_dot(cs_lnum_t n, const cs_real_t *x, const cs_real_t *y)
Return the dot product of 2 vectors: x.y.
Definition: cs_blas.cpp:1495
void cs_dot_xx_yy_xy_xz_yz(cs_lnum_t n, const cs_real_t *restrict x, const cs_real_t *restrict y, const cs_real_t *restrict z, double *xx, double *yy, double *xy, double *xz, double *yz)
Return 5 dot products of 3 vectors: x.x, y.y, x.y, x.z, and y.z.
Definition: cs_blas.cpp:1684
double cs_dot_xx(cs_lnum_t n, const cs_real_t *x)
Return dot products of a vector with itself: x.x.
Definition: cs_blas.cpp:1516
void cs_dot_xy_yz(cs_lnum_t n, const cs_real_t *restrict x, const cs_real_t *restrict y, const cs_real_t *restrict z, double *xy, double *yz)
Return 2 dot products of 3 vectors: x.y, and y.z.
Definition: cs_blas.cpp:1623
void cs_dot_xx_xy_yz(cs_lnum_t n, const cs_real_t *restrict x, const cs_real_t *restrict y, const cs_real_t *restrict z, double *xx, double *xy, double *yz)
Return 3 dot products of 3 vectors: x.x, x.y, and y.z.
Definition: cs_blas.cpp:1652
void cs_dot_xx_xy(cs_lnum_t n, const cs_real_t *restrict x, const cs_real_t *restrict y, double *xx, double *xy)
Return 2 dot products of 2 vectors: x.x, and x.y.
Definition: cs_blas.cpp:1596
#define restrict
Definition: cs_defs.h:145
#define BEGIN_C_DECLS
Definition: cs_defs.h:542
double cs_real_t
Floating-point value.
Definition: cs_defs.h:342
#define END_C_DECLS
Definition: cs_defs.h:543
int cs_lnum_t
local mesh entity id
Definition: cs_defs.h:335
struct _cs_matrix_t cs_matrix_t
Definition: cs_matrix.h:110
cs_sles_convergence_state_t
Definition: cs_sles.h:56
struct _cs_sles_it_t cs_sles_it_t
Definition: cs_sles_it.h:86
struct _cs_sles_it_convergence_t cs_sles_it_convergence_t
Definition: cs_sles_it.h:90
cs_sles_it_type_t
Definition: cs_sles_it.h:55
cs_sles_pc_state_t() cs_sles_pc_apply_t(void *context, const cs_real_t *x_in, cs_real_t *x_out)
Function pointer for application of a preconditioner.
Definition: cs_sles_pc.h:145
struct _cs_sles_pc_t cs_sles_pc_t
Definition: cs_sles_pc.h:66
struct _cs_time_plot_t cs_time_plot_t
Definition: cs_time_plot.h:48
Definition: cs_timer.h:55