Fast Auxiliary Space Preconditioning 2.7.7 Aug/28/2022
Loading...
Searching...
No Matches
SolCSR.c
Go to the documentation of this file.
1
18#include <time.h>
19
20#ifdef _OPENMP
21#include <omp.h>
22#endif
23
24#include "fasp.h"
25#include "fasp_functs.h"
26
27/*---------------------------------*/
28/*-- Declare Private Functions --*/
29/*---------------------------------*/
30
31#include "KryUtil.inl"
32
33/*---------------------------------*/
34/*-- Public Functions --*/
35/*---------------------------------*/
36
57 ITS_param* itparam)
58{
59 const SHORT prtlvl = itparam->print_level;
60 const SHORT itsolver_type = itparam->itsolver_type;
61 const SHORT stop_type = itparam->stop_type;
62 const SHORT restart = itparam->restart;
63 const INT MaxIt = itparam->maxit;
64 const REAL tol = itparam->tol;
65 const REAL abstol = itparam->abstol;
66
67 /* Local Variables */
68 REAL solve_start, solve_end;
69 INT iter;
70
71#if DEBUG_MODE > 0
72 printf("### DEBUG: [-Begin-] %s ...\n", __FUNCTION__);
73 printf("### DEBUG: rhs/sol size: %d %d\n", b->row, x->row);
74#endif
75
76 fasp_gettime(&solve_start);
77
78 /* check matrix data */
80
81 /* Safe-guard checks on parameters */
82 ITS_CHECK(MaxIt, tol);
83
84 /* Choose a desirable Krylov iterative solver */
85 switch (itsolver_type) {
86 case SOLVER_CG:
87 iter = fasp_solver_dcsr_pcg(A, b, x, pc, tol, abstol, MaxIt, stop_type,
88 prtlvl);
89 break;
90
91 case SOLVER_BiCGstab:
92 iter = fasp_solver_dcsr_pbcgs(A, b, x, pc, tol, abstol, MaxIt, stop_type,
93 prtlvl);
94 break;
95
96 case SOLVER_MinRes:
97 iter = fasp_solver_dcsr_pminres(A, b, x, pc, tol, abstol, MaxIt, stop_type,
98 prtlvl);
99 break;
100
101 case SOLVER_GMRES:
102 iter = fasp_solver_dcsr_pgmres(A, b, x, pc, tol, abstol, MaxIt, restart,
103 stop_type, prtlvl);
104 break;
105
106 case SOLVER_VGMRES:
107 iter = fasp_solver_dcsr_pvgmres(A, b, x, pc, tol, abstol, MaxIt, restart,
108 stop_type, prtlvl);
109 break;
110
111 case SOLVER_VFGMRES:
112 iter = fasp_solver_dcsr_pvfgmres(A, b, x, pc, tol, abstol, MaxIt, restart,
113 stop_type, prtlvl);
114 break;
115
116 case SOLVER_GCG:
117 iter = fasp_solver_dcsr_pgcg(A, b, x, pc, tol, abstol, MaxIt, stop_type,
118 prtlvl);
119 break;
120
121 case SOLVER_GCR:
122 iter = fasp_solver_dcsr_pgcr(A, b, x, pc, tol, abstol, MaxIt, restart,
123 stop_type, prtlvl);
124 break;
125
126 default:
127 printf("### ERROR: Unknown iterative solver type %d! [%s]\n", itsolver_type,
128 __FUNCTION__);
129 return ERROR_SOLVER_TYPE;
130 }
131
132 if ((prtlvl >= PRINT_SOME) && (iter >= 0)) {
133 fasp_gettime(&solve_end);
134 fasp_cputime("Iterative method", solve_end - solve_start);
135 }
136
137#if DEBUG_MODE > 0
138 printf("### DEBUG: [--End--] %s ...\n", __FUNCTION__);
139#endif
140
141 return iter;
142}
143
164 ITS_param* itparam)
165{
166 const SHORT prtlvl = itparam->print_level;
167 const SHORT itsolver_type = itparam->itsolver_type;
168 const SHORT stop_type = itparam->stop_type;
169 const SHORT restart = itparam->restart;
170 const INT MaxIt = itparam->maxit;
171 const REAL tol = itparam->tol;
172
173 /* Local Variables */
174 REAL solve_start, solve_end;
175 INT iter;
176
177#if DEBUG_MODE > 0
178 printf("### DEBUG: [-Begin-] %s ...\n", __FUNCTION__);
179 printf("### DEBUG: rhs/sol size: %d %d\n", b->row, x->row);
180#endif
181
182 fasp_gettime(&solve_start);
183
184 /* check matrix data */
186
187 /* Safe-guard checks on parameters */
188 ITS_CHECK(MaxIt, tol);
189
190 /* Choose a desirable Krylov iterative solver */
191 switch (itsolver_type) {
192 case SOLVER_CG:
193 iter = fasp_solver_dcsr_spcg(A, b, x, pc, tol, MaxIt, stop_type, prtlvl);
194 break;
195
196 case SOLVER_BiCGstab:
197 iter = fasp_solver_dcsr_spbcgs(A, b, x, pc, tol, MaxIt, stop_type, prtlvl);
198 break;
199
200 case SOLVER_MinRes:
201 iter =
202 fasp_solver_dcsr_spminres(A, b, x, pc, tol, MaxIt, stop_type, prtlvl);
203 break;
204
205 case SOLVER_GMRES:
206 iter = fasp_solver_dcsr_spgmres(A, b, x, pc, tol, MaxIt, restart, stop_type,
207 prtlvl);
208 break;
209
210 case SOLVER_VGMRES:
211 iter = fasp_solver_dcsr_spvgmres(A, b, x, pc, tol, MaxIt, restart,
212 stop_type, prtlvl);
213 break;
214
215 default:
216 printf("### ERROR: Unknown iterative solver type %d! [%s]\n", itsolver_type,
217 __FUNCTION__);
218 return ERROR_SOLVER_TYPE;
219 }
220
221 if ((prtlvl >= PRINT_SOME) && (iter >= 0)) {
222 fasp_gettime(&solve_end);
223 fasp_cputime("Iterative method", solve_end - solve_start);
224 }
225
226#if DEBUG_MODE > 0
227 printf("### DEBUG: [--End--] %s ...\n", __FUNCTION__);
228#endif
229
230 return iter;
231}
232
250{
251 const SHORT prtlvl = itparam->print_level;
252
253 /* Local Variables */
254 INT status = FASP_SUCCESS;
255 REAL solve_start, solve_end;
256
257#if DEBUG_MODE > 0
258 printf("### DEBUG: [-Begin-] %s ...\n", __FUNCTION__);
259 printf("### DEBUG: matrix size: %d %d %d\n", A->row, A->col, A->nnz);
260 printf("### DEBUG: rhs/sol size: %d %d\n", b->row, x->row);
261#endif
262
263 fasp_gettime(&solve_start);
264
265 status = fasp_solver_dcsr_itsolver(A, b, x, NULL, itparam);
266
267 if (prtlvl >= PRINT_MIN) {
268 fasp_gettime(&solve_end);
269 fasp_cputime("Krylov method totally", solve_end - solve_start);
270 }
271
272#if DEBUG_MODE > 0
273 printf("### DEBUG: [--End--] %s ...\n", __FUNCTION__);
274#endif
275
276 return status;
277}
278
296{
297 const SHORT prtlvl = itparam->print_level;
298
299 /* Local Variables */
300 INT status = FASP_SUCCESS;
301 REAL solve_start, solve_end;
302
303#if DEBUG_MODE > 0
304 printf("### DEBUG: [-Begin-] %s ...\n", __FUNCTION__);
305 printf("### DEBUG: matrix size: %d %d %d\n", A->row, A->col, A->nnz);
306 printf("### DEBUG: rhs/sol size: %d %d\n", b->row, x->row);
307#endif
308
309 fasp_gettime(&solve_start);
310
311 status = fasp_solver_dcsr_itsolver_s(A, b, x, NULL, itparam);
312
313 if (prtlvl >= PRINT_MIN) {
314 fasp_gettime(&solve_end);
315 fasp_cputime("Krylov method totally", solve_end - solve_start);
316 }
317
318#if DEBUG_MODE > 0
319 printf("### DEBUG: [--End--] %s ...\n", __FUNCTION__);
320#endif
321
322 return status;
323}
324
342{
343 const SHORT prtlvl = itparam->print_level;
344
345 /* Local Variables */
346 INT status = FASP_SUCCESS;
347 REAL solve_start, solve_end;
348
349#if DEBUG_MODE > 0
350 printf("### DEBUG: [-Begin-] %s ...\n", __FUNCTION__);
351 printf("### DEBUG: matrix size: %d %d %d\n", A->row, A->col, A->nnz);
352 printf("### DEBUG: rhs/sol size: %d %d\n", b->row, x->row);
353#endif
354
355 fasp_gettime(&solve_start);
356
357 // setup preconditioner
358 dvector diag;
359 fasp_dcsr_getdiag(0, A, &diag);
360
361 precond pc;
362 pc.data = &diag;
364
365 // call iterative solver
366 status = fasp_solver_dcsr_itsolver(A, b, x, &pc, itparam);
367
368 if (prtlvl >= PRINT_MIN) {
369 fasp_gettime(&solve_end);
370 fasp_cputime("Diag_Krylov method totally", solve_end - solve_start);
371 }
372
373 fasp_dvec_free(&diag);
374
375#if DEBUG_MODE > 0
376 printf("### DEBUG: [--End--] %s ...\n", __FUNCTION__);
377#endif
378
379 return status;
380}
381
402 SWZ_param* schparam)
403{
404 SWZ_param swzparam;
405 swzparam.SWZ_mmsize = schparam->SWZ_mmsize;
406 swzparam.SWZ_maxlvl = schparam->SWZ_maxlvl;
407 swzparam.SWZ_type = schparam->SWZ_type;
408 swzparam.SWZ_blksolver = schparam->SWZ_blksolver;
409
410 const SHORT prtlvl = itparam->print_level;
411
412 REAL setup_start, setup_end;
413 REAL solve_start, solve_end;
414 INT status = FASP_SUCCESS;
415
416#if DEBUG_MODE > 0
417 printf("### DEBUG: [-Begin-] %s ...\n", __FUNCTION__);
418 printf("### DEBUG: matrix size: %d %d %d\n", A->row, A->col, A->nnz);
419 printf("### DEBUG: rhs/sol size: %d %d\n", b->row, x->row);
420#endif
421
422 fasp_gettime(&solve_start);
423 fasp_gettime(&setup_start);
424
425 // setup preconditioner
427
428 // symmetrize the matrix (get rid of this later)
430
431 // construct Schwarz precondtioner
433 fasp_swz_dcsr_setup(&SWZ_data, &swzparam);
434
435 fasp_gettime(&setup_end);
436 printf("SWZ_Krylov method setup %f seconds.\n", setup_end - setup_start);
437
438 precond prec;
439 prec.data = &SWZ_data;
440 prec.fct = fasp_precond_swz;
441
442 // solver part
443 status = fasp_solver_dcsr_itsolver(A, b, x, &prec, itparam);
444
445 if (prtlvl > PRINT_NONE) {
446 fasp_gettime(&solve_end);
447 fasp_cputime("SWZ_Krylov method totally", solve_end - solve_start);
448 }
449
450#if DEBUG_MODE > 0
451 printf("### DEBUG: [--End--] %s ...\n", __FUNCTION__);
452#endif
453
455
456 return status;
457}
458
477 AMG_param* amgparam)
478{
479 const SHORT prtlvl = itparam->print_level;
480 const SHORT max_levels = amgparam->max_levels;
481 const INT nnz = A->nnz, m = A->row, n = A->col;
482
483 /* Local Variables */
484 INT status = FASP_SUCCESS;
485 REAL solve_start, solve_end;
486
487#if MULTI_COLOR_ORDER
488 A->color = 0;
489 A->IC = NULL;
490 A->ICMAP = NULL;
491#endif
492
493#if DEBUG_MODE > 0
494 printf("### DEBUG: [-Begin-] %s ...\n", __FUNCTION__);
495 printf("### DEBUG: matrix size: %d %d %d\n", A->row, A->col, A->nnz);
496 printf("### DEBUG: rhs/sol size: %d %d\n", b->row, x->row);
497#endif
498
499 fasp_gettime(&solve_start);
500
501 // initialize A, b, x for mgl[0]
502 AMG_data* mgl = fasp_amg_data_create(max_levels);
503 mgl[0].A = fasp_dcsr_create(m, n, nnz);
504 fasp_dcsr_cp(A, &mgl[0].A);
505 mgl[0].b = fasp_dvec_create(n);
506 mgl[0].x = fasp_dvec_create(n);
507
508 // setup preconditioner
509 switch (amgparam->AMG_type) {
510
511 case SA_AMG: // Smoothed Aggregation AMG
512 status = fasp_amg_setup_sa(mgl, amgparam);
513 break;
514
515 case UA_AMG: // Unsmoothed Aggregation AMG
516 status = fasp_amg_setup_ua(mgl, amgparam);
517 break;
518
519 default: // Classical AMG
520 status = fasp_amg_setup_rs(mgl, amgparam);
521 }
522
523#if DEBUG_MODE > 1
525#endif
526
527 if (status < 0) goto FINISHED;
528
529 // setup preconditioner
530 precond_data pcdata;
531 fasp_param_amg_to_prec(&pcdata, amgparam);
532 pcdata.max_levels = mgl[0].num_levels;
533 pcdata.mgl_data = mgl;
534
535 precond pc;
536 pc.data = &pcdata;
537
538 if (itparam->precond_type == PREC_FMG) {
539 pc.fct = fasp_precond_famg; // Full AMG
540 } else {
541 switch (amgparam->cycle_type) {
542 case AMLI_CYCLE: // AMLI cycle
544 break;
545 case NL_AMLI_CYCLE: // Nonlinear AMLI
547 break;
548 default: // V,W-cycles or hybrid cycles
550 }
551 }
552
553 // call iterative solver
554 status = fasp_solver_dcsr_itsolver(A, b, x, &pc, itparam);
555
556 if (prtlvl >= PRINT_MIN) {
557 fasp_gettime(&solve_end);
558 fasp_cputime("AMG_Krylov method totally", solve_end - solve_start);
559 }
560
561FINISHED:
562 fasp_amg_data_free(mgl, amgparam);
563
564#if DEBUG_MODE > 0
565 printf("### DEBUG: [--End--] %s ...\n", __FUNCTION__);
566#endif
567
568 return status;
569}
570
589 ILU_param* iluparam)
590{
591 const SHORT prtlvl = itparam->print_level;
592
593 /* Local Variables */
594 INT status = FASP_SUCCESS;
595 REAL solve_start, solve_end, solve_time;
596
597#if DEBUG_MODE > 0
598 printf("### DEBUG: [-Begin-] %s ...\n", __FUNCTION__);
599 printf("### DEBUG: matrix size: %d %d %d\n", A->row, A->col, A->nnz);
600 printf("### DEBUG: rhs/sol size: %d %d\n", b->row, x->row);
601#endif
602
603 fasp_gettime(&solve_start);
604
605 // ILU setup for whole matrix
606 ILU_data LU;
607 if ((status = fasp_ilu_dcsr_setup(A, &LU, iluparam)) < 0) goto FINISHED;
608
609 // check iludata
610 if ((status = fasp_mem_iludata_check(&LU)) < 0) goto FINISHED;
611
612 // set preconditioner
613 precond pc;
614 pc.data = &LU;
616
617 // call iterative solver
618 status = fasp_solver_dcsr_itsolver(A, b, x, &pc, itparam);
619
620 if (prtlvl >= PRINT_MIN) {
621 fasp_gettime(&solve_end);
622 solve_time = solve_end - solve_start;
623
624 switch (iluparam->ILU_type) {
625 case ILUt:
626 fasp_cputime("ILUt_Krylov method totally", solve_time);
627 break;
628 case ILUtp:
629 fasp_cputime("ILUtp_Krylov method totally", solve_time);
630 break;
631 default: // ILUk
632 fasp_cputime("ILUk_Krylov method totally", solve_time);
633 }
634 }
635
636FINISHED:
638
639#if DEBUG_MODE > 0
640 printf("### DEBUG: [--End--] %s ...\n", __FUNCTION__);
641#endif
642
643 return status;
644}
645
669 ITS_param* itparam, ILU_param* iluparam, dCSRmat* M)
670{
671 const SHORT prtlvl = itparam->print_level;
672
673 /* Local Variables */
674 REAL solve_start, solve_end, solve_time;
675 INT status = FASP_SUCCESS;
676
677#if DEBUG_MODE > 0
678 printf("### DEBUG: [-Begin-] %s ...\n", __FUNCTION__);
679 printf("### DEBUG: matrix size: %d %d %d\n", A->row, A->col, A->nnz);
680 printf("### DEBUG: rhs/sol size: %d %d\n", b->row, x->row);
681#endif
682
683 fasp_gettime(&solve_start);
684
685 // ILU setup for M
686 ILU_data LU;
687 if ((status = fasp_ilu_dcsr_setup(M, &LU, iluparam)) < 0) goto FINISHED;
688
689 // check iludata
690 if ((status = fasp_mem_iludata_check(&LU)) < 0) goto FINISHED;
691
692 // set precondtioner
693 precond pc;
694 pc.data = &LU;
696
697 // call iterative solver
698 status = fasp_solver_dcsr_itsolver(A, b, x, &pc, itparam);
699
700 if (prtlvl >= PRINT_MIN) {
701 fasp_gettime(&solve_end);
702 solve_time = solve_end - solve_start;
703
704 switch (iluparam->ILU_type) {
705 case ILUt:
706 fasp_cputime("ILUt_Krylov method", solve_time);
707 break;
708 case ILUtp:
709 fasp_cputime("ILUtp_Krylov method", solve_time);
710 break;
711 default: // ILUk
712 fasp_cputime("ILUk_Krylov method", solve_time);
713 }
714 }
715
716FINISHED:
718
719#if DEBUG_MODE > 0
720 printf("### DEBUG: [--End--] %s ...\n", __FUNCTION__);
721#endif
722
723 return status;
724}
725
749 ITS_param* itparam, AMG_param* amgparam,
750 dCSRmat* A_nk, dCSRmat* P_nk, dCSRmat* R_nk)
751{
752 const SHORT prtlvl = itparam->print_level;
753 const SHORT max_levels = amgparam->max_levels;
754 const INT nnz = A->nnz, m = A->row, n = A->col;
755
756 /* Local Variables */
757 INT status = FASP_SUCCESS;
758 REAL solve_start, solve_end, solve_time;
759
760#if MULTI_COLOR_ORDER
761 A->color = 0;
762 A->IC = NULL;
763 A->ICMAP = NULL;
764#endif
765
766#if DEBUG_MODE > 0
767 printf("### DEBUG: [-Begin-] %s ...\n", __FUNCTION__);
768 printf("### DEBUG: matrix size: %d %d %d\n", A->row, A->col, A->nnz);
769 printf("### DEBUG: rhs/sol size: %d %d\n", b->row, x->row);
770#endif
771
772 fasp_gettime(&solve_start);
773
774 // initialize A, b, x for mgl[0]
775 AMG_data* mgl = fasp_amg_data_create(max_levels);
776 mgl[0].A = fasp_dcsr_create(m, n, nnz);
777 fasp_dcsr_cp(A, &mgl[0].A);
778 mgl[0].b = fasp_dvec_create(n);
779 mgl[0].x = fasp_dvec_create(n);
780
781 // setup preconditioner
782 switch (amgparam->AMG_type) {
783
784 case SA_AMG: // Smoothed Aggregation AMG
785 status = fasp_amg_setup_sa(mgl, amgparam);
786 break;
787
788 case UA_AMG: // Unsmoothed Aggregation AMG
789 status = fasp_amg_setup_ua(mgl, amgparam);
790 break;
791
792 default: // Classical AMG
793 status = fasp_amg_setup_rs(mgl, amgparam);
794 }
795
796#if DEBUG_MODE > 1
798#endif
799
800 if (status < 0) goto FINISHED;
801
802 // setup preconditioner
803 precond_data pcdata;
804 fasp_param_amg_to_prec(&pcdata, amgparam);
805 pcdata.max_levels = mgl[0].num_levels;
806 pcdata.mgl_data = mgl;
807
808 // near kernel space
809#if WITH_UMFPACK // use UMFPACK directly
810 dCSRmat A_tran;
811 A_tran = fasp_dcsr_create(A_nk->row, A_nk->col, A_nk->nnz);
812 fasp_dcsr_transz(A_nk, NULL, &A_tran);
813 // It is equivalent to do transpose and then sort
814 // fasp_dcsr_trans(A_nk, &A_tran);
815 // fasp_dcsr_sort(&A_tran);
816 pcdata.A_nk = &A_tran;
817#else
818 pcdata.A_nk = A_nk;
819#endif
820
821 pcdata.P_nk = P_nk;
822 pcdata.R_nk = R_nk;
823
824 precond pc;
825 pc.data = &pcdata;
827
828 // call iterative solver
829 status = fasp_solver_dcsr_itsolver(A, b, x, &pc, itparam);
830
831 if (prtlvl >= PRINT_MIN) {
832 fasp_gettime(&solve_end);
833 solve_time = solve_end - solve_start;
834 fasp_cputime("AMG_NK_Krylov method", solve_time);
835 }
836
837FINISHED:
838 fasp_amg_data_free(mgl, amgparam);
839
840#if WITH_UMFPACK // use UMFPACK directly
841 fasp_dcsr_free(&A_tran);
842#endif
843
844#if DEBUG_MODE > 0
845 printf("### DEBUG: [--End--] %s ...\n", __FUNCTION__);
846#endif
847
848 return status;
849}
850
851/*---------------------------------*/
852/*-- End of File --*/
853/*---------------------------------*/
void fasp_mem_usage(void)
Show total allocated memory currently.
Definition: AuxMemory.c:183
SHORT fasp_mem_iludata_check(const ILU_data *iludata)
Check wether a ILU_data has enough work space.
Definition: AuxMemory.c:203
void fasp_cputime(const char *message, const REAL cputime)
Print CPU walltime.
Definition: AuxMessage.c:179
void fasp_param_amg_to_prec(precond_data *pcdata, const AMG_param *amgparam)
Set precond_data with AMG_param.
Definition: AuxParam.c:782
void fasp_gettime(REAL *time)
Get system time.
Definition: AuxTiming.c:37
void fasp_dvec_free(dvector *u)
Free vector data space of REAL type.
Definition: AuxVector.c:145
dvector fasp_dvec_create(const INT m)
Create dvector data space of REAL type.
Definition: AuxVector.c:62
SHORT fasp_ilu_dcsr_setup(dCSRmat *A, ILU_data *iludata, ILU_param *iluparam)
Get ILU decomposition of a CSR matrix A.
INT fasp_swz_dcsr_setup(SWZ_data *swzdata, SWZ_param *swzparam)
Setup phase for the Schwarz methods.
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_shift(dCSRmat *A, const INT offset)
Re-index a REAL matrix in CSR format to make the index starting from 0 or 1.
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
void fasp_dcsr_getdiag(INT n, const dCSRmat *A, dvector *diag)
Get first n diagonal entries of a CSR matrix A.
Definition: BlaSparseCSR.c:537
dCSRmat fasp_dcsr_sympart(dCSRmat *A)
Get symmetric part of a dCSRmat matrix.
void fasp_check_dCSRmat(const dCSRmat *A)
Check whether an dCSRmat matrix is supported or not.
INT fasp_solver_dcsr_pbcgs(dCSRmat *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 CSR matrix.
Definition: KryPbcgs.c:62
INT fasp_solver_dcsr_pcg(dCSRmat *A, dvector *b, dvector *u, precond *pc, const REAL tol, const REAL abstol, const INT MaxIt, const SHORT StopType, const SHORT PrtLvl)
Preconditioned conjugate gradient method for solving Au=b.
Definition: KryPcg.c:96
INT fasp_solver_dcsr_pgcg(dCSRmat *A, dvector *b, dvector *u, precond *pc, const REAL tol, const REAL abstol, const INT MaxIt, const SHORT StopType, const SHORT PrtLvl)
Preconditioned generilzed conjugate gradient (GCG) method for solving Au=b.
Definition: KryPgcg.c:60
INT fasp_solver_dcsr_pgcr(dCSRmat *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)
A preconditioned GCR method for solving Au=b.
Definition: KryPgcr.c:55
INT fasp_solver_dcsr_pgmres(dCSRmat *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 for solving Au=b.
Definition: KryPgmres.c:66
INT fasp_solver_dcsr_pminres(dCSRmat *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:61
INT fasp_solver_dcsr_pvfgmres(dCSRmat *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:67
INT fasp_solver_dcsr_pvgmres(dCSRmat *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:66
INT fasp_solver_dcsr_spbcgs(const dCSRmat *A, const dvector *b, dvector *u, precond *pc, const REAL tol, const INT MaxIt, const SHORT StopType, const SHORT PrtLvl)
Preconditioned BiCGstab method for solving Au=b with safety net.
Definition: KrySPbcgs.c:61
INT fasp_solver_dcsr_spcg(const dCSRmat *A, const dvector *b, dvector *u, precond *pc, const REAL tol, const INT MaxIt, const SHORT StopType, const SHORT PrtLvl)
Preconditioned conjugate gradient method for solving Au=b with safety net.
Definition: KrySPcg.c:60
INT fasp_solver_dcsr_spgmres(const dCSRmat *A, const dvector *b, dvector *x, precond *pc, const REAL tol, const INT MaxIt, SHORT restart, const SHORT StopType, const SHORT PrtLvl)
Preconditioned GMRES method for solving Au=b with safe-guard.
Definition: KrySPgmres.c:66
INT fasp_solver_dcsr_spminres(const dCSRmat *A, const dvector *b, dvector *u, precond *pc, const REAL tol, const INT MaxIt, const SHORT StopType, const SHORT PrtLvl)
A preconditioned minimal residual (Minres) method for solving Au=b with safety net.
Definition: KrySPminres.c:60
INT fasp_solver_dcsr_spvgmres(const dCSRmat *A, const dvector *b, dvector *x, precond *pc, const REAL tol, const INT MaxIt, SHORT restart, const SHORT StopType, const SHORT PrtLvl)
Solve "Ax=b" using PGMRES(right preconditioned) iterative method in which the restart parameter can b...
Definition: KrySPvgmres.c:68
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_amg(REAL *r, REAL *z, void *data)
AMG preconditioner.
Definition: PreCSR.c:416
void fasp_precond_ilu(REAL *r, REAL *z, void *data)
ILU preconditioner.
Definition: PreCSR.c:198
void fasp_precond_swz(REAL *r, REAL *z, void *data)
get z from r by Schwarz
Definition: PreCSR.c:371
void fasp_precond_namli(REAL *r, REAL *z, void *data)
Nonlinear AMLI AMG preconditioner.
Definition: PreCSR.c:515
void fasp_precond_amg_nk(REAL *r, REAL *z, void *data)
AMG with extra near kernel solve as preconditioner.
Definition: PreCSR.c:548
void fasp_precond_diag(REAL *r, REAL *z, void *data)
Diagonal preconditioner z=inv(D)*r.
Definition: PreCSR.c:172
void fasp_precond_famg(REAL *r, REAL *z, void *data)
Full AMG preconditioner.
Definition: PreCSR.c:449
void fasp_precond_amli(REAL *r, REAL *z, void *data)
AMLI AMG preconditioner.
Definition: PreCSR.c:482
void fasp_swz_data_free(SWZ_data *swzdata)
Free SWZ_data data memeory space.
Definition: PreDataInit.c:494
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
void fasp_ilu_data_free(ILU_data *iludata)
Create ILU_data sturcture.
Definition: PreDataInit.c:445
INT fasp_solver_dcsr_itsolver(dCSRmat *A, dvector *b, dvector *x, precond *pc, ITS_param *itparam)
Solve Ax=b by preconditioned Krylov methods for CSR matrices.
Definition: SolCSR.c:56
INT fasp_solver_dcsr_krylov_amg_nk(dCSRmat *A, dvector *b, dvector *x, ITS_param *itparam, AMG_param *amgparam, dCSRmat *A_nk, dCSRmat *P_nk, dCSRmat *R_nk)
Solve Ax=b by AMG preconditioned Krylov methods with an extra near kernel solve.
Definition: SolCSR.c:748
INT fasp_solver_dcsr_krylov_ilu_M(dCSRmat *A, dvector *b, dvector *x, ITS_param *itparam, ILU_param *iluparam, dCSRmat *M)
Solve Ax=b by ILUs preconditioned Krylov methods: ILU of M as preconditioner.
Definition: SolCSR.c:668
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_diag(dCSRmat *A, dvector *b, dvector *x, ITS_param *itparam)
Solve Ax=b by diagonal preconditioned Krylov methods.
Definition: SolCSR.c:341
INT fasp_solver_dcsr_krylov_s(dCSRmat *A, dvector *b, dvector *x, ITS_param *itparam)
Solve Ax=b by standard Krylov methods with safe-net for CSR matrices.
Definition: SolCSR.c:295
INT fasp_solver_dcsr_krylov(dCSRmat *A, dvector *b, dvector *x, ITS_param *itparam)
Solve Ax=b by standard Krylov methods for CSR matrices.
Definition: SolCSR.c:249
INT fasp_solver_dcsr_krylov_swz(dCSRmat *A, dvector *b, dvector *x, ITS_param *itparam, SWZ_param *schparam)
Solve Ax=b by overlapping Schwarz Krylov methods.
Definition: SolCSR.c:401
INT fasp_solver_dcsr_itsolver_s(dCSRmat *A, dvector *b, dvector *x, precond *pc, ITS_param *itparam)
Solve Ax=b by preconditioned Krylov methods with safe-net for CSR matrices.
Definition: SolCSR.c:163
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
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
#define AMLI_CYCLE
Definition: fasp_const.h:181
#define NL_AMLI_CYCLE
Definition: fasp_const.h:182
#define SOLVER_GMRES
Definition: fasp_const.h:106
#define FASP_SUCCESS
Definition of return status and error messages.
Definition: fasp_const.h:19
#define PREC_FMG
Definition: fasp_const.h:142
#define ILUtp
Definition: fasp_const.h:151
#define ERROR_SOLVER_TYPE
Definition: fasp_const.h:41
#define SOLVER_VGMRES
Definition: fasp_const.h:107
#define SOLVER_BiCGstab
Definition: fasp_const.h:104
#define ILUt
Definition: fasp_const.h:150
#define SA_AMG
Definition: fasp_const.h:164
#define PRINT_SOME
Definition: fasp_const.h:75
#define SOLVER_MinRes
Definition: fasp_const.h:105
#define SOLVER_GCR
Definition: fasp_const.h:110
#define SOLVER_VFGMRES
Definition: fasp_const.h:108
#define SOLVER_CG
Definition: fasp_const.h:103
#define PRINT_NONE
Print level for all subroutines – not including DEBUG output.
Definition: fasp_const.h:73
#define PRINT_MIN
Definition: fasp_const.h:74
#define SOLVER_GCG
Definition: fasp_const.h:109
#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
SHORT num_levels
number of levels in use <= max_levels
Definition: fasp.h:812
Parameters for AMG methods.
Definition: fasp.h:455
SHORT AMG_type
type of AMG method
Definition: fasp.h:458
SHORT cycle_type
type of AMG cycle
Definition: fasp.h:476
SHORT max_levels
max number of levels of AMG
Definition: fasp.h:470
Data for ILU setup.
Definition: fasp.h:651
Parameters for ILU.
Definition: fasp.h:404
SHORT ILU_type
ILU type for decomposition.
Definition: fasp.h:410
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
Data for Schwarz methods.
Definition: fasp.h:726
dCSRmat A
pointer to the original coefficient matrix
Definition: fasp.h:731
Parameters for Schwarz method.
Definition: fasp.h:430
INT SWZ_mmsize
maximal size of blocks
Definition: fasp.h:442
INT SWZ_maxlvl
maximal level for constructing the blocks
Definition: fasp.h:439
INT SWZ_blksolver
type of Schwarz block solver
Definition: fasp.h:445
SHORT SWZ_type
type for Schwarz method
Definition: fasp.h:436
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
Data for preconditioners.
Definition: fasp.h:894
AMG_data * mgl_data
AMG preconditioner data.
Definition: fasp.h:954
dCSRmat * A_nk
Matrix data for near kernel.
Definition: fasp.h:965
dCSRmat * R_nk
Restriction for near kernel.
Definition: fasp.h:971
SHORT max_levels
max number of AMG levels
Definition: fasp.h:906
dCSRmat * P_nk
Prolongation for near kernel.
Definition: fasp.h:968
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