PaStiX Handbook  6.3.2
bcsc_d.h
Go to the documentation of this file.
1 /**
2  * @file bcsc_d.h
3  *
4  * @copyright 2004-2023 Bordeaux INP, CNRS (LaBRI UMR 5800), Inria,
5  * Univ. Bordeaux. All rights reserved.
6  *
7  * @version 6.3.2
8  * @author Mathieu Faverge
9  * @author Pierre Ramet
10  * @author Xavier Lacoste
11  * @author Vincent Bridonneau
12  * @date 2023-07-21
13  *
14  * @generated from /builds/solverstack/pastix/bcsc/bcsc_z.h, normal z -> d, Wed Dec 13 12:09:03 2023
15  *
16  **/
17 #ifndef _bcsc_d_h_
18 #define _bcsc_d_h_
19 
20 /**
21  * @addtogroup bcsc_internal
22  * @{
23  *
24  * @name PastixDouble initialization functions
25  * @{
26  */
27 void bcsc_dinit( const spmatrix_t *spm,
28  const pastix_order_t *ord,
29  const SolverMatrix *solvmtx,
30  int initAt,
31  pastix_bcsc_t *bcsc,
32  pastix_int_t valuesize );
33 
34 #if defined(PASTIX_WITH_MPI)
35 void bcsc_dstore_data( const spmatrix_t *spm,
36  const pastix_order_t *ord,
37  const pastix_int_t *col2cblk,
38  bcsc_handle_comm_t *bcsc_comm );
39 #endif
40 
41 /**
42  * @}
43  * @}
44  *
45  * @addtogroup bcsc
46  * @{
47  *
48  * @name PastixDouble vector(s) operations
49  * @{
50  */
51 void bvec_daxpy_seq( pastix_data_t *pastix_data,
52  pastix_int_t n,
53  double alpha,
54  const double *x,
55  double *y );
56 void bvec_daxpy_smp( pastix_data_t *pastix_data,
57  pastix_int_t n,
58  double alpha,
59  const double *x,
60  double *y );
61 
62 void bvec_dcopy_seq( pastix_data_t *pastix_data,
63  pastix_int_t n,
64  const double *x,
65  double *y );
66 void bvec_dcopy_smp( pastix_data_t *pastix_data,
67  pastix_int_t n,
68  const double *x,
69  double *y );
70 
71 #if defined(PRECISION_z) || defined(PRECISION_c)
72 double bvec_ddot_seq( pastix_data_t *pastix_data,
73  pastix_int_t n,
74  const double *x,
75  const double *y );
76 double bvec_ddot_smp( pastix_data_t *pastix_data,
77  pastix_int_t n,
78  const double *x,
79  const double *y );
80 #endif
81 
82 double bvec_ddot_seq( pastix_data_t *pastix_data,
83  pastix_int_t n,
84  const double *x,
85  const double *y );
86 double bvec_ddot_smp( pastix_data_t *pastix_data,
87  pastix_int_t n,
88  const double *x,
89  const double *y );
90 
91 void bvec_dgemv_seq( pastix_data_t *pastix_data,
92  pastix_int_t m,
93  pastix_int_t n,
94  double alpha,
95  const double *A,
96  pastix_int_t lda,
97  const double *x,
98  double beta,
99  double *y );
100 void bvec_dgemv_smp( pastix_data_t *pastix_data,
101  pastix_int_t m,
102  pastix_int_t n,
103  double alpha,
104  const double *A,
105  pastix_int_t lda,
106  const double *x,
107  double beta,
108  double *y );
109 
110 double bvec_dnrm2_seq( pastix_data_t *pastix_data,
111  pastix_int_t n,
112  const double *x );
113 double bvec_dnrm2_smp( pastix_data_t *pastix_data,
114  pastix_int_t n,
115  const double *x );
116 
117 void bvec_dscal_seq( pastix_data_t *pastix_data,
118  pastix_int_t n,
119  double alpha,
120  double *x );
121 void bvec_dscal_smp( pastix_data_t *pastix_data,
122  pastix_int_t n,
123  double alpha,
124  double *x );
125 
126 #if defined( PASTIX_WITH_MPI )
127 int bvec_dexchange_data_rep( pastix_data_t *pastix_data,
128  pastix_int_t nrhs,
129  double *b,
130  pastix_int_t ldb,
131  pastix_rhs_t Pb );
132 int bvec_dallocate_buf_dst( bvec_handle_comm_t *rhs_comm );
133 int bvec_dexchange_data_dst( pastix_data_t *pastix_data,
134  pastix_dir_t dir,
135  pastix_int_t nrhs,
136  double *b,
137  pastix_int_t ldb,
138  pastix_rhs_t Pb,
139  const pastix_int_t *glob2loc );
140 #endif
141 
142 int bvec_dlapmr( pastix_data_t *pastix_data,
143  pastix_dir_t dir,
144  pastix_int_t m,
145  pastix_int_t n,
146  double *A,
147  pastix_int_t lda,
148  pastix_rhs_t PA );
149 
150 /**
151  * @}
152  *
153  * @name PastixDouble matrix operations
154  * @{
155  */
156 double bcsc_dnorm( pastix_normtype_t ntype,
157  const pastix_bcsc_t *bcsc );
158 
159 void bcsc_dspsv( pastix_data_t *pastix_data,
160  double *b,
161  float *work );
162 
163 void bcsc_dspmv( const pastix_data_t *pastix_data,
164  pastix_trans_t trans,
165  double alpha,
166  const double *x,
167  double beta,
168  double *y );
169 
170 void bcsc_dspmv_seq( const pastix_data_t *pastix_data,
171  pastix_trans_t trans,
172  double alpha,
173  const double *x,
174  double beta,
175  double *y );
176 void bcsc_dspmv_smp( const pastix_data_t *pastix_data,
177  pastix_trans_t trans,
178  double alpha,
179  const double *x,
180  double beta,
181  double *y );
182 
183 /**
184  * @}
185  *
186  * @name PastixDouble MPI vector operations
187  * @{
188  */
189 const double *bvec_dgather_remote( const pastix_data_t *pastix_data,
190  const double *y );
191 void bvec_dnullify_remote( const pastix_data_t *pastix_data,
192  double *y );
193 void bvec_dallreduce( const pastix_data_t *pastix_data,
194  double *y );
195 /**
196  * @}
197  * @}
198  */
199 #endif /* _bcsc_d_h_ */
BEGIN_C_DECLS typedef int pastix_int_t
Definition: datatypes.h:51
void bcsc_dinit(const spmatrix_t *spm, const pastix_order_t *ord, const SolverMatrix *solvmtx, int initAt, pastix_bcsc_t *bcsc, pastix_int_t valuesize)
Initializes a centralize double block csc.
Definition: bcsc_dinit.c:1374
void bvec_dgemv_seq(pastix_data_t *pastix_data, pastix_int_t m, pastix_int_t n, double alpha, const double *A, pastix_int_t lda, const double *x, double beta, double *y)
Compute.
void bvec_daxpy_seq(pastix_data_t *pastix_data, pastix_int_t n, double alpha, const double *x, double *y)
Compute y <- alpha * x + y. (Sequential version)
double bvec_dnrm2_seq(pastix_data_t *pastix_data, pastix_int_t n, const double *x)
Compute the norm 2 of a vector. (Sequential version)
Definition: bvec_dcompute.c:71
void bcsc_dspmv_seq(const pastix_data_t *pastix_data, pastix_trans_t trans, double alpha, const double *x, double beta, double *y)
Compute the matrix-vector product y = alpha * A * x + beta * y (Sequential version)
Definition: bcsc_dspmv.c:246
double bvec_ddot_smp(pastix_data_t *pastix_data, pastix_int_t n, const double *x, const double *y)
Compute a regular scalar product x.y (Parallel version)
void bvec_dgemv_smp(pastix_data_t *pastix_data, pastix_int_t m, pastix_int_t n, double alpha, const double *A, pastix_int_t lda, const double *x, double beta, double *y)
Compute.
void bvec_daxpy_smp(pastix_data_t *pastix_data, pastix_int_t n, double alpha, const double *x, double *y)
Perform y = alpha * x + y (Parallel version)
void bcsc_dspmv_smp(const pastix_data_t *pastix_data, pastix_trans_t trans, double alpha, const double *x, double beta, double *y)
Perform y = alpha A x + beta y (Parallel version)
Definition: bcsc_dspmv.c:552
void bvec_dcopy_smp(pastix_data_t *pastix_data, pastix_int_t n, const double *x, double *y)
Copy a vector y = x (parallel version)
const double * bvec_dgather_remote(const pastix_data_t *pastix_data, const double *y)
Gather a distributed right hand side (bvec storage) on all nodes.
void bvec_dcopy_seq(pastix_data_t *pastix_data, pastix_int_t n, const double *x, double *y)
Copy a vector y = x (Sequential version)
void bcsc_dspsv(pastix_data_t *pastix_data, double *b, float *work)
Solve A x = b with A the sparse matrix.
void bvec_dscal_smp(pastix_data_t *pastix_data, pastix_int_t n, double alpha, double *x)
Scale a vector (Parallel version)
double bcsc_dnorm(pastix_normtype_t ntype, const pastix_bcsc_t *bcsc)
Compute the norm of an bcsc matrix.
Definition: bcsc_dnorm.c:252
int bvec_dlapmr(pastix_data_t *pastix_data, pastix_dir_t dir, pastix_int_t m, pastix_int_t n, double *A, pastix_int_t lda, pastix_rhs_t PA)
Apply a row permutation to a right hand side A (LAPACK xlatmr)
Definition: bvec_dlapmr.c:816
void bvec_dallreduce(const pastix_data_t *pastix_data, double *y)
Apply an all reduce of the vector on all nodes.
void bcsc_dspmv(const pastix_data_t *pastix_data, pastix_trans_t trans, double alpha, const double *x, double beta, double *y)
Compute the matrix-vector product y = alpha * op(A) * x + beta * y.
Definition: bcsc_dspmv.c:629
double bvec_ddot_seq(pastix_data_t *pastix_data, pastix_int_t n, const double *x, const double *y)
Compute the scalar product x.y. (Sequential version)
double bvec_dnrm2_smp(pastix_data_t *pastix_data, pastix_int_t n, const double *x)
Compute the norm 2 of a vector. (Parallel version)
void bvec_dnullify_remote(const pastix_data_t *pastix_data, double *y)
Set to 0 remote coefficients.
void bvec_dscal_seq(pastix_data_t *pastix_data, pastix_int_t n, double alpha, double *x)
Scale a vector by the scalar alpha. (Sequential version)
Structure to manage communications with distributed spm.
Definition: bcsc.h:92
Structure to manage communications with distributed rhs.
Definition: bvec.h:66
enum pastix_normtype_e pastix_normtype_t
Norms.
enum pastix_dir_e pastix_dir_t
Direction.
enum pastix_trans_e pastix_trans_t
Transpostion.
Order structure.
Definition: order.h:47
Main PaStiX data structure.
Definition: pastixdata.h:67
Main PaStiX RHS structure.
Definition: pastixdata.h:150
Solver column block structure.
Definition: solver.h:200