Fast Auxiliary Space Preconditioning 2.7.7 Aug/28/2022
Loading...
Searching...
No Matches
PreMGCycle.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_functs.h"
22
23/*---------------------------------*/
24/*-- Declare Private Functions --*/
25/*---------------------------------*/
26
27#include "PreMGSmoother.inl"
28#include "PreMGUtil.inl"
29
30/*---------------------------------*/
31/*-- Public Functions --*/
32/*---------------------------------*/
33
49{
50 const SHORT prtlvl = param->print_level;
51 const SHORT amg_type = param->AMG_type;
52 const SHORT smoother = param->smoother;
53 const SHORT smooth_order = param->smooth_order;
54 const SHORT cycle_type = param->cycle_type;
55 const SHORT coarse_solver = param->coarse_solver;
56 const SHORT nl = mgl[0].num_levels;
57 const REAL relax = param->relaxation;
58 const REAL tol = param->tol * 1e-4;
59 const SHORT ndeg = param->polynomial_degree;
60
61 // Schwarz parameters
62 SWZ_param swzparam;
63 if (param->SWZ_levels > 0) {
64 swzparam.SWZ_blksolver = param->SWZ_blksolver;
65 }
66
67 // local variables
68 REAL alpha = 1.0;
69 INT num_lvl[MAX_AMG_LVL] = {0}, l = 0;
70
71 // more general cycling types on each level --zcs 05/07/2020
72 INT ncycles[MAX_AMG_LVL] = {1};
73 SHORT i;
74 for (i = 0; i < MAX_AMG_LVL; ++i) ncycles[i] = 1; // initially V-cycle
75 switch (cycle_type) {
76 case 12:
77 for (i = MAX_AMG_LVL - 2; i > 0; i -= 2) ncycles[i] = 2;
78 break;
79 case 21:
80 for (i = MAX_AMG_LVL - 1; i > 0; i -= 2) ncycles[i] = 2;
81 break;
82 default:
83 for (i = 0; i < MAX_AMG_LVL; i += 1) ncycles[i] = cycle_type;
84 }
85
86#if DEBUG_MODE > 0
87 printf("### DEBUG: [-Begin-] %s ...\n", __FUNCTION__);
88 printf("### DEBUG: n=%d, nnz=%d\n", mgl[0].A.row, mgl[0].A.nnz);
89#endif
90
91#if DEBUG_MODE > 1
92 printf("### DEBUG: AMG_level = %d, ILU_level = %d\n", nl, mgl->ILU_levels);
93#endif
94
95ForwardSweep:
96 while (l < nl - 1) {
97
98 num_lvl[l]++;
99
100 // pre-smoothing with ILU method
101 if (l < mgl->ILU_levels) {
102 fasp_smoother_dcsr_ilu(&mgl[l].A, &mgl[l].b, &mgl[l].x, &mgl[l].LU);
103 }
104
105 // or pre-smoothing with Schwarz method
106 else if (l < mgl->SWZ_levels) {
107 switch (mgl[l].Schwarz.SWZ_type) {
109 fasp_dcsr_swz_forward(&mgl[l].Schwarz, &swzparam, &mgl[l].x,
110 &mgl[l].b);
111 fasp_dcsr_swz_backward(&mgl[l].Schwarz, &swzparam, &mgl[l].x,
112 &mgl[l].b);
113 break;
114 default:
115 fasp_dcsr_swz_forward(&mgl[l].Schwarz, &swzparam, &mgl[l].x,
116 &mgl[l].b);
117 break;
118 }
119 }
120
121 // or pre-smoothing with standard smoother
122 else {
123#if MULTI_COLOR_ORDER
124 // printf("fasp_smoother_dcsr_gs_multicolor, %s, %d\n", __FUNCTION__,
125 // __LINE__);
126 fasp_smoother_dcsr_gs_multicolor(&mgl[l].x, &mgl[l].A, &mgl[l].b,
127 param->presmooth_iter, 1);
128#else
129 fasp_dcsr_presmoothing(smoother, &mgl[l].A, &mgl[l].b, &mgl[l].x,
130 param->presmooth_iter, 0, mgl[l].A.row - 1, 1, relax,
131 ndeg, smooth_order, mgl[l].cfmark.val);
132#endif
133 }
134
135 // form residual r = b - A x
136 fasp_darray_cp(mgl[l].A.row, mgl[l].b.val, mgl[l].w.val);
137 fasp_blas_dcsr_aAxpy(-1.0, &mgl[l].A, mgl[l].x.val, mgl[l].w.val);
138
139 // restriction r1 = R*r0
140 switch (amg_type) {
141 case UA_AMG:
142 fasp_blas_dcsr_mxv_agg(&mgl[l].R, mgl[l].w.val, mgl[l + 1].b.val);
143 break;
144 default:
145 fasp_blas_dcsr_mxv(&mgl[l].R, mgl[l].w.val, mgl[l + 1].b.val);
146 break;
147 }
148
149 // prepare for the next level
150 ++l;
151 fasp_dvec_set(mgl[l].A.row, &mgl[l].x, 0.0);
152 }
153
154 // If AMG only has one level or we have arrived at the coarsest level,
155 // call the coarse space solver:
156 switch (coarse_solver) {
157
158#if WITH_PARDISO
159 case SOLVER_PARDISO:
160 {
161 /* use Intel MKL PARDISO direct solver on the coarsest level */
162 fasp_pardiso_solve(&mgl[nl - 1].A, &mgl[nl - 1].b, &mgl[nl - 1].x,
163 &mgl[nl - 1].pdata, 0);
164 break;
165 }
166#endif
167
168#if WITH_MUMPS
169 case SOLVER_MUMPS:
170 {
171 // use MUMPS direct solver on the coarsest level
172 mgl[nl - 1].mumps.job = 2;
173 fasp_solver_mumps_steps(&mgl[nl - 1].A, &mgl[nl - 1].b, &mgl[nl - 1].x,
174 &mgl[nl - 1].mumps);
175 break;
176 }
177#endif
178
179#if WITH_UMFPACK
180 case SOLVER_UMFPACK:
181 {
182 // use UMFPACK direct solver on the coarsest level
183 fasp_umfpack_solve(&mgl[nl - 1].A, &mgl[nl - 1].b, &mgl[nl - 1].x,
184 mgl[nl - 1].Numeric, 0);
185 break;
186 }
187#endif
188
189#if WITH_SuperLU
190 case SOLVER_SUPERLU:
191 {
192 // use SuperLU direct solver on the coarsest level
193 fasp_solver_superlu(&mgl[nl - 1].A, &mgl[nl - 1].b, &mgl[nl - 1].x, 0);
194 break;
195 }
196#endif
197
198 default:
199 // use iterative solver on the coarsest level
200 fasp_coarse_itsolver(&mgl[nl - 1].A, &mgl[nl - 1].b, &mgl[nl - 1].x, tol,
201 prtlvl);
202 }
203
204 // BackwardSweep:
205 while (l > 0) {
206
207 --l;
208
209 // find the optimal scaling factor alpha
210 if (param->coarse_scaling == ON) {
211 alpha =
212 fasp_blas_darray_dotprod(mgl[l + 1].A.row, mgl[l + 1].x.val,
213 mgl[l + 1].b.val) /
214 fasp_blas_dcsr_vmv(&mgl[l + 1].A, mgl[l + 1].x.val, mgl[l + 1].x.val);
215 alpha = MIN(alpha, 1.0); // Add this for safety! --Chensong on 10/04/2014
216 }
217
218 // prolongation u = u + alpha*P*e1
219 switch (amg_type) {
220 case UA_AMG:
221 fasp_blas_dcsr_aAxpy_agg(alpha, &mgl[l].P, mgl[l + 1].x.val,
222 mgl[l].x.val);
223 break;
224 default:
225 fasp_blas_dcsr_aAxpy(alpha, &mgl[l].P, mgl[l + 1].x.val, mgl[l].x.val);
226 break;
227 }
228
229 // post-smoothing with ILU method
230 if (l < mgl->ILU_levels) {
231 fasp_smoother_dcsr_ilu(&mgl[l].A, &mgl[l].b, &mgl[l].x, &mgl[l].LU);
232 }
233
234 // post-smoothing with Schwarz method
235 else if (l < mgl->SWZ_levels) {
236 switch (mgl[l].Schwarz.SWZ_type) {
238 fasp_dcsr_swz_backward(&mgl[l].Schwarz, &swzparam, &mgl[l].x,
239 &mgl[l].b);
240 fasp_dcsr_swz_forward(&mgl[l].Schwarz, &swzparam, &mgl[l].x,
241 &mgl[l].b);
242 break;
243 default:
244 fasp_dcsr_swz_backward(&mgl[l].Schwarz, &swzparam, &mgl[l].x,
245 &mgl[l].b);
246 break;
247 }
248 }
249
250 // post-smoothing with standard methods
251 else {
252#if MULTI_COLOR_ORDER
253 fasp_smoother_dcsr_gs_multicolor(&mgl[l].x, &mgl[l].A, &mgl[l].b,
254 param->postsmooth_iter, -1);
255#else
256 fasp_dcsr_postsmoothing(smoother, &mgl[l].A, &mgl[l].b, &mgl[l].x,
257 param->postsmooth_iter, 0, mgl[l].A.row - 1, -1,
258 relax, ndeg, smooth_order, mgl[l].cfmark.val);
259#endif
260 }
261
262 // General cycling on each level --zcs
263 if (num_lvl[l] < ncycles[l])
264 break;
265 else
266 num_lvl[l] = 0;
267 }
268
269 if (l > 0) goto ForwardSweep;
270
271#if DEBUG_MODE > 0
272 printf("### DEBUG: [--End--] %s ...\n", __FUNCTION__);
273#endif
274}
275
288{
289 const SHORT prtlvl = param->print_level;
290 const SHORT nl = mgl[0].num_levels;
291 const SHORT smoother = param->smoother;
292 const SHORT cycle_type = param->cycle_type;
293 const SHORT coarse_solver = param->coarse_solver;
294 const REAL relax = param->relaxation;
295 INT steps = param->presmooth_iter;
296
297 // local variables
298 INT nu_l[MAX_AMG_LVL] = {0}, l = 0;
299 REAL alpha = 1.0;
300 INT i;
301
302 dvector r_nk, z_nk;
303
304 if (mgl[0].A_nk != NULL) {
305 fasp_dvec_alloc(mgl[0].A_nk->row, &r_nk);
306 fasp_dvec_alloc(mgl[0].A_nk->row, &z_nk);
307 }
308
309#if DEBUG_MODE > 0
310 printf("### DEBUG: [-Begin-] %s ...\n", __FUNCTION__);
311#endif
312
313#if DEBUG_MODE > 1
314 printf("### DEBUG: AMG_level = %d, ILU_level = %d\n", nl, mgl->ILU_levels);
315#endif
316
317ForwardSweep:
318 while (l < nl - 1) {
319 nu_l[l]++;
320 // pre smoothing
321 if (l < mgl->ILU_levels) {
322 fasp_smoother_dbsr_ilu(&mgl[l].A, &mgl[l].b, &mgl[l].x, &mgl[l].LU);
323 for (i = 0; i < steps; i++)
324 fasp_smoother_dbsr_gs_ascend(&mgl[l].A, &mgl[l].b, &mgl[l].x,
325 mgl[l].diaginv.val);
326 } else {
327 if (steps > 0) {
328 switch (smoother) {
329 case SMOOTHER_JACOBI:
330 for (i = 0; i < steps; i++)
331 fasp_smoother_dbsr_jacobi1(&mgl[l].A, &mgl[l].b, &mgl[l].x,
332 mgl[l].diaginv.val);
333 break;
334 case SMOOTHER_GS:
335 for (i = 0; i < steps; i++)
336 fasp_smoother_dbsr_gs_ascend(&mgl[l].A, &mgl[l].b,
337 &mgl[l].x, mgl[l].diaginv.val);
338 break;
339 case SMOOTHER_SGS:
340 for (i = 0; i < steps; i++) {
341 fasp_smoother_dbsr_gs_ascend(&mgl[l].A, &mgl[l].b,
342 &mgl[l].x, mgl[l].diaginv.val);
344 &mgl[l].A, &mgl[l].b, &mgl[l].x, mgl[l].diaginv.val);
345 }
346 break;
347 case SMOOTHER_SOR:
348 for (i = 0; i < steps; i++)
349 fasp_smoother_dbsr_sor_ascend(&mgl[l].A, &mgl[l].b,
350 &mgl[l].x, mgl[l].diaginv.val,
351 relax);
352 break;
353 case SMOOTHER_SSOR:
354 for (i = 0; i < steps; i++) {
355 fasp_smoother_dbsr_sor_ascend(&mgl[l].A, &mgl[l].b,
356 &mgl[l].x, mgl[l].diaginv.val,
357 relax);
358 }
359 fasp_smoother_dbsr_sor_descend(&mgl[l].A, &mgl[l].b, &mgl[l].x,
360 mgl[l].diaginv.val, relax);
361 break;
362 default:
363 printf("### ERROR: Unknown smoother type %d!\n", smoother);
364 fasp_chkerr(ERROR_SOLVER_TYPE, __FUNCTION__);
365 }
366 }
367 }
368
369 // extra kernel solve
370 if (mgl[l].A_nk != NULL) {
371
372 //--------------------------------------------
373 // extra kernel solve
374 //--------------------------------------------
375 // form residual r = b - A x
376 fasp_darray_cp(mgl[l].A.ROW * mgl[l].A.nb, mgl[l].b.val, mgl[l].w.val);
377 fasp_blas_dbsr_aAxpy(-1.0, &mgl[l].A, mgl[l].x.val, mgl[l].w.val);
378
379 // r_nk = R_nk*r
380 fasp_blas_dcsr_mxv(mgl[l].R_nk, mgl[l].w.val, r_nk.val);
381
382 // z_nk = A_nk^{-1}*r_nk
383#if WITH_UMFPACK // use UMFPACK directly
384 fasp_solver_umfpack(mgl[l].A_nk, &r_nk, &z_nk, 0);
385#else
386 fasp_coarse_itsolver(mgl[l].A_nk, &r_nk, &z_nk, 1e-12, 0);
387#endif
388
389 // z = z + P_nk*z_nk;
390 fasp_blas_dcsr_aAxpy(1.0, mgl[l].P_nk, z_nk.val, mgl[l].x.val);
391 }
392
393 // form residual r = b - A x
394 fasp_darray_cp(mgl[l].A.ROW * mgl[l].A.nb, mgl[l].b.val, mgl[l].w.val);
395 fasp_blas_dbsr_aAxpy(-1.0, &mgl[l].A, mgl[l].x.val, mgl[l].w.val);
396
397 // restriction r1 = R*r0
398 fasp_blas_dbsr_mxv(&mgl[l].R, mgl[l].w.val, mgl[l + 1].b.val);
399
400 // prepare for the next level
401 ++l;
402 fasp_dvec_set(mgl[l].A.ROW * mgl[l].A.nb, &mgl[l].x, 0.0);
403 }
404
405 // If AMG only has one level or we have arrived at the coarsest level,
406 // call the coarse space solver:
407 switch (coarse_solver) {
408
409#if WITH_PARDISO
410 case SOLVER_PARDISO:
411 {
412 /* use Intel MKL PARDISO direct solver on the coarsest level */
413 fasp_pardiso_solve(&mgl[nl - 1].Ac, &mgl[nl - 1].b, &mgl[nl - 1].x,
414 &mgl[nl - 1].pdata, 0);
415 break;
416 }
417#endif
418
419#if WITH_MUMPS
420 case SOLVER_MUMPS:
421 /* use MUMPS direct solver on the coarsest level */
422 mgl[nl - 1].mumps.job = 2;
423 fasp_solver_mumps_steps(&mgl[nl - 1].Ac, &mgl[nl - 1].b, &mgl[nl - 1].x,
424 &mgl[nl - 1].mumps);
425 break;
426#endif
427
428#if WITH_UMFPACK
429 case SOLVER_UMFPACK:
430 /* use UMFPACK direct solver on the coarsest level */
431 fasp_umfpack_solve(&mgl[nl - 1].Ac, &mgl[nl - 1].b, &mgl[nl - 1].x,
432 mgl[nl - 1].Numeric, 0);
433 break;
434#endif
435
436#if WITH_SuperLU
437 case SOLVER_SUPERLU:
438 /* use SuperLU direct solver on the coarsest level */
439 fasp_solver_superlu(&mgl[nl - 1].Ac, &mgl[nl - 1].b, &mgl[nl - 1].x, 0);
440 break;
441#endif
442
443 default:
444 {
445 /* use iterative solver on the coarsest level */
446 const INT csize = mgl[nl - 1].A.ROW * mgl[nl - 1].A.nb;
447 const INT cmaxit = MIN(csize * csize, 200);
448 const REAL ctol = param->tol; // coarse level tolerance
449 const REAL atol = ctol * 1e-8; // coarse level tolerance
450 if (fasp_solver_dbsr_pvgmres(&mgl[nl - 1].A, &mgl[nl - 1].b,
451 &mgl[nl - 1].x, NULL, ctol, atol, cmaxit,
452 25, 1, 0) < 0) {
453 if (prtlvl > PRINT_MIN) {
454 printf("### WARNING: Coarse level solver did not converge!\n");
455 printf("### WARNING: Consider to increase maxit to %d!\n",
456 2 * cmaxit);
457 }
458 }
459 }
460 }
461
462 // BackwardSweep:
463 while (l > 0) {
464 --l;
465
466 // prolongation u = u + alpha*P*e1
467 if (param->coarse_scaling == ON) {
468 dvector PeH, Aeh;
469 PeH.row = Aeh.row = mgl[l].b.row;
470 PeH.val = mgl[l].w.val + mgl[l].b.row;
471 Aeh.val = PeH.val + mgl[l].b.row;
472
473 fasp_blas_dbsr_mxv(&mgl[l].P, mgl[l + 1].x.val, PeH.val);
474 fasp_blas_dbsr_mxv(&mgl[l].A, PeH.val, Aeh.val);
475
476 alpha = (fasp_blas_darray_dotprod(mgl[l].b.row, Aeh.val, mgl[l].w.val)) /
477 (fasp_blas_darray_dotprod(mgl[l].b.row, Aeh.val, Aeh.val));
478 alpha = MIN(alpha, 1.0); // Add this for safety! --Chensong on 10/04/2014
479 fasp_blas_darray_axpy(mgl[l].b.row, alpha, PeH.val, mgl[l].x.val);
480 } else {
481 fasp_blas_dbsr_aAxpy(alpha, &mgl[l].P, mgl[l + 1].x.val, mgl[l].x.val);
482 }
483
484 // extra kernel solve
485 if (mgl[l].A_nk != NULL) {
486 //--------------------------------------------
487 // extra kernel solve
488 //--------------------------------------------
489 // form residual r = b - A x
490 fasp_darray_cp(mgl[l].A.ROW * mgl[l].A.nb, mgl[l].b.val, mgl[l].w.val);
491 fasp_blas_dbsr_aAxpy(-1.0, &mgl[l].A, mgl[l].x.val, mgl[l].w.val);
492
493 // r_nk = R_nk*r
494 fasp_blas_dcsr_mxv(mgl[l].R_nk, mgl[l].w.val, r_nk.val);
495
496 // z_nk = A_nk^{-1}*r_nk
497#if WITH_UMFPACK // use UMFPACK directly
498 fasp_solver_umfpack(mgl[l].A_nk, &r_nk, &z_nk, 0);
499#else
500 fasp_coarse_itsolver(mgl[l].A_nk, &r_nk, &z_nk, 1e-12, 0);
501#endif
502
503 // z = z + P_nk*z_nk;
504 fasp_blas_dcsr_aAxpy(1.0, mgl[l].P_nk, z_nk.val, mgl[l].x.val);
505 }
506
507 // post-smoothing
508 if (l < mgl->ILU_levels) {
509 fasp_smoother_dbsr_ilu(&mgl[l].A, &mgl[l].b, &mgl[l].x, &mgl[l].LU);
510 for (i = 0; i < steps; i++)
511 fasp_smoother_dbsr_gs_descend(&mgl[l].A, &mgl[l].b, &mgl[l].x,
512 mgl[l].diaginv.val);
513 } else {
514 if (steps > 0) {
515 switch (smoother) {
516 case SMOOTHER_JACOBI:
517 for (i = 0; i < steps; i++)
518 fasp_smoother_dbsr_jacobi1(&mgl[l].A, &mgl[l].b, &mgl[l].x,
519 mgl[l].diaginv.val);
520 break;
521 case SMOOTHER_GS:
522 for (i = 0; i < steps; i++)
524 &mgl[l].A, &mgl[l].b, &mgl[l].x, mgl[l].diaginv.val);
525 break;
526 case SMOOTHER_SGS:
527 for (i = 0; i < steps; i++) {
528 fasp_smoother_dbsr_gs_ascend(&mgl[l].A, &mgl[l].b,
529 &mgl[l].x, mgl[l].diaginv.val);
531 &mgl[l].A, &mgl[l].b, &mgl[l].x, mgl[l].diaginv.val);
532 }
533 break;
534 case SMOOTHER_SOR:
535 for (i = 0; i < steps; i++)
536 fasp_smoother_dbsr_sor_descend(&mgl[l].A, &mgl[l].b,
537 &mgl[l].x,
538 mgl[l].diaginv.val, relax);
539 break;
540 case SMOOTHER_SSOR:
541 for (i = 0; i < steps; i++)
542 fasp_smoother_dbsr_sor_ascend(&mgl[l].A, &mgl[l].b,
543 &mgl[l].x, mgl[l].diaginv.val,
544 relax);
545 fasp_smoother_dbsr_sor_descend(&mgl[l].A, &mgl[l].b, &mgl[l].x,
546 mgl[l].diaginv.val, relax);
547 break;
548 default:
549 printf("### ERROR: Unknown smoother type %d!\n", smoother);
550 fasp_chkerr(ERROR_SOLVER_TYPE, __FUNCTION__);
551 }
552 }
553 }
554
555 if (nu_l[l] < cycle_type)
556 break;
557 else
558 nu_l[l] = 0;
559 }
560
561 if (l > 0) goto ForwardSweep;
562
563#if DEBUG_MODE > 0
564 printf("### DEBUG: [--End--] %s ...\n", __FUNCTION__);
565#endif
566}
567
568/*---------------------------------*/
569/*-- End of File --*/
570/*---------------------------------*/
void fasp_darray_cp(const INT n, const REAL *x, REAL *y)
Copy an array to the other y=x.
Definition: AuxArray.c:210
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_dvec_set(INT n, dvector *x, const REAL val)
Initialize dvector x[i]=val for i=0:n-1.
Definition: AuxVector.c:222
void fasp_dvec_alloc(const INT m, dvector *u)
Create dvector data space of REAL type.
Definition: AuxVector.c:105
REAL fasp_blas_darray_dotprod(const INT n, const REAL *x, const REAL *y)
Inner product of two arraies x and y.
Definition: BlaArray.c:771
void fasp_blas_darray_axpy(const INT n, const REAL a, const REAL *x, REAL *y)
y = a*x + y
Definition: BlaArray.c:90
void fasp_dcsr_swz_forward(SWZ_data *swzdata, SWZ_param *swzparam, dvector *x, dvector *b)
Schwarz smoother: forward sweep.
void fasp_dcsr_swz_backward(SWZ_data *swzdata, SWZ_param *swzparam, dvector *x, dvector *b)
Schwarz smoother: backward sweep.
void fasp_blas_dbsr_aAxpy(const REAL alpha, const dBSRmat *A, const REAL *x, REAL *y)
Compute y := alpha*A*x + y.
Definition: BlaSpmvBSR.c:514
void fasp_blas_dbsr_mxv(const dBSRmat *A, const REAL *x, REAL *y)
Compute y := A*x.
Definition: BlaSpmvBSR.c:1055
void fasp_blas_dcsr_mxv_agg(const dCSRmat *A, const REAL *x, REAL *y)
Matrix-vector multiplication y = A*x (nonzeros of A = 1)
Definition: BlaSpmvCSR.c:438
void fasp_blas_dcsr_mxv(const dCSRmat *A, const REAL *x, REAL *y)
Matrix-vector multiplication y = A*x.
Definition: BlaSpmvCSR.c:242
void fasp_blas_dcsr_aAxpy_agg(const REAL alpha, const dCSRmat *A, const REAL *x, REAL *y)
Matrix-vector multiplication y = alpha*A*x + y (nonzeros of A = 1)
Definition: BlaSpmvCSR.c:727
REAL fasp_blas_dcsr_vmv(const dCSRmat *A, const REAL *x, const REAL *y)
vector-Matrix-vector multiplication alpha = y'*A*x
Definition: BlaSpmvCSR.c:839
void fasp_blas_dcsr_aAxpy(const REAL alpha, const dCSRmat *A, const REAL *x, REAL *y)
Matrix-vector multiplication y = alpha*A*x + y.
Definition: BlaSpmvCSR.c:494
void fasp_smoother_dbsr_gs_descend(dBSRmat *A, dvector *b, dvector *u, REAL *diaginv)
Gauss-Seidel relaxation in the descending order.
void fasp_smoother_dbsr_jacobi1(dBSRmat *A, dvector *b, dvector *u, REAL *diaginv)
Jacobi relaxation.
void fasp_smoother_dbsr_gs_ascend(dBSRmat *A, dvector *b, dvector *u, REAL *diaginv)
Gauss-Seidel relaxation in the ascending order.
void fasp_smoother_dbsr_sor_ascend(dBSRmat *A, dvector *b, dvector *u, REAL *diaginv, REAL weight)
SOR relaxation in the ascending order.
void fasp_smoother_dbsr_ilu(dBSRmat *A, dvector *b, dvector *x, void *data)
ILU method as the smoother in solving Au=b with multigrid method.
void fasp_smoother_dbsr_sor_descend(dBSRmat *A, dvector *b, dvector *u, REAL *diaginv, REAL weight)
SOR relaxation in the descending order.
void fasp_smoother_dcsr_ilu(dCSRmat *A, dvector *b, dvector *x, void *data)
ILU method as a smoother.
INT fasp_solver_dbsr_pvgmres(dBSRmat *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:416
void fasp_solver_mgcycle_bsr(AMG_data_bsr *mgl, AMG_param *param)
Solve Ax=b with non-recursive multigrid cycle.
Definition: PreMGCycle.c:287
void fasp_solver_mgcycle(AMG_data *mgl, AMG_param *param)
Solve Ax=b with non-recursive multigrid cycle.
Definition: PreMGCycle.c:48
int fasp_solver_mumps_steps(dCSRmat *ptrA, dvector *b, dvector *u, Mumps_data *mumps)
Solve Ax=b by MUMPS in three steps.
Definition: XtrMumps.c:196
int fasp_solver_superlu(dCSRmat *ptrA, dvector *b, dvector *u, const SHORT prtlvl)
Solve Au=b by SuperLU.
Definition: XtrSuperlu.c:47
INT fasp_solver_umfpack(dCSRmat *ptrA, dvector *b, dvector *u, const SHORT prtlvl)
Solve Au=b by UMFpack.
Definition: XtrUmfpack.c:44
Main header file for the FASP project.
#define MIN(a, b)
Definition: fasp.h:83
#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 SMOOTHER_SOR
Definition: fasp_const.h:195
#define SOLVER_PARDISO
Definition: fasp_const.h:126
#define MAX_AMG_LVL
Definition: fasp_const.h:259
#define SMOOTHER_JACOBI
Definition: fasp_const.h:190
#define SOLVER_MUMPS
Definition: fasp_const.h:125
#define ERROR_SOLVER_TYPE
Definition: fasp_const.h:41
#define SOLVER_SUPERLU
Definition: fasp_const.h:123
#define SCHWARZ_SYMMETRIC
Definition: fasp_const.h:158
#define SMOOTHER_GS
Definition: fasp_const.h:191
#define SOLVER_UMFPACK
Definition: fasp_const.h:124
#define SMOOTHER_SGS
Definition: fasp_const.h:193
#define ON
Definition of switch.
Definition: fasp_const.h:67
#define SMOOTHER_SSOR
Definition: fasp_const.h:196
#define PRINT_MIN
Definition: fasp_const.h:74
#define UA_AMG
Definition: fasp_const.h:165
Data for multigrid levels in dBSRmat format.
Definition: fasp_block.h:146
dBSRmat A
pointer to the matrix at level level_num
Definition: fasp_block.h:155
Mumps_data mumps
data for MUMPS
Definition: fasp_block.h:245
dvector b
pointer to the right-hand side at level level_num
Definition: fasp_block.h:164
INT ILU_levels
number of levels use ILU smoother
Definition: fasp_block.h:217
INT num_levels
number of levels in use <= max_levels
Definition: fasp_block.h:152
dvector x
pointer to the iterative solution at level level_num
Definition: fasp_block.h:167
dvector w
temporary work space
Definition: fasp_block.h:242
Data for AMG methods.
Definition: fasp.h:804
dCSRmat A
pointer to the matrix at level level_num
Definition: fasp.h:817
Mumps_data mumps
data for MUMPS
Definition: fasp.h:866
dvector b
pointer to the right-hand side at level level_num
Definition: fasp.h:826
INT ILU_levels
number of levels use ILU smoother
Definition: fasp.h:843
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
dvector w
temporary work space
Definition: fasp.h:863
ivector cfmark
pointer to the CF marker at level level_num
Definition: fasp.h:840
Parameters for AMG methods.
Definition: fasp.h:455
SHORT print_level
print level for AMG
Definition: fasp.h:461
SHORT polynomial_degree
degree of the polynomial smoother
Definition: fasp.h:497
SHORT coarse_scaling
switch of scaling of the coarse grid correction
Definition: fasp.h:503
SHORT AMG_type
type of AMG method
Definition: fasp.h:458
REAL tol
stopping tolerance for AMG solver
Definition: fasp.h:467
SHORT coarse_solver
coarse solver type
Definition: fasp.h:500
REAL relaxation
relaxation parameter for Jacobi and SOR smoother
Definition: fasp.h:494
SHORT smoother
smoother type
Definition: fasp.h:482
SHORT cycle_type
type of AMG cycle
Definition: fasp.h:476
INT SWZ_blksolver
type of Schwarz block solver
Definition: fasp.h:590
SHORT postsmooth_iter
number of postsmoothers
Definition: fasp.h:491
INT SWZ_levels
number of levels use Schwarz smoother
Definition: fasp.h:578
SHORT presmooth_iter
number of presmoothers
Definition: fasp.h:488
SHORT smooth_order
smoother order
Definition: fasp.h:485
INT job
work for MUMPS
Definition: fasp.h:615
Parameters for Schwarz method.
Definition: fasp.h:430
INT SWZ_blksolver
type of Schwarz block solver
Definition: fasp.h:445
INT nb
dimension of each sub-block
Definition: fasp_block.h:46
INT ROW
number of rows of sub-blocks in matrix A, M
Definition: fasp_block.h:37
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
REAL * val
actual vector entries
Definition: fasp.h:360
INT row
number of rows
Definition: fasp.h:357
INT * val
actual vector entries
Definition: fasp.h:374