OpenCAEPoro  0.2.0 Sep/22/2022
A simulator for multicomponent porous media flow
FaspSolver.cpp
Go to the documentation of this file.
1 
12 #include "FaspSolver.hpp"
13 
14 void FaspSolver::SetupParam(const string& dir, const string& file)
15 {
16  solveDir = dir;
17  solveFile = file;
18 
19  string myfile = solveDir + solveFile;
20  InitParam(); // Set default solver parameters
21  ifstream ifs(myfile);
22  if (!ifs.is_open()) {
23  cout << "The input file " << myfile << " is missing!" << endl;
24  myfile = solveDir + "../conf/csr.fasp";
25  ifs.open(myfile);
26  if (!ifs.is_open()) {
27  cout << "The input file " << myfile << " is missing!" << endl;
28  cout << "Using the default parameters of FASP" << endl;
29  } else {
30  ifs.close();
31  cout << "Using the input file " << myfile << endl;
32  fasp_param_input(myfile.data(), &inParam);
33  }
34  } else {
35  ifs.close(); // if file has been opened, close it first
36  fasp_param_input(myfile.data(), &inParam);
37  }
38  fasp_param_init(&inParam, &itParam, &amgParam, &iluParam, &swzParam);
39 }
40 
41 void ScalarFaspSolver::Allocate(const vector<USI>& rowCapacity, const OCP_USI& maxDim,
42  const USI& blockDim)
43 {
44  OCP_USI nnz = 0;
45  for (OCP_USI n = 0; n < maxDim; n++) {
46  nnz += rowCapacity[n];
47  }
48  A = fasp_dcsr_create(maxDim, maxDim, nnz);
49 }
50 
51 void ScalarFaspSolver::InitParam()
52 {
53  // Input/output
54  inParam.print_level = PRINT_MIN;
55  inParam.output_type = 0;
56 
57  // Problem information
58  inParam.solver_type = SOLVER_VFGMRES;
59  inParam.precond_type = PREC_AMG;
60  inParam.stop_type = STOP_REL_RES;
61 
62  // LinearSolver parameters
63  inParam.itsolver_tol = 1e-4;
64  inParam.itsolver_maxit = 100;
65  inParam.restart = 30;
66 
67  // ILU method parameters
68  inParam.ILU_type = ILUk;
69  inParam.ILU_lfil = 0;
70  inParam.ILU_droptol = 0.001;
71  inParam.ILU_relax = 0;
72  inParam.ILU_permtol = 0.0;
73 
74  // Schwarz method parameters
75  inParam.SWZ_mmsize = 200;
76  inParam.SWZ_maxlvl = 2;
77  inParam.SWZ_type = 1;
78  inParam.SWZ_blksolver = SOLVER_DEFAULT;
79 
80  // AMG method parameters
81  inParam.AMG_type = CLASSIC_AMG;
82  inParam.AMG_levels = 20;
83  inParam.AMG_cycle_type = V_CYCLE;
84  inParam.AMG_smoother = SMOOTHER_GS;
85  inParam.AMG_smooth_order = CF_ORDER;
86  inParam.AMG_presmooth_iter = 1;
87  inParam.AMG_postsmooth_iter = 1;
88  inParam.AMG_relaxation = 1.0;
89  inParam.AMG_coarse_dof = 500;
90  inParam.AMG_coarse_solver = 0;
91  inParam.AMG_tol = 1e-6;
92  inParam.AMG_maxit = 1;
93  inParam.AMG_ILU_levels = 0;
94  inParam.AMG_SWZ_levels = 0;
95  inParam.AMG_coarse_scaling = OFF;
96  inParam.AMG_amli_degree = 1;
97  inParam.AMG_nl_amli_krylov_type = 2;
98 
99  // Classical AMG specific
100  inParam.AMG_coarsening_type = 1;
101  inParam.AMG_interpolation_type = 1;
102  inParam.AMG_max_row_sum = 0.9;
103  inParam.AMG_strong_threshold = 0.3;
104  inParam.AMG_truncation_threshold = 0.2;
105  inParam.AMG_aggressive_level = 0;
106  inParam.AMG_aggressive_path = 1;
107 
108  // Aggregation AMG specific
109  inParam.AMG_aggregation_type = PAIRWISE;
110  inParam.AMG_quality_bound = 8.0;
111  inParam.AMG_pair_number = 2;
112  inParam.AMG_strong_coupled = 0.25;
113  inParam.AMG_max_aggregation = 9;
114  inParam.AMG_tentative_smooth = 0.67;
115  inParam.AMG_smooth_filter = ON;
116  inParam.AMG_smooth_restriction = ON;
117 }
118 
119 void ScalarFaspSolver::AssembleMat(const vector<vector<USI>>& colId,
120  const vector<vector<OCP_DBL>>& val,
121  const OCP_USI& dim, const USI& blockDim,
122  vector<OCP_DBL>& rhs, vector<OCP_DBL>& u)
123 {
124  // b & x
125  b.row = dim;
126  b.val = rhs.data();
127  x.row = dim;
128  x.val = u.data();
129  // A
130  OCP_USI nnz = 0;
131  for (OCP_USI i = 0; i < dim; i++) {
132  nnz += colId[i].size();
133  }
134 
135  A.row = dim;
136  A.col = dim;
137  A.nnz = nnz;
138 
139  // IA
140  OCP_USI count = 0;
141  A.IA[0] = 0;
142  for (OCP_USI i = 1; i < dim + 1; i++) {
143  USI nnz_Row = colId[i - 1].size();
144  A.IA[i] = A.IA[i - 1] + nnz_Row;
145 
146  for (USI j = 0; j < nnz_Row; j++) {
147  A.JA[count] = colId[i - 1][j];
148  A.val[count] = val[i - 1][j];
149  count++;
150  }
151  }
152 }
153 
154 OCP_INT ScalarFaspSolver::Solve()
155 {
156  OCP_INT status = FASP_SUCCESS;
157 
158  const OCP_INT print_level = inParam.print_level;
159  const OCP_INT solver_type = inParam.solver_type;
160  const OCP_INT precond_type = inParam.precond_type;
161  const OCP_INT output_type = inParam.output_type;
162 
163  if (output_type) {
164  const char* outputfile = "out/Solver.out";
165  printf("Redirecting outputs to file: %s ...\n", outputfile);
166  freopen(outputfile, "w", stdout); // open a file for stdout
167  }
168 
169  // Preconditioned Krylov methods
170  if (solver_type >= 1 && solver_type <= 20) {
171 
172  // Using no preconditioner for Krylov iterative methods
173  if (precond_type == PREC_NULL) {
174  status = fasp_solver_dcsr_krylov(&A, &b, &x, &itParam);
175  }
176 
177  // Using diag(A) as preconditioner for Krylov iterative methods
178  else if (precond_type == PREC_DIAG) {
179  status = fasp_solver_dcsr_krylov_diag(&A, &b, &x, &itParam);
180  }
181 
182  // Using AMG as preconditioner for Krylov iterative methods
183  else if (precond_type == PREC_AMG || precond_type == PREC_FMG) {
184  if (print_level > PRINT_NONE) fasp_param_amg_print(&amgParam);
185  status = fasp_solver_dcsr_krylov_amg(&A, &b, &x, &itParam, &amgParam);
186  }
187 
188  // Using ILU as preconditioner for Krylov iterative methods
189  else if (precond_type == PREC_ILU) {
190  if (print_level > PRINT_NONE) fasp_param_ilu_print(&iluParam);
191  status = fasp_solver_dcsr_krylov_ilu(&A, &b, &x, &itParam, &iluParam);
192  }
193 
194  // Undefined iterative methods
195  else {
196  printf("### ERROR: Wrong preconditioner type %d!!!\n", precond_type);
197  status = ERROR_SOLVER_PRECTYPE;
198  }
199 
200  }
201 
202  // AMG as the iterative solver
203  else if (solver_type == SOLVER_AMG) {
204  if (print_level > PRINT_NONE) fasp_param_amg_print(&amgParam);
205  fasp_solver_amg(&A, &b, &x, &amgParam);
206  }
207 
208  // Full AMG as the iterative solver
209  else if (solver_type == SOLVER_FMG) {
210  if (print_level > PRINT_NONE) fasp_param_amg_print(&amgParam);
211  fasp_solver_famg(&A, &b, &x, &amgParam);
212  }
213 
214 #if WITH_MUMPS // use MUMPS directly
215  else if (solver_type == SOLVER_MUMPS) {
216  status = fasp_solver_mumps(&A, &b, &x, print_level);
217  if (status >= 0) status = 1; // Direct solver returns 1
218  }
219 #endif
220 
221 #if WITH_SuperLU // use SuperLU directly
222  else if (solver_type == SOLVER_SUPERLU) {
223  status = fasp_solver_superlu(&A, &b, &x, print_level);
224  if (status >= 0) status = 1; // Direct solver returns 1
225  }
226 #endif
227 
228 #if WITH_UMFPACK // use UMFPACK directly
229  else if (solver_type == SOLVER_UMFPACK) {
230  // Need to sort the matrix A for UMFPACK to work
231  dCSRmat A_trans = fasp_dcsr_create(A.row, A.col, A.nnz);
232  fasp_dcsr_transz(&A, NULL, &A_trans);
233  fasp_dcsr_sort(&A_trans);
234  status = fasp_solver_umfpack(&A_trans, &b, &x, print_level);
235  fasp_dcsr_free(&A_trans);
236  if (status >= 0) status = 1; // Direct solver returns 1
237  }
238 #endif
239 
240 #ifdef WITH_PARDISO // use PARDISO directly
241  else if (solver_type == SOLVER_PARDISO) {
242  fasp_dcsr_sort(&A);
243  status = fasp_solver_pardiso(&A, &b, &x, print_level);
244  if (status >= 0) status = 1; // Direct solver returns 1
245  }
246 #endif
247 
248  else {
249  printf("### ERROR: Wrong solver type %d!!!\n", solver_type);
250  status = ERROR_SOLVER_TYPE;
251  }
252 
253  if (status < 0) {
254  printf("### ERROR: Solver failed! Exit status = %d.\n\n", status);
255  }
256 
257  if (output_type) fclose(stdout);
258  return status;
259 }
260 
261 void VectorFaspSolver::Allocate(const vector<USI>& rowCapacity, const OCP_USI& maxDim,
262  const USI& blockDim)
263 {
264  OCP_USI nnz = 0;
265  for (OCP_USI n = 0; n < maxDim; n++) {
266  nnz += rowCapacity[n];
267  }
268  A = fasp_dbsr_create(maxDim, maxDim, nnz, blockDim, 0);
269  Asc = fasp_dbsr_create(maxDim, maxDim, nnz, blockDim, 0);
270  fsc = fasp_dvec_create(maxDim * blockDim);
271  order = fasp_ivec_create(maxDim);
272  Dmat.resize(maxDim * blockDim * blockDim);
273 }
274 
275 void VectorFaspSolver::InitParam()
276 {
277  // Input/output
278  inParam.print_level = PRINT_MIN;
279  inParam.output_type = 0;
280 
281  // Problem information
282  inParam.solver_type = SOLVER_VFGMRES;
283  inParam.decoup_type = 1;
284  inParam.precond_type = 64;
285  inParam.stop_type = STOP_REL_RES;
286 
287  // Solver parameters
288  inParam.itsolver_tol = 1e-3;
289  inParam.itsolver_maxit = 100;
290  inParam.restart = 30;
291 
292  // ILU method parameters
293  inParam.ILU_type = ILUk;
294  inParam.ILU_lfil = 0;
295  inParam.ILU_droptol = 0.001;
296  inParam.ILU_relax = 0;
297  inParam.ILU_permtol = 0.0;
298 
299  // Schwarz method parameters
300  inParam.SWZ_mmsize = 200;
301  inParam.SWZ_maxlvl = 2;
302  inParam.SWZ_type = 1;
303  inParam.SWZ_blksolver = SOLVER_DEFAULT;
304 
305  // AMG method parameters
306  inParam.AMG_type = CLASSIC_AMG;
307  inParam.AMG_levels = 20;
308  inParam.AMG_cycle_type = V_CYCLE;
309  inParam.AMG_smoother = SMOOTHER_GS;
310  inParam.AMG_smooth_order = CF_ORDER;
311  inParam.AMG_presmooth_iter = 1;
312  inParam.AMG_postsmooth_iter = 1;
313  inParam.AMG_relaxation = 1.0;
314  inParam.AMG_coarse_dof = 500;
315  inParam.AMG_coarse_solver = 0;
316  inParam.AMG_tol = 1e-6;
317  inParam.AMG_maxit = 1;
318  inParam.AMG_ILU_levels = 0;
319  inParam.AMG_SWZ_levels = 0;
320  inParam.AMG_coarse_scaling = OFF; // Require investigation --Chensong
321  inParam.AMG_amli_degree = 1;
322  inParam.AMG_nl_amli_krylov_type = 2;
323 
324  // Classical AMG specific
325  inParam.AMG_coarsening_type = 1;
326  inParam.AMG_interpolation_type = 1;
327  inParam.AMG_max_row_sum = 0.9;
328  inParam.AMG_strong_threshold = 0.3;
329  inParam.AMG_truncation_threshold = 0.2;
330  inParam.AMG_aggressive_level = 0;
331  inParam.AMG_aggressive_path = 1;
332 
333  // Aggregation AMG specific
334  inParam.AMG_aggregation_type = PAIRWISE;
335  inParam.AMG_quality_bound = 8.0;
336  inParam.AMG_pair_number = 2;
337  inParam.AMG_strong_coupled = 0.25;
338  inParam.AMG_max_aggregation = 9;
339  inParam.AMG_tentative_smooth = 0.67;
340  inParam.AMG_smooth_filter = ON;
341  inParam.AMG_smooth_restriction = ON;
342 }
343 
344 void VectorFaspSolver::AssembleMat(const vector<vector<USI>>& colId,
345  const vector<vector<OCP_DBL>>& val,
346  const OCP_USI& dim, const USI& blockDim,
347  vector<OCP_DBL>& rhs, vector<OCP_DBL>& u)
348 {
349  OCP_USI nrow = dim * blockDim;
350  // b & x
351  b.row = nrow;
352  b.val = rhs.data();
353  x.row = nrow;
354  x.val = u.data();
355 
356  // fsc & order
357  fsc.row = nrow;
358  order.row = nrow;
359 
360  // nnz
361  OCP_USI nnz = 0;
362  for (OCP_USI i = 0; i < dim; i++) {
363  nnz += colId[i].size();
364  }
365 
366  // Asc
367  Asc.ROW = dim;
368  Asc.COL = dim;
369  Asc.nb = blockDim;
370  Asc.NNZ = nnz;
371 
372  // A
373  A.ROW = dim;
374  A.COL = dim;
375  A.nb = blockDim;
376  A.NNZ = nnz;
377 
378  OCP_USI count = 0;
379  OCP_USI count2 = 0;
380  OCP_USI size_row;
381  A.IA[0] = 0;
382  for (OCP_USI i = 1; i < dim + 1; i++) {
383  USI nnb_Row = colId[i - 1].size();
384  A.IA[i] = A.IA[i - 1] + nnb_Row;
385 
386  for (USI j = 0; j < nnb_Row; j++) {
387  A.JA[count] = colId[i - 1][j];
388  count++;
389  }
390  size_row = nnb_Row * blockDim * blockDim;
391  for (USI k = 0; k < size_row; k++) {
392  A.val[count2] = val[i - 1][k];
393  count2++;
394  }
395  }
396 
397 #ifdef _DEBUG
398  // check x and b ---- for test
399  for (OCP_USI i = 0; i < nrow; i++) {
400  if (!isfinite(b.val[i]))
401  OCP_ABORT("vFasp b is infinite!");
402  if (!isfinite(x.val[i]))
403  OCP_ABORT("vFasp x is infinite!");
404  }
405  // check A ---- for test
406  for (OCP_USI i = 0; i < A.NNZ; i++) {
407  if (!isfinite(A.val[i]))
408  OCP_ABORT("vFasp A is infinite!");
409  }
410 #endif // DEBUG
411 
412 }
413 
414 OCP_INT VectorFaspSolver::Solve()
415 {
416  OCP_INT status = FASP_SUCCESS;
417 
418  // Set local parameters
419  const OCP_INT print_level = inParam.print_level;
420  const OCP_INT solver_type = inParam.solver_type;
421  const OCP_INT decoup_type = inParam.decoup_type;
422  const OCP_INT precond_type = inParam.precond_type;
423  const OCP_INT output_type = inParam.output_type;
424 
425  if (output_type) {
426  const char* outputfile = "../output/test.out";
427  printf("Redirecting outputs to file: %s ...\n", outputfile);
428  freopen(outputfile, "w", stdout); // open a file for stdout
429  }
430 
431  fasp_dvec_set(x.row, &x, 0);
432 
433  // Preconditioned Krylov methods
434  if (solver_type >= 1 && solver_type <= 10) {
435 
436  // Preconditioned Krylov methods in BSR format
437  switch (precond_type) {
438  case PC_NULL:
439  status = fasp_solver_dbsr_krylov(&A, &b, &x, &itParam);
440  break;
441  case PC_DIAG:
442  status = fasp_solver_dbsr_krylov_diag(&A, &b, &x, &itParam);
443  break;
444  case PC_BILU:
445  status = fasp_solver_dbsr_krylov_ilu(&A, &b, &x, &itParam, &iluParam);
446  break;
447 #if WITH_FASP4BLKOIL
448  case PC_FASP1:
449  Decoupling(&A, &b, &Asc, &fsc, &order, Dmat.data(), decoup_type);
450 #if WITH_FASP4CUDA // zhaoli 2022.04.04
451  status = fasp_solver_dbsr_krylov_FASP1_cuda_interface(
452  &Asc, &fsc, &x, &itParam, &iluParam, &amgParam, NULL, &order);
453 #else
454  status = fasp_solver_dbsr_krylov_FASP1a(
455  &Asc, &fsc, &x, &itParam, &iluParam, &amgParam, NULL, &order);
456 #endif
457  break;
458  case PC_FASP1_SHARE: // zhaoli 2021.03.24
459  Decoupling(&A, &b, &Asc, &fsc, &order, Dmat.data(), decoup_type);
460 #if WITH_FASP4CUDA
461  status = fasp_solver_dbsr_krylov_FASP1_cuda_share_interface(
462  &Asc, &fsc, &x, &itParam, &iluParam, &amgParam, NULL, &order, RESET_CONST);
463 #else
464  status = fasp_solver_dbsr_krylov_FASP1a_share_interface(
465  &Asc, &fsc, &x, &itParam, &iluParam, &amgParam, NULL, &order, RESET_CONST);
466 #endif
467  break;
468  case PC_FASP2:
469  Decoupling(&A, &b, &Asc, &fsc, &order, Dmat.data(), decoup_type);
470  status = fasp_solver_dbsr_krylov_FASP2(
471  &Asc, &fsc, &x, &itParam, &iluParam, &amgParam, NULL, &order);
472  break;
473  case PC_FASP3:
474  Decoupling(&A, &b, &Asc, &fsc, &order, Dmat.data(), decoup_type);
475  status = fasp_solver_dbsr_krylov_FASP3(
476  &Asc, &fsc, &x, &itParam, &iluParam, &amgParam, NULL, &order);
477  break;
478  case PC_FASP4:
479  Decoupling(&A, &b, &Asc, &fsc, &order, Dmat.data(), decoup_type);
480 #if WITH_FASP4CUDA
481  status = fasp_solver_dbsr_krylov_FASP4_cuda(
482  &Asc, &fsc, &x, &itParam, &iluParam, &amgParam, NULL, &order);
483 #else
484  status = fasp_solver_dbsr_krylov_FASP4(
485  &Asc, &fsc, &x, &itParam, &iluParam, &amgParam, NULL, &order);
486 #endif
487  break;
488  case PC_FASP4_SHARE: // zhaoli 2021.04.24
489  Decoupling(&A, &b, &Asc, &fsc, &order, Dmat.data(), decoup_type);
490 #if WITH_FASP4CUDA // zhaoli 2022.08.03
491  status = fasp_solver_dbsr_krylov_FASP4_cuda_share_interface(
492  &Asc, &fsc, &x, &itParam, &iluParam, &amgParam, NULL, &order,
493  RESET_CONST);
494 #else
495  status = fasp_solver_dbsr_krylov_FASP4_share_interface(
496  &Asc, &fsc, &x, &itParam, &iluParam, &amgParam, NULL, &order,
497  RESET_CONST);
498 #endif
499  break;
500  case PC_FASP5:
501  Decoupling(&A, &b, &Asc, &fsc, &order, Dmat.data(), decoup_type);
502  status = fasp_solver_dbsr_krylov_FASP5(
503  &Asc, &fsc, &x, &itParam, &iluParam, &amgParam, NULL, &order);
504  break;
505 #endif
506  default:
507  OCP_ABORT("Preconditioner type " + to_string(precond_type) +
508  " not supported!");
509  }
510  fill(Dmat.begin(), Dmat.end(), 0.0);
511  }
512 
513 #if WITH_MUMPS // use MUMPS directly
514  else if (solver_type == SOLVER_MUMPS) {
515  dCSRmat Acsr = fasp_format_dbsr_dcsr(&A);
516  status = fasp_solver_mumps(&Acsr, &b, &x, print_level);
517  fasp_dcsr_free(&Acsr);
518  if (status >= 0) status = 1; // Direct solver returns 1
519  }
520 #endif
521 
522 #if WITH_SuperLU // use SuperLU directly
523  else if (solver_type == SOLVER_SUPERLU) {
524  dCSRmat Acsr = fasp_format_dbsr_dcsr(&A);
525  status = fasp_solver_superlu(&Acsr, &b, &x, print_level);
526  fasp_dcsr_free(&Acsr);
527  if (status >= 0) status = 1; // Direct solver returns 1
528  }
529 #endif
530 
531 #if WITH_UMFPACK // use UMFPACK directly
532  else if (solver_type == SOLVER_UMFPACK) {
533  // Need to sort the matrix A for UMFPACK to work
534  dCSRmat Acsr = fasp_format_dbsr_dcsr(&A);
535  dCSRmat A_trans = fasp_dcsr_create(Acsr.row, Acsr.col, Acsr.nnz);
536  fasp_dcsr_transz(&Acsr, NULL, &A_trans);
537  fasp_dcsr_sort(&A_trans);
538  status = fasp_solver_umfpack(&A_trans, &b, &x, print_level);
539  fasp_dcsr_free(&A_trans);
540  fasp_dcsr_free(&Acsr);
541  if (status >= 0) status = 1; // Direct solver returns 1
542  }
543 #endif
544 
545 #ifdef WITH_PARDISO // use PARDISO directly
546  else if (solver_type == SOLVER_PARDISO) {
547  dCSRmat Acsr = fasp_format_dbsr_dcsr(&A);
548  fasp_dcsr_sort(&Acsr);
549  status = fasp_solver_pardiso(&Acsr, &b, &x, print_level);
550  fasp_dcsr_free(&Acsr);
551  if (status >= 0) status = 1; // Direct solver returns 1
552  }
553 #endif
554 
555  else {
556  printf("### ERROR: Wrong solver type %d!!!\n", solver_type);
557  status = ERROR_SOLVER_TYPE;
558  }
559 
560  if (print_level > PRINT_MIN) {
561  if (status < 0) {
562  cout << "\n### WARNING: Solver does not converge!\n" << endl;
563  } else {
564  cout << "\nSolver converges successfully!\n" << endl;
565  }
566  }
567 
568  if (output_type) fclose(stdout);
569 
570  return status;
571 }
572 
573 /*----------------------------------------------------------------------------*/
574 /* Brief Change History of This File */
575 /*----------------------------------------------------------------------------*/
576 /* Author Date Actions */
577 /*----------------------------------------------------------------------------*/
578 /* Shizhe Li Nov/22/2021 Create file */
579 /* Chensong Zhang Jan/19/2022 Set FASP4BLKOIL as optional */
580 /* Li Zhao Apr/04/2022 Set FASP4CUDA as optional */
581 /*----------------------------------------------------------------------------*/
Declaration of classes interfacing to the FASP solvers.
#define PC_BILU
BILU: block ILU preconditioner.
Definition: FaspSolver.hpp:56
#define PC_FASP2
FASP2: MSP, experimental only.
Definition: FaspSolver.hpp:51
#define PC_FASP3
FASP3: MSP, monolithic preconditioner.
Definition: FaspSolver.hpp:52
#define PC_FASP5
FASP5: MSP, experimental only.
Definition: FaspSolver.hpp:54
#define RESET_CONST
Sharing threshold for PC_FASP1_SHARE, PC_FASP4_SHARE.
Definition: FaspSolver.hpp:61
#define PC_NULL
None: no preconditioner.
Definition: FaspSolver.hpp:49
#define PC_FASP4
FASP4: MSP, default for FIM from 2015.
Definition: FaspSolver.hpp:53
#define PC_DIAG
DIAG: diagonal preconditioner.
Definition: FaspSolver.hpp:55
#define PC_FASP1_SHARE
Sharing setup stage for PC_FASP1, use with caution.
Definition: FaspSolver.hpp:59
#define PC_FASP4_SHARE
Sharing setup stage for PC_FASP4, use with caution.
Definition: FaspSolver.hpp:60
#define PC_FASP1
FASP1: MSP, default for FIM from 2020.
Definition: FaspSolver.hpp:50
unsigned int USI
Generic unsigned integer.
Definition: OCPConst.hpp:22
unsigned int OCP_USI
Long unsigned integer.
Definition: OCPConst.hpp:24
int OCP_INT
Long integer.
Definition: OCPConst.hpp:25
#define OCP_ABORT(msg)
Abort if critical error happens.
Definition: UtilError.hpp:47
AMG_param amgParam
Parameters for AMG method.
Definition: FaspSolver.hpp:78
input_param inParam
Parameters from input files.
Definition: FaspSolver.hpp:76
ILU_param iluParam
Parameters for ILU method.
Definition: FaspSolver.hpp:79
string solveFile
Relative path of fasp file.
Definition: FaspSolver.hpp:75
SWZ_param swzParam
Parameters for Schwarz method.
Definition: FaspSolver.hpp:80
string solveDir
Current work dir.
Definition: FaspSolver.hpp:74
ITS_param itParam
Parameters for iterative method.
Definition: FaspSolver.hpp:77
void SetupParam(const string &dir, const string &file) override
Set FASP parameters.
Definition: FaspSolver.cpp:14
virtual void InitParam()=0
Initialize the params for linear solvers.