Actual source code: ex6.c

  1: static char help[] = "Reads a PETSc matrix and vector from a file and solves a linear system.\n\
  2: Input arguments are:\n\
  3:   -f <input_file> : file to load. For example see $PETSC_DIR/share/petsc/datafiles/matrices\n\n";

  5: #include <petscksp.h>
  6: #include <petsclog.h>

  8: static PetscErrorCode KSPTestResidualMonitor(KSP ksp, PetscInt i, PetscReal r, void* ctx)
  9: {
 10:   Vec            *t,*v;
 11:   PetscReal      err;

 14:   KSPCreateVecs(ksp,2,&t,2,&v);
 15:   KSPBuildResidualDefault(ksp,t[0],v[0],&v[0]);
 16:   KSPBuildResidual(ksp,t[1],v[1],&v[1]);
 17:   VecAXPY(v[1],-1.0,v[0]);
 18:   VecNorm(v[1],NORM_INFINITY,&err);
 20:   VecDestroyVecs(2,&t);
 21:   VecDestroyVecs(2,&v);
 22:   return 0;
 23: }

 25: int main(int argc,char **args)
 26: {
 27:   PetscInt       its;
 28: #if defined(PETSC_USE_LOG)
 29:   PetscLogStage  stage1,stage2;
 30: #endif
 31:   PetscReal      norm;
 32:   Vec            x,b,u;
 33:   Mat            A;
 34:   char           file[PETSC_MAX_PATH_LEN];
 35:   PetscViewer    fd;
 36:   PetscBool      table = PETSC_FALSE,flg,test_residual = PETSC_FALSE,b_in_f = PETSC_TRUE;
 37:   KSP            ksp;

 39:   PetscInitialize(&argc,&args,(char*)0,help);
 40:   PetscOptionsGetBool(NULL,NULL,"-table",&table,NULL);
 41:   PetscOptionsGetBool(NULL,NULL,"-test_residual",&test_residual,NULL);
 42:   PetscOptionsGetBool(NULL,NULL,"-b_in_f",&b_in_f,NULL);

 44:   /* Read matrix and RHS */
 45:   PetscOptionsGetString(NULL,NULL,"-f",file,sizeof(file),&flg);
 47:   PetscViewerBinaryOpen(PETSC_COMM_WORLD,file,FILE_MODE_READ,&fd);
 48:   MatCreate(PETSC_COMM_WORLD,&A);
 49:   MatLoad(A,fd);
 50:   if (b_in_f) {
 51:     VecCreate(PETSC_COMM_WORLD,&b);
 52:     VecLoad(b,fd);
 53:   } else {
 54:     MatCreateVecs(A,NULL,&b);
 55:     VecSetRandom(b,NULL);
 56:   }
 57:   PetscViewerDestroy(&fd);

 59:   /*
 60:    If the load matrix is larger then the vector, due to being padded
 61:    to match the blocksize then create a new padded vector
 62:   */
 63:   {
 64:     PetscInt    m,n,j,mvec,start,end,indx;
 65:     Vec         tmp;
 66:     PetscScalar *bold;

 68:     MatGetLocalSize(A,&m,&n);
 69:     VecCreate(PETSC_COMM_WORLD,&tmp);
 70:     VecSetSizes(tmp,m,PETSC_DECIDE);
 71:     VecSetFromOptions(tmp);
 72:     VecGetOwnershipRange(b,&start,&end);
 73:     VecGetLocalSize(b,&mvec);
 74:     VecGetArray(b,&bold);
 75:     for (j=0; j<mvec; j++) {
 76:       indx = start+j;
 77:       VecSetValues(tmp,1,&indx,bold+j,INSERT_VALUES);
 78:     }
 79:     VecRestoreArray(b,&bold);
 80:     VecDestroy(&b);
 81:     VecAssemblyBegin(tmp);
 82:     VecAssemblyEnd(tmp);
 83:     b    = tmp;
 84:   }
 85:   VecDuplicate(b,&x);
 86:   VecDuplicate(b,&u);

 88:   VecSet(x,0.0);
 89:   PetscBarrier((PetscObject)A);

 91:   PetscLogStageRegister("mystage 1",&stage1);
 92:   PetscLogStagePush(stage1);
 93:   KSPCreate(PETSC_COMM_WORLD,&ksp);
 94:   KSPSetOperators(ksp,A,A);
 95:   KSPSetFromOptions(ksp);
 96:   if (test_residual) {
 97:     KSPMonitorSet(ksp,KSPTestResidualMonitor,NULL,NULL);
 98:   }
 99:   KSPSetUp(ksp);
100:   KSPSetUpOnBlocks(ksp);
101:   PetscLogStagePop();
102:   PetscBarrier((PetscObject)A);

104:   PetscLogStageRegister("mystage 2",&stage2);
105:   PetscLogStagePush(stage2);
106:   KSPSolve(ksp,b,x);
107:   PetscLogStagePop();

109:   /* Show result */
110:   MatMult(A,x,u);
111:   VecAXPY(u,-1.0,b);
112:   VecNorm(u,NORM_2,&norm);
113:   KSPGetIterationNumber(ksp,&its);
114:   /*  matrix PC   KSP   Options       its    residual  */
115:   if (table) {
116:     char        *matrixname,kspinfo[120];
117:     PetscViewer viewer;
118:     PetscViewerStringOpen(PETSC_COMM_WORLD,kspinfo,sizeof(kspinfo),&viewer);
119:     KSPView(ksp,viewer);
120:     PetscStrrchr(file,'/',&matrixname);
121:     PetscPrintf(PETSC_COMM_WORLD,"%-8.8s %3D %2.0e %s \n",matrixname,its,norm,kspinfo);
122:     PetscViewerDestroy(&viewer);
123:   } else {
124:     PetscPrintf(PETSC_COMM_WORLD,"Number of iterations = %3D\n",its);
125:     PetscPrintf(PETSC_COMM_WORLD,"Residual norm = %g\n",(double)norm);
126:   }

128:   /* Cleanup */
129:   KSPDestroy(&ksp);
130:   VecDestroy(&x);
131:   VecDestroy(&b);
132:   VecDestroy(&u);
133:   MatDestroy(&A);
134:   PetscFinalize();
135:   return 0;
136: }

138: /*TEST

140:     test:
141:       args: -ksp_type preonly  -pc_type lu -options_left no  -f ${DATAFILESPATH}/matrices/arco1
142:       requires: datafilespath double !complex !defined(PETSC_USE_64BIT_INDICES)

144:     test:
145:       suffix: 2
146:       args: -sub_pc_type ilu -options_left no  -f ${DATAFILESPATH}/matrices/arco1 -ksp_gmres_restart 100 -ksp_gmres_cgs_refinement_type refine_always -sub_ksp_type preonly -pc_type bjacobi -pc_bjacobi_blocks 8 -sub_pc_factor_in_place -ksp_monitor_short
147:       requires: datafilespath double  !complex !defined(PETSC_USE_64BIT_INDICES)

149:     test:
150:       suffix: 7
151:       args: -ksp_gmres_cgs_refinement_type refine_always -pc_type asm -pc_asm_blocks 6 -f ${DATAFILESPATH}/matrices/small -matload_block_size 6  -ksp_monitor_short
152:       requires: datafilespath double  !complex !defined(PETSC_USE_64BIT_INDICES)

154:     test:
155:       requires: double !complex !defined(PETSC_USE_64BIT_INDICES)
156:       suffix: 3
157:       filter: sed -e "s/CONVERGED_RTOL/CONVERGED_ATOL/g"
158:       args: -f ${wPETSC_DIR}/share/petsc/datafiles/matrices/spd-real-int32-float64 -pc_type none -ksp_type {{cg groppcg pipecg pipecgrr pipelcg pipeprcg cgne nash stcg gltr fcg pipefcg gmres pipefgmres fgmres lgmres dgmres pgmres tcqmr bcgs ibcgs qmrcgs fbcgs fbcgsr bcgsl pipebcgs cgs tfqmr cr pipecr lsqr qcg bicg minres symmlq lcd gcr pipegcr cgls}} -ksp_max_it 20 -ksp_error_if_not_converged -ksp_converged_reason -test_residual

160:     test:
161:       requires: double !complex !defined(PETSC_USE_64BIT_INDICES)
162:       suffix: 3_maxits
163:       output_file: output/ex6_maxits.out
164:       args: -f ${wPETSC_DIR}/share/petsc/datafiles/matrices/spd-real-int32-float64 -pc_type none -ksp_type {{chebyshev cg groppcg pipecg pipecgrr pipelcg pipeprcg cgne nash stcg gltr fcg pipefcg gmres pipefgmres fgmres lgmres dgmres pgmres tcqmr bcgs ibcgs qmrcgs fbcgs fbcgsr bcgsl pipebcgs cgs tfqmr cr pipecr qcg bicg minres symmlq lcd gcr pipegcr cgls richardson}} -ksp_max_it 4 -ksp_error_if_not_converged -ksp_converged_maxits -ksp_converged_reason -test_residual -ksp_norm_type none

166:     testset:
167:       requires: double !complex !defined(PETSC_USE_64BIT_INDICES)
168:       output_file: output/ex6_skip.out
169:       args: -f ${wPETSC_DIR}/share/petsc/datafiles/matrices/spd-real-int32-float64 -pc_type none -ksp_max_it 8 -ksp_error_if_not_converged -ksp_convergence_test skip -ksp_converged_reason -test_residual
170:       #SYMMLQ converges in 4 iterations and then generate nans
171:       test:
172:         suffix: 3_skip
173:         args: -ksp_type {{chebyshev cg groppcg pipecg pipecgrr pipelcg pipeprcg cgne nash stcg gltr fcg pipefcg gmres fgmres lgmres dgmres pgmres tcqmr bcgs ibcgs qmrcgs fbcgs fbcgsr bcgsl pipebcgs cgs tfqmr cr pipecr qcg bicg minres lcd gcr cgls richardson}}
174:       test:
175:         requires: !pgf90_compiler
176:         suffix: 3_skip_pipefgmres
177:         args: -ksp_type pipefgmres
178:       #PIPEGCR generates nans on linux-knl
179:       test:
180:         requires: !defined(PETSC_USE_AVX512_KERNELS)
181:         suffix: 3_skip_pipegcr
182:         args: -ksp_type pipegcr
183:       test:
184:         requires: hpddm
185:         suffix: 3_skip_hpddm
186:         args: -ksp_type hpddm -ksp_hpddm_type {{cg gmres bgmres bcg bfbcg gcrodr bgcrodr}}

188:     test:
189:       requires: double !complex !defined(PETSC_USE_64BIT_INDICES) hpddm
190:       suffix: 3_hpddm
191:       output_file: output/ex6_3.out
192:       filter: sed -e "s/CONVERGED_RTOL/CONVERGED_ATOL/g"
193:       args: -f ${wPETSC_DIR}/share/petsc/datafiles/matrices/spd-real-int32-float64 -pc_type none -ksp_type hpddm -ksp_hpddm_type {{cg gmres bgmres bcg bfbcg gcrodr bgcrodr}} -ksp_max_it 20 -ksp_error_if_not_converged -ksp_converged_reason -test_residual

195:     # test CG shortcut for residual access
196:     test:
197:       suffix: 4
198:       args: -ksp_converged_reason -ksp_max_it 20 -ksp_converged_maxits -ksp_type {{cg pipecg groppcg}} -ksp_norm_type {{preconditioned unpreconditioned natural}separate output} -pc_type {{bjacobi none}separate output} -f ${DATAFILESPATH}/matrices/poisson_2d13p -b_in_f 0 -test_residual
199:       requires: datafilespath double !complex !defined(PETSC_USE_64BIT_INDICES)

201: TEST*/