Fast Auxiliary Space Preconditioning 2.7.7 Aug/28/2022
Loading...
Searching...
No Matches
SolAMG.c
Go to the documentation of this file.
1
16#include <time.h>
17
18#include "fasp.h"
19#include "fasp_functs.h"
20
21/*---------------------------------*/
22/*-- Public Functions --*/
23/*---------------------------------*/
24
50{
51 const REAL tol = param->tol;
52 const SHORT max_levels = param->max_levels;
53 const SHORT prtlvl = param->print_level;
54 const SHORT amg_type = param->AMG_type;
55 const SHORT cycle_type = param->cycle_type;
56 const INT maxit = param->maxit;
57 const INT nnz = A->nnz, m = A->row, n = A->col;
58
59 // local variables
60 SHORT status;
61 INT iter = 0;
62 AMG_data* mgl = fasp_amg_data_create(max_levels);
63 REAL AMG_start = 0, AMG_end;
64
65#if MULTI_COLOR_ORDER
66 A->color = 0;
67 A->IC = NULL;
68 A->ICMAP = NULL;
69#endif
70
71#if DEBUG_MODE > 0
72 printf("### DEBUG: [-Begin-] %s ...\n", __FUNCTION__);
73#endif
74
75 if (prtlvl > PRINT_NONE) fasp_gettime(&AMG_start);
76
77 // check matrix data
79
80 // Step 0: initialize mgl[0] with A, b and x
81 mgl[0].A = fasp_dcsr_create(m, n, nnz);
82 fasp_dcsr_cp(A, &mgl[0].A);
83
84 mgl[0].b = fasp_dvec_create(n);
85 fasp_dvec_cp(b, &mgl[0].b);
86
87 mgl[0].x = fasp_dvec_create(n);
88 fasp_dvec_cp(x, &mgl[0].x);
89
90 // Step 1: AMG setup phase
91 switch (amg_type) {
92
93 case SA_AMG: // Smoothed Aggregation AMG setup
94 status = fasp_amg_setup_sa(mgl, param);
95 break;
96
97 case UA_AMG: // Unsmoothed Aggregation AMG setup
98 status = fasp_amg_setup_ua(mgl, param);
99 break;
100
101 default: // Classical AMG setup
102 status = fasp_amg_setup_rs(mgl, param);
103 break;
104 }
105
106 // Step 2: AMG solve phase
107 if (status == FASP_SUCCESS) { // call a multilevel cycle
108
109 switch (cycle_type) {
110
111 case AMLI_CYCLE: // AMLI-cycle
112 iter = fasp_amg_solve_amli(mgl, param);
113 break;
114
115 case NL_AMLI_CYCLE: // Nonlinear AMLI-cycle
116 iter = fasp_amg_solve_namli(mgl, param);
117 break;
118
119 default: // V,W-cycles or hybrid cycles (determined by param)
120 iter = fasp_amg_solve(mgl, param);
121 break;
122 }
123
124 fasp_dvec_cp(&mgl[0].x, x);
125
126 }
127
128 else { // call a backup solver
129
130 if (prtlvl > PRINT_MIN) {
131 printf("### WARNING: AMG setup failed!\n");
132 printf("### WARNING: Use a backup solver instead!\n");
133 }
134 fasp_solver_dcsr_spgmres(A, b, x, NULL, tol, maxit, 20, 1, prtlvl);
135 }
136
137 // clean-up memory
138 fasp_amg_data_free(mgl, param);
139
140 // print out CPU time if needed
141 if (prtlvl > PRINT_NONE) {
142 fasp_gettime(&AMG_end);
143 fasp_cputime("AMG totally", AMG_end - AMG_start);
144 }
145
146#if DEBUG_MODE > 0
147 printf("### DEBUG: [--End--] %s ...\n", __FUNCTION__);
148#endif
149
150 return iter;
151}
152
153/*---------------------------------*/
154/*-- End of File --*/
155/*---------------------------------*/
void fasp_cputime(const char *message, const REAL cputime)
Print CPU walltime.
Definition: AuxMessage.c:179
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
void fasp_dvec_cp(const dvector *x, dvector *y)
Copy dvector x to dvector y.
Definition: AuxVector.c:334
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_cp(const dCSRmat *A, dCSRmat *B)
copy a dCSRmat to a new one B=A
Definition: BlaSparseCSR.c:851
void fasp_check_dCSRmat(const dCSRmat *A)
Check whether an dCSRmat matrix is supported or not.
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
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
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_amg_solve_amli(AMG_data *mgl, AMG_param *param)
AMLI – SOLVE phase.
Definition: PreMGSolve.c:142
INT fasp_amg_solve_namli(AMG_data *mgl, AMG_param *param)
Nonlinear AMLI – SOLVE phase.
Definition: PreMGSolve.c:230
INT fasp_amg_solve(AMG_data *mgl, AMG_param *param)
AMG – SOLVE phase.
Definition: PreMGSolve.c:49
INT fasp_solver_amg(dCSRmat *A, dvector *b, dvector *x, AMG_param *param)
Solve Ax = b by algebraic multigrid methods.
Definition: SolAMG.c:49
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 FASP_SUCCESS
Definition of return status and error messages.
Definition: fasp_const.h:19
#define SA_AMG
Definition: fasp_const.h:164
#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 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 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
SHORT cycle_type
type of AMG cycle
Definition: fasp.h:476
SHORT max_levels
max number of levels of AMG
Definition: fasp.h:470
INT maxit
max number of iterations of AMG
Definition: fasp.h:464
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