Fast Auxiliary Space Preconditioning 2.7.7 Aug/28/2022
Loading...
Searching...
No Matches
SolBLC.c
Go to the documentation of this file.
1
17#include <math.h>
18#include <time.h>
19
20#include "fasp.h"
21#include "fasp_block.h"
22#include "fasp_functs.h"
23
24/*---------------------------------*/
25/*-- Declare Private Functions --*/
26/*---------------------------------*/
27
28#include "KryUtil.inl"
29
30/*---------------------------------*/
31/*-- Public Functions --*/
32/*---------------------------------*/
33
55 ITS_param* itparam)
56{
57 const SHORT prtlvl = itparam->print_level;
58 const SHORT itsolver_type = itparam->itsolver_type;
59 const SHORT stop_type = itparam->stop_type;
60 const SHORT restart = itparam->restart;
61 const INT MaxIt = itparam->maxit;
62 const REAL tol = itparam->tol;
63 const REAL abstol = itparam->abstol;
64
65 REAL solve_start, solve_end;
67
68#if DEBUG_MODE > 0
69 printf("### DEBUG: [-Begin-] %s ...\n", __FUNCTION__);
70 printf("### DEBUG: rhs/sol size: %d %d\n", b->row, x->row);
71#endif
72
73 fasp_gettime(&solve_start);
74
75 /* Safe-guard checks on parameters */
76 ITS_CHECK(MaxIt, tol);
77
78 switch (itsolver_type) {
79
80 case SOLVER_BiCGstab:
81 iter = fasp_solver_dblc_pbcgs(A, b, x, pc, tol, abstol, MaxIt, stop_type,
82 prtlvl);
83 break;
84
85 case SOLVER_MinRes:
86 iter = fasp_solver_dblc_pminres(A, b, x, pc, tol, abstol, MaxIt, stop_type,
87 prtlvl);
88 break;
89
90 case SOLVER_GMRES:
91 iter = fasp_solver_dblc_pgmres(A, b, x, pc, tol, abstol, MaxIt, restart,
92 stop_type, prtlvl);
93 break;
94
95 case SOLVER_VGMRES:
96 iter = fasp_solver_dblc_pvgmres(A, b, x, pc, tol, abstol, MaxIt, restart,
97 stop_type, prtlvl);
98 break;
99
100 case SOLVER_VFGMRES:
101 iter = fasp_solver_dblc_pvfgmres(A, b, x, pc, tol, abstol, MaxIt, restart,
102 stop_type, prtlvl);
103 break;
104
105 default:
106 printf("### ERROR: Unknown iterative solver type %d! [%s]\n", itsolver_type,
107 __FUNCTION__);
108 return ERROR_SOLVER_TYPE;
109 }
110
111 if ((prtlvl >= PRINT_MIN) && (iter >= 0)) {
112 fasp_gettime(&solve_end);
113 fasp_cputime("Iterative method", solve_end - solve_start);
114 }
115
116#if DEBUG_MODE > 0
117 printf("### DEBUG: [--End--] %s ...\n", __FUNCTION__);
118#endif
119
120 return iter;
121}
122
140{
141 const SHORT prtlvl = itparam->print_level;
142
143 INT status = FASP_SUCCESS;
144 REAL solve_start, solve_end;
145
146#if DEBUG_MODE > 0
147 printf("### DEBUG: [-Begin-] %s ...\n", __FUNCTION__);
148#endif
149
150 // solver part
151 fasp_gettime(&solve_start);
152
153 status = fasp_solver_dblc_itsolver(A, b, x, NULL, itparam);
154
155 fasp_gettime(&solve_end);
156
157 if (prtlvl >= PRINT_MIN)
158 fasp_cputime("Krylov method totally", solve_end - solve_start);
159
160#if DEBUG_MODE > 0
161 printf("### DEBUG: [--End--] %s ...\n", __FUNCTION__);
162#endif
163
164 return status;
165}
166
189 ITS_param* itparam, AMG_param* amgparam,
190 dCSRmat* A_diag)
191{
192 const SHORT prtlvl = itparam->print_level;
193 const SHORT precond_type = itparam->precond_type;
194
195 INT status = FASP_SUCCESS;
196 REAL setup_start, setup_end;
197 REAL solve_start, solve_end;
198
199 const SHORT max_levels = amgparam->max_levels;
200 INT m, n, nnz, i;
201
202 AMG_data** mgl = NULL;
203
204#if WITH_UMFPACK
205 void** LU_diag = (void**)fasp_mem_calloc(3, sizeof(void*));
206#endif
207
208#if DEBUG_MODE > 0
209 printf("### DEBUG: [-Begin-] %s ...\n", __FUNCTION__);
210#endif
211
212 /* setup preconditioner */
213 fasp_gettime(&solve_start);
214 fasp_gettime(&setup_start);
215
216 /* diagonal blocks are solved exactly */
217 if (precond_type > 20 && precond_type < 30) {
218#if WITH_UMFPACK
219 // Need to sort the diagonal blocks for UMFPACK format
220 dCSRmat A_tran;
221
222 for (i = 0; i < 3; i++) {
223
224 A_tran = fasp_dcsr_create(A_diag[i].row, A_diag[i].col, A_diag[i].nnz);
225 fasp_dcsr_transz(&A_diag[i], NULL, &A_tran);
226 fasp_dcsr_cp(&A_tran, &A_diag[i]);
227
228 printf("Factorization for %d-th diagnol: \n", i);
229 LU_diag[i] = fasp_umfpack_factorize(&A_diag[i], prtlvl);
230 }
231
232 fasp_dcsr_free(&A_tran);
233#endif
234 }
235
236 /* diagonal blocks are solved by AMG */
237 else if (precond_type > 30 && precond_type < 40) {
238
239 mgl = (AMG_data**)fasp_mem_calloc(3, sizeof(AMG_data*));
240
241 for (i = 0; i < 3; i++) {
242
243 mgl[i] = fasp_amg_data_create(max_levels);
244 m = A_diag[i].row;
245 n = A_diag[i].col;
246 nnz = A_diag[i].nnz;
247 mgl[i][0].A = fasp_dcsr_create(m, n, nnz);
248 fasp_dcsr_cp(&A_diag[i], &mgl[i][0].A);
249 mgl[i][0].b = fasp_dvec_create(n);
250 mgl[i][0].x = fasp_dvec_create(n);
251
252 switch (amgparam->AMG_type) {
253 case SA_AMG: // Smoothed Aggregation AMG
254 status = fasp_amg_setup_sa(mgl[i], amgparam);
255 break;
256 case UA_AMG: // Unsmoothed Aggregation AMG
257 status = fasp_amg_setup_ua(mgl[i], amgparam);
258 break;
259 default: // Classical AMG
260 status = fasp_amg_setup_rs(mgl[i], amgparam);
261 break;
262 }
263
264 fasp_chkerr(status, __FUNCTION__);
265 }
266
267 }
268
269 else {
270 fasp_chkerr(ERROR_SOLVER_PRECTYPE, __FUNCTION__);
271 }
272
273 precond_data_blc precdata;
274 precdata.Ablc = A;
275 precdata.A_diag = A_diag;
276 precdata.r = fasp_dvec_create(b->row);
277
278 /* diagonal blocks are solved exactly */
279 if (precond_type > 20 && precond_type < 30) {
280#if WITH_UMFPACK
281 precdata.LU_diag = LU_diag;
282#endif
283 }
284 /* diagonal blocks are solved by AMG */
285 else if (precond_type > 30 && precond_type < 40) {
286 precdata.amgparam = amgparam;
287 precdata.mgl = mgl;
288 } else {
289 fasp_chkerr(ERROR_SOLVER_PRECTYPE, __FUNCTION__);
290 }
291
292 precond prec;
293 prec.data = &precdata;
294
295 switch (precond_type) {
296 case 21:
298 break;
299
300 case 22:
302 break;
303
304 case 23:
306 break;
307
308 case 24:
310 break;
311
312 case 31:
314 break;
315
316 case 32:
318 break;
319
320 case 33:
322 break;
323
324 case 34:
326 break;
327
328 default:
329 fasp_chkerr(ERROR_SOLVER_PRECTYPE, __FUNCTION__);
330 break;
331 }
332
333 if (prtlvl >= PRINT_MIN) {
334 fasp_gettime(&setup_end);
335 fasp_cputime("Setup totally", setup_end - setup_start);
336 }
337
338 // solve part
339 status = fasp_solver_dblc_itsolver(A, b, x, &prec, itparam);
340
341 fasp_gettime(&solve_end);
342
343 if (prtlvl >= PRINT_MIN)
344 fasp_cputime("Krylov method totally", solve_end - solve_start);
345
346 // clean up
347 /* diagonal blocks are solved exactly */
348 if (precond_type > 20 && precond_type < 30) {
349#if WITH_UMFPACK
350 for (i = 0; i < 3; i++) fasp_umfpack_free_numeric(LU_diag[i]);
351#endif
352 }
353 /* diagonal blocks are solved by AMG */
354 else if (precond_type > 30 && precond_type < 40) {
355 for (i = 0; i < 3; i++) fasp_amg_data_free(mgl[i], amgparam);
356 fasp_mem_free(mgl);
357 mgl = NULL;
358 } else {
359 fasp_chkerr(ERROR_SOLVER_PRECTYPE, __FUNCTION__);
360 }
361
362#if DEBUG_MODE > 0
363 printf("### DEBUG: [--End--] %s ...\n", __FUNCTION__);
364#endif
365
366 return status;
367}
368
391 ITS_param* itparam, AMG_param* amgparam,
392 dCSRmat* A_diag)
393{
394 const SHORT prtlvl = itparam->print_level;
395 const SHORT precond_type = itparam->precond_type;
396
397 INT status = FASP_SUCCESS;
398 REAL setup_start, setup_end;
399 REAL solve_start, solve_end;
400
401#if DEBUG_MODE > 0
402 printf("### DEBUG: [-Begin-] %s ...\n", __FUNCTION__);
403#endif
404
405 /* setup preconditioner */
406 fasp_gettime(&solve_start);
407 fasp_gettime(&setup_start);
408
409#if WITH_UMFPACK
410 void** LU_diag = (void**)fasp_mem_calloc(4, sizeof(void*));
411 INT i;
412#endif
413
414 /* diagonal blocks are solved exactly */
415 if (precond_type > 20 && precond_type < 30) {
416
417#if WITH_UMFPACK
418 // Need to sort the matrices local_A for UMFPACK format
419 dCSRmat A_tran;
420
421 for (i = 0; i < 4; i++) {
422
423 A_tran = fasp_dcsr_create(A_diag[i].row, A_diag[i].col, A_diag[i].nnz);
424 fasp_dcsr_transz(&A_diag[i], NULL, &A_tran);
425 fasp_dcsr_cp(&A_tran, &A_diag[i]);
426
427 printf("Factorization for %d-th diagnol: \n", i);
428 LU_diag[i] = fasp_umfpack_factorize(&A_diag[i], prtlvl);
429 }
430
431 fasp_dcsr_free(&A_tran);
432#endif
433
434 } else {
435 fasp_chkerr(ERROR_SOLVER_PRECTYPE, __FUNCTION__);
436 }
437
438 precond_data_blc precdata;
439
440 precdata.Ablc = A;
441 precdata.A_diag = A_diag;
442#if WITH_UMFPACK
443 precdata.LU_diag = LU_diag;
444#endif
445 precdata.r = fasp_dvec_create(b->row);
446
447 precond prec;
448 prec.data = &precdata;
449
450 switch (precond_type) {
451 case 21:
453 break;
454
455 case 22:
457 break;
458 }
459
460 if (prtlvl >= PRINT_MIN) {
461 fasp_gettime(&setup_end);
462 fasp_cputime("Setup totally", setup_end - setup_start);
463 }
464
465 // solver part
466 status = fasp_solver_dblc_itsolver(A, b, x, &prec, itparam);
467
468 fasp_gettime(&solve_end);
469
470 if (prtlvl >= PRINT_MIN)
471 fasp_cputime("Krylov method totally", solve_end - solve_start);
472
473 // clean
474#if WITH_UMFPACK
475 for (i = 0; i < 4; i++) fasp_umfpack_free_numeric(LU_diag[i]);
476#endif
477
478#if DEBUG_MODE > 0
479 printf("### DEBUG: [--End--] %s ...\n", __FUNCTION__);
480#endif
481
482 return status;
483}
484
508 ITS_param* itparam, INT NumLayers, dBLCmat* Ai,
509 dCSRmat* local_A, ivector* local_index)
510{
511 const SHORT prtlvl = itparam->print_level;
512
513 INT status = FASP_SUCCESS;
514 REAL setup_start, setup_end;
515 REAL solve_start, solve_end;
516
517 void** local_LU = NULL;
518
519#if DEBUG_MODE > 0
520 printf("### DEBUG: [-Begin-] %s ...\n", __FUNCTION__);
521#endif
522
523 /* setup preconditioner */
524 fasp_gettime(&solve_start);
525 fasp_gettime(&setup_start);
526
527#if WITH_UMFPACK
528 // Need to sort the matrices local_A for UMFPACK format
529 INT l;
530 dCSRmat A_tran;
531 local_LU = (void**)fasp_mem_calloc(NumLayers, sizeof(void*));
532
533 for (l = 0; l < NumLayers; l++) {
534
535 A_tran = fasp_dcsr_create(local_A[l].row, local_A[l].col, local_A[l].nnz);
536 fasp_dcsr_transz(&local_A[l], NULL, &A_tran);
537 fasp_dcsr_cp(&A_tran, &local_A[l]);
538
539 printf("Factorization for layer %d: \n", l);
540 local_LU[l] = fasp_umfpack_factorize(&local_A[l], prtlvl);
541 }
542
543 fasp_dcsr_free(&A_tran);
544#endif
545
546 precond_data_sweeping precdata;
547 precdata.NumLayers = NumLayers;
548 precdata.A = A;
549 precdata.Ai = Ai;
550 precdata.local_A = local_A;
551 precdata.local_LU = local_LU;
552 precdata.local_index = local_index;
553 precdata.r = fasp_dvec_create(b->row);
554 precdata.w = (REAL*)fasp_mem_calloc(10 * b->row, sizeof(REAL));
555
556 precond prec;
557 prec.data = &precdata;
559
560 if (prtlvl >= PRINT_MIN) {
561 fasp_gettime(&setup_end);
562 fasp_cputime("Setup totally", setup_end - setup_start);
563 }
564
565 /* solver part */
566 status = fasp_solver_dblc_itsolver(A, b, x, &prec, itparam);
567
568 fasp_gettime(&solve_end);
569
570 if (prtlvl >= PRINT_MIN)
571 fasp_cputime("Krylov method totally", solve_end - solve_start);
572
573 // clean
574#if WITH_UMFPACK
575 for (l = 0; l < NumLayers; l++) fasp_umfpack_free_numeric(local_LU[l]);
576#endif
577
578#if DEBUG_MODE > 0
579 printf("### DEBUG: [--End--] %s ...\n", __FUNCTION__);
580#endif
581
582 return status;
583}
584
585/*---------------------------------*/
586/*-- End of File --*/
587/*---------------------------------*/
void fasp_mem_free(void *mem)
Free up previous allocated memory body and set pointer to NULL.
Definition: AuxMemory.c:152
void * fasp_mem_calloc(const unsigned int size, const unsigned int type)
Allocate, initiate, and check memory.
Definition: AuxMemory.c:65
void fasp_cputime(const char *message, const REAL cputime)
Print CPU walltime.
Definition: AuxMessage.c:179
void fasp_chkerr(const SHORT status, const char *fctname)
Check error status and print out error messages before quit.
Definition: AuxMessage.c:213
void fasp_gettime(REAL *time)
Get system time.
Definition: AuxTiming.c:37
dvector fasp_dvec_create(const INT m)
Create dvector data space of REAL type.
Definition: AuxVector.c:62
dCSRmat fasp_dcsr_create(const INT m, const INT n, const INT nnz)
Create CSR sparse matrix data memory space.
Definition: BlaSparseCSR.c:47
void fasp_dcsr_free(dCSRmat *A)
Free CSR sparse matrix data memory space.
Definition: BlaSparseCSR.c:184
void fasp_dcsr_transz(dCSRmat *A, INT *p, dCSRmat *AT)
Generalized transpose of A: (n x m) matrix given in dCSRmat format.
void fasp_dcsr_cp(const dCSRmat *A, dCSRmat *B)
copy a dCSRmat to a new one B=A
Definition: BlaSparseCSR.c:851
INT fasp_solver_dblc_pbcgs(dBLCmat *A, dvector *b, dvector *u, precond *pc, const REAL tol, const REAL abstol, const INT MaxIt, const SHORT StopType, const SHORT PrtLvl)
Preconditioned BiCGstab method for solving Au=b for BLC matrix.
Definition: KryPbcgs.c:705
INT fasp_solver_dblc_pgmres(dBLCmat *A, dvector *b, dvector *x, precond *pc, const REAL tol, const REAL abstol, const INT MaxIt, const SHORT restart, const SHORT StopType, const SHORT PrtLvl)
Preconditioned GMRES method for solving Au=b.
Definition: KryPgmres.c:687
INT fasp_solver_dblc_pminres(dBLCmat *A, dvector *b, dvector *u, precond *pc, const REAL tol, const REAL abstol, const INT MaxIt, const SHORT StopType, const SHORT PrtLvl)
A preconditioned minimal residual (Minres) method for solving Au=b.
Definition: KryPminres.c:470
INT fasp_solver_dblc_pvfgmres(dBLCmat *A, dvector *b, dvector *x, precond *pc, const REAL tol, const REAL abstol, const INT MaxIt, const SHORT restart, const SHORT StopType, const SHORT PrtLvl)
Solve "Ax=b" using PFGMRES (right preconditioned) iterative method in which the restart parameter can...
Definition: KryPvfgmres.c:707
INT fasp_solver_dblc_pvgmres(dBLCmat *A, dvector *b, dvector *x, precond *pc, const REAL tol, const REAL abstol, const INT MaxIt, const SHORT restart, const SHORT StopType, const SHORT PrtLvl)
Right preconditioned GMRES method in which the restart parameter can be adaptively modified during it...
Definition: KryPvgmres.c:765
SHORT fasp_amg_setup_rs(AMG_data *mgl, AMG_param *param)
Setup phase of Ruge and Stuben's classic AMG.
Definition: PreAMGSetupRS.c:52
SHORT fasp_amg_setup_sa(AMG_data *mgl, AMG_param *param)
Set up phase of smoothed aggregation AMG.
Definition: PreAMGSetupSA.c:63
SHORT fasp_amg_setup_ua(AMG_data *mgl, AMG_param *param)
Set up phase of unsmoothed aggregation AMG.
Definition: PreAMGSetupUA.c:55
void fasp_precond_dblc_upper_3(REAL *r, REAL *z, void *data)
block upper triangular preconditioning (3x3 blocks)
Definition: PreBLC.c:557
void fasp_precond_dblc_SGS_3(REAL *r, REAL *z, void *data)
Block symmetric GS preconditioning (3x3 blocks)
Definition: PreBLC.c:725
void fasp_precond_dblc_diag_4(REAL *r, REAL *z, void *data)
Block diagonal preconditioning (4x4 blocks)
Definition: PreBLC.c:191
void fasp_precond_dblc_lower_3(REAL *r, REAL *z, void *data)
block lower triangular preconditioning (3x3 blocks)
Definition: PreBLC.c:291
void fasp_precond_dblc_SGS_3_amg(REAL *r, REAL *z, void *data)
Block symmetric GS preconditioning (3x3 blocks)
Definition: PreBLC.c:838
void fasp_precond_dblc_diag_3_amg(REAL *r, REAL *z, void *data)
Block diagonal preconditioning (3x3 blocks)
Definition: PreBLC.c:126
void fasp_precond_dblc_lower_4(REAL *r, REAL *z, void *data)
block lower triangular preconditioning (4x4 blocks)
Definition: PreBLC.c:453
void fasp_precond_dblc_upper_3_amg(REAL *r, REAL *z, void *data)
block upper triangular preconditioning (3x3 blocks)
Definition: PreBLC.c:645
void fasp_precond_dblc_lower_3_amg(REAL *r, REAL *z, void *data)
block lower triangular preconditioning (3x3 blocks)
Definition: PreBLC.c:379
void fasp_precond_dblc_diag_3(REAL *r, REAL *z, void *data)
Block diagonal preconditioner (3x3 blocks)
Definition: PreBLC.c:38
void fasp_precond_dblc_sweeping(REAL *r, REAL *z, void *data)
Sweeping preconditioner for Maxwell equations.
Definition: PreBLC.c:939
AMG_data * fasp_amg_data_create(SHORT max_levels)
Create and initialize AMG_data for classical and SA AMG.
Definition: PreDataInit.c:64
void fasp_amg_data_free(AMG_data *mgl, AMG_param *param)
Free AMG_data data memeory space.
Definition: PreDataInit.c:101
INT fasp_solver_dblc_krylov(dBLCmat *A, dvector *b, dvector *x, ITS_param *itparam)
Solve Ax = b by standard Krylov methods.
Definition: SolBLC.c:139
INT fasp_solver_dblc_itsolver(dBLCmat *A, dvector *b, dvector *x, precond *pc, ITS_param *itparam)
Solve Ax = b by standard Krylov methods.
Definition: SolBLC.c:54
INT fasp_solver_dblc_krylov_block4(dBLCmat *A, dvector *b, dvector *x, ITS_param *itparam, AMG_param *amgparam, dCSRmat *A_diag)
Solve Ax = b by standard Krylov methods.
Definition: SolBLC.c:390
INT fasp_solver_dblc_krylov_block3(dBLCmat *A, dvector *b, dvector *x, ITS_param *itparam, AMG_param *amgparam, dCSRmat *A_diag)
Solve Ax = b by standard Krylov methods.
Definition: SolBLC.c:188
INT fasp_solver_dblc_krylov_sweeping(dBLCmat *A, dvector *b, dvector *x, ITS_param *itparam, INT NumLayers, dBLCmat *Ai, dCSRmat *local_A, ivector *local_index)
Solve Ax = b by standard Krylov methods.
Definition: SolBLC.c:507
Main header file for the FASP project.
#define REAL
Definition: fasp.h:75
#define SHORT
FASP integer and floating point numbers.
Definition: fasp.h:71
#define INT
Definition: fasp.h:72
Header file for FASP block matrices.
#define SOLVER_GMRES
Definition: fasp_const.h:106
#define FASP_SUCCESS
Definition of return status and error messages.
Definition: fasp_const.h:19
#define ERROR_SOLVER_TYPE
Definition: fasp_const.h:41
#define SOLVER_VGMRES
Definition: fasp_const.h:107
#define ERROR_SOLVER_PRECTYPE
Definition: fasp_const.h:42
#define SOLVER_BiCGstab
Definition: fasp_const.h:104
#define SA_AMG
Definition: fasp_const.h:164
#define SOLVER_MinRes
Definition: fasp_const.h:105
#define SOLVER_VFGMRES
Definition: fasp_const.h:108
#define PRINT_MIN
Definition: fasp_const.h:74
#define UA_AMG
Definition: fasp_const.h:165
Data for AMG methods.
Definition: fasp.h:804
dCSRmat A
pointer to the matrix at level level_num
Definition: fasp.h:817
dvector b
pointer to the right-hand side at level level_num
Definition: fasp.h:826
dvector x
pointer to the iterative solution at level level_num
Definition: fasp.h:829
Parameters for AMG methods.
Definition: fasp.h:455
SHORT AMG_type
type of AMG method
Definition: fasp.h:458
SHORT max_levels
max number of levels of AMG
Definition: fasp.h:470
Parameters for iterative solvers.
Definition: fasp.h:386
SHORT itsolver_type
Definition: fasp.h:389
SHORT print_level
Definition: fasp.h:388
SHORT precond_type
Definition: fasp.h:391
REAL tol
Definition: fasp.h:395
INT restart
Definition: fasp.h:393
SHORT stop_type
Definition: fasp.h:392
REAL abstol
Definition: fasp.h:396
INT maxit
Definition: fasp.h:394
Block REAL CSR matrix format.
Definition: fasp_block.h:74
Sparse matrix of REAL type in CSR format.
Definition: fasp.h:151
INT col
column of matrix A, n
Definition: fasp.h:157
INT row
row number of matrix A, m
Definition: fasp.h:154
INT nnz
number of nonzero entries
Definition: fasp.h:160
Vector with n entries of REAL type.
Definition: fasp.h:354
INT row
number of rows
Definition: fasp.h:357
Vector with n entries of INT type.
Definition: fasp.h:368
Data for block preconditioners in dBLCmat format.
Definition: fasp_block.h:363
dBLCmat * Ablc
Definition: fasp_block.h:368
dCSRmat * A_diag
Definition: fasp_block.h:370
AMG_data ** mgl
Definition: fasp_block.h:382
AMG_param * amgparam
Definition: fasp_block.h:384
Data for sweeping preconditioner.
Definition: fasp_block.h:398
Preconditioner data and action.
Definition: fasp.h:1095
void * data
data for preconditioner, void pointer
Definition: fasp.h:1098
void(* fct)(REAL *, REAL *, void *)
action for preconditioner, void function pointer
Definition: fasp.h:1101