PaStiX Handbook  6.3.2
pastix.h
Go to the documentation of this file.
1 /**
2  *
3  * @file pastix.h
4  *
5  * PaStiX main header file.
6  *
7  * @copyright 2004-2023 Bordeaux INP, CNRS (LaBRI UMR 5800), Inria,
8  * Univ. Bordeaux. All rights reserved.
9  *
10  * @version 6.3.2
11  * @author David Goudin
12  * @author Francois Pellegrini
13  * @author Gregoire Pichon
14  * @author Mathieu Faverge
15  * @author Pascal Henon
16  * @author Pierre Ramet
17  * @author Xavier Lacoste
18  * @author Theophile Terraz
19  * @author Tony Delarue
20  * @date 2023-07-21
21  *
22  **/
23 #ifndef _pastix_h_
24 #define _pastix_h_
25 
26 #include "pastix/config.h"
27 #include <spm.h>
28 #include <stdio.h>
29 #include <stdlib.h>
30 #include <string.h>
31 #include <stdarg.h>
32 #include <assert.h>
33 #include <math.h>
34 
35 #ifndef DOXYGEN_SHOULD_SKIP_THIS
36 #if defined(PASTIX_WITH_MPI)
37 #include <mpi.h>
38 typedef MPI_Comm PASTIX_Comm;
39 #else
40 typedef uintptr_t PASTIX_Comm;
41 #ifndef MPI_COMM_WORLD
42 #define MPI_COMM_WORLD 0
43 #endif
44 #endif
45 #endif /* DOXYGEN_SHOULD_SKIP_THIS */
46 
47 #include "pastix/api.h"
48 #include "pastix/datatypes.h"
49 #include "pastix/order.h"
50 
51 BEGIN_C_DECLS
52 
53 /*
54  * Main function for compatibility with former versions
55  */
56 int pastix( pastix_data_t **pastix_data,
57  PASTIX_Comm pastix_comm,
58  pastix_int_t n,
59  pastix_int_t *colptr,
60  pastix_int_t *rowptr,
61  void *values,
62  pastix_int_t *perm,
63  pastix_int_t *invp,
64  void *B,
65  pastix_int_t nrhs,
66  pastix_int_t *iparm,
67  double *dparm );
68 
69 /*
70  * Solver initialization
71  */
72 void pastixInitParam( pastix_int_t *iparm,
73  double *dparm );
74 void pastixInit ( pastix_data_t **pastix_data,
75  PASTIX_Comm pastix_comm,
76  pastix_int_t *iparm,
77  double *dparm );
78 void pastixInitWithAffinity( pastix_data_t **pastix_data,
79  PASTIX_Comm pastix_comm,
80  pastix_int_t *iparm,
81  double *dparm,
82  const int *bindtab );
83 void pastixFinalize ( pastix_data_t **pastix_data );
84 
85 /*
86  * Main steps of the solver
87  */
88 int pastix_task_analyze( pastix_data_t *pastix_data,
89  const spmatrix_t *spm );
90 int pastix_task_numfact( pastix_data_t *pastix_data,
91  spmatrix_t *spm );
92 int pastix_task_solve ( pastix_data_t *pastix_data,
93  pastix_int_t m,
94  pastix_int_t nrhs,
95  void *B,
96  pastix_int_t ldb );
97 int pastix_task_refine ( pastix_data_t *pastix_data,
98  pastix_int_t n,
99  pastix_int_t nrhs,
100  void *B,
101  pastix_int_t ldb,
102  void *X,
103  pastix_int_t ldx );
105  pastix_int_t n,
106  pastix_int_t nrhs,
107  void *B,
108  pastix_int_t ldb,
109  void *X,
110  pastix_int_t ldx );
111 
112 /*
113  * Analyze subtasks
114  */
115 int pastix_subtask_order ( pastix_data_t *pastix_data,
116  const spmatrix_t *spm,
117  pastix_order_t *myorder );
118 int pastix_subtask_symbfact ( pastix_data_t *pastix_data );
119 int pastix_subtask_reordering( pastix_data_t *pastix_data );
120 int pastix_subtask_blend ( pastix_data_t *pastix_data );
121 
122 /*
123  * Numerical factorization subtasks
124  */
125 int pastix_subtask_spm2bcsc ( pastix_data_t *pastix_data,
126  spmatrix_t *spm );
127 int pastix_subtask_bcsc2ctab ( pastix_data_t *pastix_data );
128 int pastix_subtask_sopalin ( pastix_data_t *pastix_data );
129 
130 /*
131  * Numerical solve and refinement subtasks
132  */
133 int pastix_subtask_applyorder( pastix_data_t *pastix_data,
134  pastix_dir_t dir,
135  pastix_int_t m,
136  pastix_int_t n,
137  void *B,
138  pastix_int_t ldb,
139  pastix_rhs_t Bp );
140 int pastix_subtask_trsm( pastix_data_t *pastix_data,
141  pastix_side_t side,
142  pastix_uplo_t uplo,
143  pastix_trans_t trans,
144  pastix_diag_t diag,
145  pastix_rhs_t b );
146 int pastix_subtask_diag( pastix_data_t *pastix_data,
147  pastix_rhs_t b );
148 int pastix_subtask_solve( pastix_data_t *pastix_data,
149  pastix_rhs_t b );
150 int pastix_subtask_refine( pastix_data_t *pastix_data,
151  pastix_rhs_t b,
152  pastix_rhs_t x );
153 int pastix_subtask_solve_adv( pastix_data_t *pastix_data,
154  pastix_trans_t transA,
155  pastix_rhs_t b );
156 
157 /*
158  * Schur complement manipulation routines.
159  */
160 void pastixSetSchurUnknownList( pastix_data_t *pastix_data,
161  pastix_int_t n,
162  const pastix_int_t *list );
163 int pastixGetSchur ( const pastix_data_t *pastix_data,
164  void *S,
165  pastix_int_t lds );
166 
167 /*
168  * Right hand side routines.
169  */
170 int pastixRhsInit( pastix_rhs_t *rhs );
173  pastix_rhs_t sB );
175  pastix_rhs_t dB );
176 
177 int pastixRhsSchurGet( const pastix_data_t *pastix_data,
178  pastix_int_t m,
179  pastix_int_t n,
180  pastix_rhs_t rhsB,
181  void *B,
182  pastix_int_t ldb );
183 int pastixRhsSchurSet( const pastix_data_t *pastix_data,
184  pastix_int_t m,
185  pastix_int_t n,
186  void *B,
187  pastix_int_t ldb,
188  pastix_rhs_t rhsB );
189 
190 /*
191  * DoF subroutine to expand the problem.
192  */
193 void pastixExpand( const pastix_data_t *pastix_data,
194  spmatrix_t *spm );
195 
196 /*
197  * Function to provide access to the diagonal
198  */
199 int pastixGetDiag( const pastix_data_t *pastix_data,
200  void *x,
201  pastix_int_t incx );
202 
203 /*
204  * Function to provide a common way to read binary options in examples/testings
205  */
206 void pastixGetOptions( int argc,
207  char **argv,
208  pastix_int_t *iparm,
209  double *dparm,
210  int *check,
211  int *scatter,
212  spm_driver_t *driver,
213  char **filename );
214 
215 /*
216  * Function to provide a common way to output
217  * the iparm/dparm parameters in a CSV file.
218  */
219 void pastixDumpParam ( const pastix_data_t *pastix_data );
220 int pastixCheckParam( const pastix_int_t *iparm,
221  const double *dparm );
222 
223 END_C_DECLS
224 
225 #endif /* _pastix_h_ */
BEGIN_C_DECLS typedef int pastix_int_t
Definition: datatypes.h:51
int pastixRhsInit(pastix_rhs_t *rhs)
Initialize an RHS data structure.
Definition: pastix_rhs.c:41
int pastixRhsFinalize(pastix_rhs_t rhs)
Cleanup an RHS data structure.
Definition: pastix_rhs.c:86
int pastix_subtask_symbfact(pastix_data_t *pastix_data)
Computes the symbolic factorization step.
int pastix_subtask_order(pastix_data_t *pastix_data, const spmatrix_t *spm, pastix_order_t *myorder)
Computes the ordering of the given graph in parameters.
int pastix_subtask_blend(pastix_data_t *pastix_data)
Compute the proportional mapping and the final solver structure.
int pastix_subtask_reordering(pastix_data_t *pastix_data)
Apply the reordering step to compact off-diagonal blocks.
void pastixDumpParam(const pastix_data_t *pastix_data)
Dump the iparm and dparm parameters in a CSV file.
Definition: api.c:1031
void pastixFinalize(pastix_data_t **pastix_data)
Finalize the solver instance.
Definition: api.c:919
enum pastix_diag_e pastix_diag_t
Diagonal.
void pastixInitParam(pastix_int_t *iparm, double *dparm)
Initialize the iparm and dparm arrays to their default values.
Definition: api.c:411
enum pastix_dir_e pastix_dir_t
Direction.
void pastixInit(pastix_data_t **pastix_data, PASTIX_Comm pastix_comm, pastix_int_t *iparm, double *dparm)
Initialize the solver instance.
Definition: api.c:896
enum pastix_uplo_e pastix_uplo_t
Upper/Lower part.
int pastixCheckParam(const pastix_int_t *iparm, const double *dparm)
Check the values of iparm and dparm arrays.
Definition: api.c:1079
BEGIN_C_DECLS int pastix(pastix_data_t **pastix_data, PASTIX_Comm pastix_comm, pastix_int_t n, pastix_int_t *colptr, pastix_int_t *rowptr, void *values, pastix_int_t *perm, pastix_int_t *invp, void *B, pastix_int_t nrhs, pastix_int_t *iparm, double *dparm)
Main function for compatibility with former releases.
Definition: pastix.c:103
void pastixInitWithAffinity(pastix_data_t **pastix_data, PASTIX_Comm pastix_comm, pastix_int_t *iparm, double *dparm, const int *bindtab)
Initialize the solver instance with a bintab array to specify the thread binding.
Definition: api.c:700
enum pastix_side_e pastix_side_t
Side of the operation.
enum pastix_trans_e pastix_trans_t
Transpostion.
void pastixGetOptions(int argc, char **argv, pastix_int_t *iparm, double *dparm, int *check, int *scatter, spm_driver_t *driver, char **filename)
PaStiX helper function to read command line options in examples.
Definition: get_options.c:149
int pastix_subtask_sopalin(pastix_data_t *pastix_data)
Factorize the given problem using Cholesky or LU decomposition.
int pastix_subtask_spm2bcsc(pastix_data_t *pastix_data, spmatrix_t *spm)
Fill the internal block CSC structure.
int pastix_subtask_bcsc2ctab(pastix_data_t *pastix_data)
Fill the internal solver matrix structure.
Order structure.
Definition: order.h:47
int pastix_subtask_refine(pastix_data_t *pastix_data, pastix_rhs_t b, pastix_rhs_t x)
Perform the iterative refinement without apply the permutations.
void pastixSetSchurUnknownList(pastix_data_t *pastix_data, pastix_int_t n, const pastix_int_t *list)
Set the list of unknowns that belongs to the schur complement.
Definition: schur.c:50
int pastixGetDiag(const pastix_data_t *pastix_data, void *x, pastix_int_t incx)
Return the diagonal of the matrix.
Definition: diag.c:54
int pastixGetSchur(const pastix_data_t *pastix_data, void *S, pastix_int_t lds)
Return the Schur complement.
Definition: schur.c:90
int pastix_subtask_solve_adv(pastix_data_t *pastix_data, pastix_trans_t transA, pastix_rhs_t b)
Solve the given problem without applying the permutation.
int pastixRhsDoubletoSingle(const pastix_rhs_t dB, pastix_rhs_t sB)
Reduces the precision of an RHS.
Definition: pastix_rhs.c:146
int pastixRhsSchurSet(const pastix_data_t *pastix_data, pastix_int_t m, pastix_int_t n, void *B, pastix_int_t ldb, pastix_rhs_t rhsB)
Set the vector in an RHS data structure.
Definition: schur.c:284
int pastix_subtask_applyorder(pastix_data_t *pastix_data, pastix_dir_t dir, pastix_int_t m, pastix_int_t n, void *B, pastix_int_t ldb, pastix_rhs_t Bp)
Apply a permutation on the right-and-side vector before the solve step.
int pastix_subtask_diag(pastix_data_t *pastix_data, pastix_rhs_t b)
Apply a diagonal operation on the right-and-side vectors.
int pastix_subtask_solve(pastix_data_t *pastix_data, pastix_rhs_t b)
Solve the given problem without applying the permutation.
int pastix_subtask_trsm(pastix_data_t *pastix_data, pastix_side_t side, pastix_uplo_t uplo, pastix_trans_t trans, pastix_diag_t diag, pastix_rhs_t b)
Apply a triangular solve on the right-and-side vectors.
int pastixRhsSchurGet(const pastix_data_t *pastix_data, pastix_int_t m, pastix_int_t n, pastix_rhs_t rhsB, void *B, pastix_int_t ldb)
Get the vector in an RHS data structure.
Definition: schur.c:182
int pastixRhsSingleToDouble(const pastix_rhs_t sB, pastix_rhs_t dB)
Increases the precision of an RHS.
Definition: pastix_rhs.c:226
int pastix_task_refine(pastix_data_t *pastix_data, pastix_int_t n, pastix_int_t nrhs, void *B, pastix_int_t ldb, void *X, pastix_int_t ldx)
Perform iterative refinement.
int pastix_task_analyze(pastix_data_t *pastix_data, const spmatrix_t *spm)
Perform all the preprocessing steps: ordering, symbolic factorization, reordering,...
int pastix_task_solve(pastix_data_t *pastix_data, pastix_int_t m, pastix_int_t nrhs, void *B, pastix_int_t ldb)
Solve the given problem.
int pastix_task_solve_and_refine(pastix_data_t *pastix_data, pastix_int_t n, pastix_int_t nrhs, void *B, pastix_int_t ldb, void *X, pastix_int_t ldx)
Performs solve and iterative refinement without unnecessary permutations.
int pastix_task_numfact(pastix_data_t *pastix_data, spmatrix_t *spm)
Perform all the numerical factorization steps: fill the internal block CSC and the solver matrix stru...
Main PaStiX data structure.
Definition: pastixdata.h:67
Main PaStiX RHS structure.
Definition: pastixdata.h:150
void pastixExpand(const pastix_data_t *pastix_data, spmatrix_t *spm)
Expand an spm structure and the already computed data structure associated if any.
Definition: pastix.c:506