Fast Auxiliary Space Preconditioning 2.7.7 Aug/28/2022
Loading...
Searching...
No Matches
SolWrapper.c
Go to the documentation of this file.
1
19#include "fasp.h"
20#include "fasp_block.h"
21#include "fasp_functs.h"
22
23/*---------------------------------*/
24/*-- Public Functions --*/
25/*---------------------------------*/
26
45#if WITH_PARDISO
46void fasp_fwrapper_dcsr_pardiso_(
47 INT* n, INT* nnz, INT* ia, INT* ja, REAL* a, REAL* b, REAL* u, INT* ptrlvl)
48{
49 dCSRmat mat; // coefficient matrix
50 dvector rhs, sol; // right-hand-side, solution
51
52 // set up coefficient matrix
53 mat.row = *n;
54 mat.col = *n;
55 mat.nnz = *nnz;
56 mat.IA = ia;
57 mat.JA = ja;
58 mat.val = a;
59
60 rhs.row = *n;
61 rhs.val = b;
62 sol.row = *n;
63 sol.val = u;
64
65 fasp_dcsr_sort(&mat);
66
67 fasp_solver_pardiso(&mat, &rhs, &sol, *ptrlvl);
68}
69#endif
70
89#if WITH_UMFPACK
90void fasp_fwrapper_dcsr_strumpack_(
91 INT* n, INT* nnz, INT* ia, INT* ja, REAL* a, REAL* b, REAL* u, INT* ptrlvl)
92{
93 dCSRmat mat; // coefficient matrix
94 dvector rhs, sol; // right-hand-side, solution
95
96 // set up coefficient matrix
97 mat.row = *n;
98 mat.col = *n;
99 mat.nnz = *nnz;
100 mat.IA = ia;
101 mat.JA = ja;
102 mat.val = a;
103
104 rhs.row = *n;
105 rhs.val = b;
106 sol.row = *n;
107 sol.val = u;
108
109 fasp_dcsr_sort(&mat);
110
111 fasp_solver_strumpack(&mat, &rhs, &sol, *ptrlvl);
112}
113#endif
114
137 INT* nnz,
138 INT* ia,
139 INT* ja,
140 REAL* a,
141 REAL* b,
142 REAL* u,
143 REAL* tol,
144 INT* maxit,
145 INT* ptrlvl)
146{
147 dCSRmat mat; // coefficient matrix
148 dvector rhs, sol; // right-hand-side, solution
149 AMG_param amgparam; // parameters for AMG
150
151 // setup AMG parameters
152 fasp_param_amg_init(&amgparam);
153
154 amgparam.tol = *tol;
155 amgparam.print_level = *ptrlvl;
156 amgparam.maxit = *maxit;
157
158 // set up coefficient matrix
159 mat.row = *n;
160 mat.col = *n;
161 mat.nnz = *nnz;
162 mat.IA = ia;
163 mat.JA = ja;
164 mat.val = a;
165
166 rhs.row = *n;
167 rhs.val = b;
168 sol.row = *n;
169 sol.val = u;
170
171 fasp_solver_amg(&mat, &rhs, &sol, &amgparam);
172}
173
196 INT* nnz,
197 INT* ia,
198 INT* ja,
199 REAL* a,
200 REAL* b,
201 REAL* u,
202 REAL* tol,
203 INT* maxit,
204 INT* ptrlvl)
205{
206 dCSRmat mat; // coefficient matrix
207 dvector rhs, sol; // right-hand-side, solution
208 ILU_param iluparam; // parameters for ILU
209 ITS_param itsparam; // parameters for itsolver
210
211 // setup ILU parameters
212 fasp_param_ilu_init(&iluparam);
213
214 iluparam.print_level = *ptrlvl;
215
216 // setup Krylov method parameters
217 fasp_param_solver_init(&itsparam);
218
219 itsparam.itsolver_type = SOLVER_VFGMRES;
220 itsparam.tol = *tol;
221 itsparam.maxit = *maxit;
222 itsparam.print_level = *ptrlvl;
223
224 // set up coefficient matrix
225 mat.row = *n;
226 mat.col = *n;
227 mat.nnz = *nnz;
228 mat.IA = ia;
229 mat.JA = ja;
230 mat.val = a;
231
232 rhs.row = *n;
233 rhs.val = b;
234 sol.row = *n;
235 sol.val = u;
236
237 fasp_solver_dcsr_krylov_ilu(&mat, &rhs, &sol, &itsparam, &iluparam);
238}
239
262 INT* nnz,
263 INT* ia,
264 INT* ja,
265 REAL* a,
266 REAL* b,
267 REAL* u,
268 REAL* tol,
269 INT* maxit,
270 INT* ptrlvl)
271{
272 dCSRmat mat; // coefficient matrix
273 dvector rhs, sol; // right-hand-side, solution
274 input_param inparam; // parameters from input files
275 AMG_param amgparam; // parameters for AMG
276 ITS_param itsparam; // parameters for itsolver
277 ILU_param iluparam; // parameters for ILU
278
280 char* inputfile = "ini/amg.dat"; // Added for fasp4ns 2022.04.08 --zcs
281 fasp_param_input(inputfile, &inparam);
282 fasp_param_init(&inparam, &itsparam, &amgparam, &iluparam, NULL);
283
284 itsparam.tol = *tol;
285 itsparam.maxit = *maxit;
286 itsparam.print_level = *ptrlvl;
287
288 // set up coefficient matrix
289 mat.row = *n;
290 mat.col = *n;
291 mat.nnz = *nnz;
292 mat.IA = ia;
293 mat.JA = ja;
294 mat.val = a;
295
296 rhs.row = *n;
297 rhs.val = b;
298 sol.row = *n;
299 sol.val = u;
300
301 fasp_solver_dcsr_krylov_amg(&mat, &rhs, &sol, &itsparam, &amgparam);
302}
303
327 INT* nnz,
328 INT* nb,
329 INT* ia,
330 INT* ja,
331 REAL* a,
332 REAL* b,
333 REAL* u,
334 REAL* tol,
335 INT* maxit,
336 INT* ptrlvl)
337{
338 dBSRmat mat; // coefficient matrix in BSR format
339 dvector rhs, sol; // right-hand-side, solution
340
341 ILU_param iluparam; // parameters for ILU
342 ITS_param itsparam; // parameters for itsolver
343
344 // setup ILU parameters
345 fasp_param_ilu_init(&iluparam);
346 iluparam.ILU_lfil = 0;
347 iluparam.print_level = *ptrlvl;
348
349 // setup Krylov method parameters
350 fasp_param_solver_init(&itsparam);
351
352 itsparam.itsolver_type = SOLVER_VFGMRES;
353 itsparam.tol = *tol;
354 itsparam.maxit = *maxit;
355 itsparam.print_level = *ptrlvl;
356
357 // set up coefficient matrix
358 mat.ROW = *n;
359 mat.COL = *n;
360 mat.NNZ = *nnz;
361 mat.nb = *nb;
362 mat.IA = ia;
363 mat.JA = ja;
364 mat.val = a;
365
366 rhs.row = *n * *nb;
367 rhs.val = b;
368 sol.row = *n * *nb;
369 sol.val = u;
370
371 // solve
372 fasp_solver_dbsr_krylov_ilu(&mat, &rhs, &sol, &itsparam, &iluparam);
373}
374
398 INT* nnz,
399 INT* nb,
400 INT* ia,
401 INT* ja,
402 REAL* a,
403 REAL* b,
404 REAL* u,
405 REAL* tol,
406 INT* maxit,
407 INT* ptrlvl)
408{
409 dBSRmat mat; // coefficient matrix in CSR format
410 dvector rhs, sol; // right-hand-side, solution
411
412 AMG_param amgparam; // parameters for AMG
413 ITS_param itsparam; // parameters for itsolver
414
415 // setup AMG parameters
416 fasp_param_amg_init(&amgparam);
417 amgparam.AMG_type = UA_AMG;
418 amgparam.print_level = *ptrlvl;
419
420 // setup Krylov method parameters
421 fasp_param_solver_init(&itsparam);
422 itsparam.tol = *tol;
423 itsparam.print_level = *ptrlvl;
424 itsparam.maxit = *maxit;
425 itsparam.itsolver_type = SOLVER_VFGMRES;
426
427 // set up coefficient matrix
428 mat.ROW = *n;
429 mat.COL = *n;
430 mat.NNZ = *nnz;
431 mat.nb = *nb;
432 mat.IA = ia;
433 mat.JA = ja;
434 mat.val = a;
435
436 rhs.row = *n * *nb;
437 rhs.val = b;
438 sol.row = *n * *nb;
439 sol.val = u;
440
441 // solve
442 fasp_solver_dbsr_krylov_amg(&mat, &rhs, &sol, &itsparam, &amgparam);
443}
444
445/*---------------------------------*/
446/*-- End of File --*/
447/*---------------------------------*/
void fasp_param_input(const char *fname, input_param *inparam)
Read input parameters from disk file.
Definition: AuxInput.c:86
void fasp_param_ilu_init(ILU_param *iluparam)
Initialize ILU parameters.
Definition: AuxParam.c:595
void fasp_param_init(const input_param *iniparam, ITS_param *itsparam, AMG_param *amgparam, ILU_param *iluparam, SWZ_param *swzparam)
Initialize parameters, global variables, etc.
Definition: AuxParam.c:306
void fasp_param_solver_init(ITS_param *itsparam)
Initialize ITS_param.
Definition: AuxParam.c:572
void fasp_param_amg_init(AMG_param *amgparam)
Initialize AMG parameters.
Definition: AuxParam.c:431
void fasp_dcsr_sort(dCSRmat *A)
Sort each row of A in ascending order w.r.t. column indices.
Definition: BlaSparseCSR.c:385
INT fasp_solver_amg(dCSRmat *A, dvector *b, dvector *x, AMG_param *param)
Solve Ax = b by algebraic multigrid methods.
Definition: SolAMG.c:49
INT fasp_solver_dbsr_krylov_ilu(dBSRmat *A, dvector *b, dvector *x, ITS_param *itparam, ILU_param *iluparam)
Solve Ax=b by ILUs preconditioned Krylov methods.
Definition: SolBSR.c:286
INT fasp_solver_dbsr_krylov_amg(dBSRmat *A, dvector *b, dvector *x, ITS_param *itparam, AMG_param *amgparam)
Solve Ax=b by AMG preconditioned Krylov methods.
Definition: SolBSR.c:349
INT fasp_solver_dcsr_krylov_ilu(dCSRmat *A, dvector *b, dvector *x, ITS_param *itparam, ILU_param *iluparam)
Solve Ax=b by ILUs preconditioned Krylov methods.
Definition: SolCSR.c:588
INT fasp_solver_dcsr_krylov_amg(dCSRmat *A, dvector *b, dvector *x, ITS_param *itparam, AMG_param *amgparam)
Solve Ax=b by AMG preconditioned Krylov methods.
Definition: SolCSR.c:476
void fasp_fwrapper_dcsr_amg_(INT *n, INT *nnz, INT *ia, INT *ja, REAL *a, REAL *b, REAL *u, REAL *tol, INT *maxit, INT *ptrlvl)
Solve Ax=b by Ruge and Stuben's classic AMG.
Definition: SolWrapper.c:136
void fasp_fwrapper_dbsr_krylov_ilu_(INT *n, INT *nnz, INT *nb, INT *ia, INT *ja, REAL *a, REAL *b, REAL *u, REAL *tol, INT *maxit, INT *ptrlvl)
Solve Ax=b by Krylov method preconditioned by block ILU in BSR format.
Definition: SolWrapper.c:326
void fasp_fwrapper_dcsr_krylov_ilu_(INT *n, INT *nnz, INT *ia, INT *ja, REAL *a, REAL *b, REAL *u, REAL *tol, INT *maxit, INT *ptrlvl)
Solve Ax=b by Krylov method preconditioned by ILUk.
Definition: SolWrapper.c:195
void fasp_fwrapper_dcsr_krylov_amg_(INT *n, INT *nnz, INT *ia, INT *ja, REAL *a, REAL *b, REAL *u, REAL *tol, INT *maxit, INT *ptrlvl)
Solve Ax=b by Krylov method preconditioned by classic AMG.
Definition: SolWrapper.c:261
void fasp_fwrapper_dbsr_krylov_amg_(INT *n, INT *nnz, INT *nb, INT *ia, INT *ja, REAL *a, REAL *b, REAL *u, REAL *tol, INT *maxit, INT *ptrlvl)
Solve Ax=b by Krylov method preconditioned by block AMG in BSR format.
Definition: SolWrapper.c:397
INT fasp_solver_pardiso(dCSRmat *ptrA, dvector *b, dvector *u, const SHORT prtlvl)
Solve Ax=b by PARDISO directly.
Definition: XtrPardiso.c:45
INT fasp_solver_strumpack(dCSRmat *ptrA, dvector *b, dvector *u, const SHORT prtlvl)
Solve Au=b by UMFpack.
Definition: XtrStrumpack.c:41
Main header file for the FASP project.
#define REAL
Definition: fasp.h:75
#define INT
Definition: fasp.h:72
Header file for FASP block matrices.
#define SOLVER_VFGMRES
Definition: fasp_const.h:108
#define UA_AMG
Definition: fasp_const.h:165
Parameters for AMG methods.
Definition: fasp.h:455
SHORT print_level
print level for AMG
Definition: fasp.h:461
SHORT AMG_type
type of AMG method
Definition: fasp.h:458
REAL tol
stopping tolerance for AMG solver
Definition: fasp.h:467
INT maxit
max number of iterations of AMG
Definition: fasp.h:464
Parameters for ILU.
Definition: fasp.h:404
INT ILU_lfil
level of fill-in for ILUk
Definition: fasp.h:413
SHORT print_level
print level
Definition: fasp.h:407
Parameters for iterative solvers.
Definition: fasp.h:386
SHORT itsolver_type
Definition: fasp.h:389
SHORT print_level
Definition: fasp.h:388
REAL tol
Definition: fasp.h:395
INT maxit
Definition: fasp.h:394
Block sparse row storage matrix of REAL type.
Definition: fasp_block.h:34
INT COL
number of cols of sub-blocks in matrix A, N
Definition: fasp_block.h:40
INT NNZ
number of nonzero sub-blocks in matrix A, NNZ
Definition: fasp_block.h:43
REAL * val
Definition: fasp_block.h:57
INT nb
dimension of each sub-block
Definition: fasp_block.h:46
INT * IA
integer array of row pointers, the size is ROW+1
Definition: fasp_block.h:60
INT ROW
number of rows of sub-blocks in matrix A, M
Definition: fasp_block.h:37
INT * JA
Definition: fasp_block.h:64
Sparse matrix of REAL type in CSR format.
Definition: fasp.h:151
INT col
column of matrix A, n
Definition: fasp.h:157
REAL * val
nonzero entries of A
Definition: fasp.h:169
INT row
row number of matrix A, m
Definition: fasp.h:154
INT * IA
integer array of row pointers, the size is m+1
Definition: fasp.h:163
INT nnz
number of nonzero entries
Definition: fasp.h:160
INT * JA
integer array of column indexes, the size is nnz
Definition: fasp.h:166
Vector with n entries of REAL type.
Definition: fasp.h:354
REAL * val
actual vector entries
Definition: fasp.h:360
INT row
number of rows
Definition: fasp.h:357
Input parameters.
Definition: fasp.h:1125