Actual source code: ex75.c
petsc-master 2019-12-03
1: #include <petsc.h>
3: static char help[] = "Solves a series of linear systems using KSPHPDDM.\n\n";
5: int main(int argc,char **args)
6: {
7: Vec x,b; /* computed solution and RHS */
8: Mat A; /* linear system matrix */
9: KSP ksp; /* linear solver context */
10: PetscInt i,j,nmat = 10;
11: PetscViewer viewer;
12: char name[256];
13: char dir[PETSC_MAX_PATH_LEN];
16: PetscInitialize(&argc,&args,NULL,help);if (ierr) return ierr;
17: PetscStrcpy(dir,".");
18: PetscOptionsGetString(NULL,NULL,"-load_dir",dir,sizeof(dir),NULL);
19: PetscOptionsGetInt(NULL,NULL,"-nmat",&nmat,NULL);
20: MatCreate(PETSC_COMM_WORLD,&A);
21: KSPCreate(PETSC_COMM_WORLD,&ksp);
22: KSPSetOperators(ksp,A,A);
23: for (i=0; i<nmat; i++) {
24: j = i+400;
25: PetscSNPrintf(name,sizeof(name),"%s/A_%d.dat",dir,j);
26: PetscViewerBinaryOpen(PETSC_COMM_WORLD,name,FILE_MODE_READ,&viewer);
27: MatLoad(A,viewer);
28: PetscViewerDestroy(&viewer);
29: MatCreateVecs(A,&x,&b);
30: PetscSNPrintf(name,sizeof(name),"%s/rhs_%d.dat",dir,j);
31: PetscViewerBinaryOpen(PETSC_COMM_WORLD,name,FILE_MODE_READ,&viewer);
32: VecLoad(b,viewer);
33: PetscViewerDestroy(&viewer);
34: KSPSetOperators(ksp,A,A);
35: KSPSetFromOptions(ksp);
36: KSPSolve(ksp,b,x);
37: VecDestroy(&x);
38: VecDestroy(&b);
39: }
40: MatDestroy(&A);
41: KSPDestroy(&ksp);
42: PetscFinalize();
43: return ierr;
44: }
46: /*TEST
48: test:
49: suffix: 1
50: nsize: 1
51: requires: hpddm datafilespath double !complex !define(PETSC_USE_64BIT_INDICES)
52: args: -nmat 1 -pc_type none -ksp_converged_reason -ksp_type {{gmres hpddm}shared ouput} -ksp_max_it 1000 -ksp_gmres_restart 1000 -ksp_rtol 1e-10 -load_dir ${DATAFILESPATH}/matrices/hpddm/GCRODR
54: test:
55: requires: hpddm datafilespath double !complex !define(PETSC_USE_64BIT_INDICES)
56: suffix: 1_icc
57: nsize: 1
58: args: -nmat 1 -pc_type icc -ksp_converged_reason -ksp_type {{gmres hpddm}shared ouput} -ksp_max_it 1000 -ksp_gmres_restart 1000 -ksp_rtol 1e-10 -load_dir ${DATAFILESPATH}/matrices/hpddm/GCRODR
60: testset:
61: requires: hpddm datafilespath double !complex !define(PETSC_USE_64BIT_INDICES)
62: args: -nmat 3 -pc_type none -ksp_converged_reason -ksp_type hpddm -ksp_max_it 1000 -ksp_gmres_restart 40 -ksp_rtol 1e-10 -ksp_hpddm_krylov_method gcrodr -ksp_hpddm_recycle 20 -load_dir ${DATAFILESPATH}/matrices/hpddm/GCRODR
63: test:
64: nsize: 1
65: suffix: 2_seq
66: output_file: output/ex75_2.out
67: test:
68: nsize: 2
69: suffix: 2_par
70: output_file: output/ex75_2.out
72: test:
73: requires: hpddm datafilespath double !complex !define(PETSC_USE_64BIT_INDICES)
74: suffix: 2_icc
75: nsize: 1
76: args: -nmat 3 -pc_type icc -ksp_converged_reason -ksp_type hpddm -ksp_max_it 1000 -ksp_gmres_restart 40 -ksp_rtol 1e-10 -ksp_hpddm_krylov_method gcrodr -ksp_hpddm_recycle 20 -load_dir ${DATAFILESPATH}/matrices/hpddm/GCRODR
78: TEST*/