Actual source code: ex18.c
petsc-3.12.1 2019-10-22
2: static char help[] = "Reads a PETSc matrix and vector from a file and solves a linear system.\n\
3: Input arguments are:\n\
4: -f <input_file> : file to load. For example see $PETSC_DIR/share/petsc/datafiles/matrices\n\n";
6: #include <petscmat.h>
7: #include <petscksp.h>
9: int main(int argc,char **args)
10: {
12: PetscInt its,m,n,mvec;
13: PetscReal norm;
14: Vec x,b,u;
15: Mat A;
16: KSP ksp;
17: char file[PETSC_MAX_PATH_LEN];
18: PetscViewer fd;
19: PetscLogStage stage1;
21: PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr;
23: /* Read matrix and RHS */
24: PetscOptionsGetString(NULL,NULL,"-f",file,PETSC_MAX_PATH_LEN,NULL);
25: PetscViewerBinaryOpen(PETSC_COMM_WORLD,file,FILE_MODE_READ,&fd);
26: MatCreate(PETSC_COMM_WORLD,&A);
27: MatSetType(A,MATSEQAIJ);
28: MatLoad(A,fd);
29: VecCreate(PETSC_COMM_WORLD,&b);
30: VecLoad(b,fd);
31: PetscViewerDestroy(&fd);
33: /*
34: If the load matrix is larger then the vector, due to being padded
35: to match the blocksize then create a new padded vector
36: */
37: MatGetSize(A,&m,&n);
38: VecGetSize(b,&mvec);
39: if (m > mvec) {
40: Vec tmp;
41: PetscScalar *bold,*bnew;
42: /* create a new vector b by padding the old one */
43: VecCreate(PETSC_COMM_WORLD,&tmp);
44: VecSetSizes(tmp,PETSC_DECIDE,m);
45: VecSetFromOptions(tmp);
46: VecGetArray(tmp,&bnew);
47: VecGetArray(b,&bold);
48: PetscArraycpy(bnew,bold,mvec);
49: VecDestroy(&b);
50: b = tmp;
51: }
53: /* Set up solution */
54: VecDuplicate(b,&x);
55: VecDuplicate(b,&u);
56: VecSet(x,0.0);
58: /* Solve system */
59: PetscLogStageRegister("Stage 1",&stage1);
60: PetscLogStagePush(stage1);
61: KSPCreate(PETSC_COMM_WORLD,&ksp);
62: KSPSetOperators(ksp,A,A);
63: KSPSetFromOptions(ksp);
64: KSPSolve(ksp,b,x);
65: PetscLogStagePop();
67: /* Show result */
68: MatMult(A,x,u);
69: VecAXPY(u,-1.0,b);
70: VecNorm(u,NORM_2,&norm);
71: KSPGetIterationNumber(ksp,&its);
72: PetscPrintf(PETSC_COMM_WORLD,"Number of iterations = %3D\n",its);
73: PetscPrintf(PETSC_COMM_WORLD,"Residual norm %g\n",(double)norm);
75: /* Cleanup */
76: KSPDestroy(&ksp);
77: VecDestroy(&x);
78: VecDestroy(&b);
79: VecDestroy(&u);
80: MatDestroy(&A);
82: PetscFinalize();
83: return ierr;
84: }
86: /*TEST
88: test:
89: args: -ksp_gmres_cgs_refinement_type refine_always -f ${DATAFILESPATH}/matrices/arco1 -ksp_monitor_short
90: requires: datafilespath double !complex !define(PETSC_USE_64BIT_INDICES)
92: TEST*/