Actual source code: ex77.c

petsc-master 2020-08-25
Report Typos and Errors
  1:  #include <petsc.h>

  3: static char help[] = "Solves a linear system with a block of right-hand sides using KSPHPDDM.\n\n";

  5: int main(int argc,char **args)
  6: {
  7:   Mat                X,B;         /* computed solutions and RHS */
  8:   Vec                cx,cb;       /* columns of X and B */
  9:   Mat                A,KA = NULL; /* linear system matrix */
 10:   KSP                ksp;         /* linear solver context */
 11:   PC                 pc;          /* preconditioner context */
 12:   Mat                F;           /* factored matrix from the preconditioner context */
 13:   PetscScalar        *x,*S = NULL,*T = NULL;
 14:   PetscReal          norm;
 15:   PetscInt           m,M,N = 5,i;
 16:   PetscMPIInt        rank,size;
 17:   const char         *deft = MATAIJ;
 18:   PetscViewer        viewer;
 19:   char               name[PETSC_MAX_PATH_LEN],type[256];
 20:   PetscBool          breakdown = PETSC_FALSE,flg;
 21:   KSPConvergedReason reason;
 22:   PetscErrorCode     ierr;

 24:   PetscInitialize(&argc,&args,NULL,help);if (ierr) return ierr;
 25:   MPI_Comm_rank(PETSC_COMM_WORLD,&rank);
 26:   MPI_Comm_size(PETSC_COMM_WORLD,&size);
 27:   PetscOptionsGetString(NULL,NULL,"-f",name,sizeof(name),&flg);
 28:   if (!flg) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"Must provide a binary file for the matrix");
 29:   PetscOptionsGetInt(NULL,NULL,"-N",&N,NULL);
 30:   PetscOptionsGetBool(NULL,NULL,"-breakdown",&breakdown,NULL);
 31:   MatCreate(PETSC_COMM_WORLD,&A);
 32:   KSPCreate(PETSC_COMM_WORLD,&ksp);
 33:   KSPSetOperators(ksp,A,A);
 34:   PetscViewerBinaryOpen(PETSC_COMM_WORLD,name,FILE_MODE_READ,&viewer);
 35:   MatLoad(A,viewer);
 36:   PetscViewerDestroy(&viewer);
 37:   PetscOptionsBegin(PETSC_COMM_WORLD,"","","");
 38:   PetscOptionsFList("-mat_type","Matrix type","MatSetType",MatList,deft,type,256,&flg);
 39:   PetscOptionsEnd();
 40:   if (flg) {
 41:     PetscStrcmp(type,MATKAIJ,&flg);
 42:     if (!flg) {
 43:       MatSetOption(A,MAT_SYMMETRIC,PETSC_TRUE);
 44:       MatConvert(A,type,MAT_INPLACE_MATRIX,&A);
 45:     } else {
 46:       if (size > 2) {
 47:         MatGetSize(A,&M,NULL);
 48:         MatCreate(PETSC_COMM_WORLD,&B);
 49:         if (rank > 1) {
 50:           MatSetSizes(B,0,0,M,M);
 51:         } else {
 52:           MatSetSizes(B,rank?M-M/2:M/2,rank?M-M/2:M/2,M,M);
 53:         }
 54:         PetscViewerBinaryOpen(PETSC_COMM_WORLD,name,FILE_MODE_READ,&viewer);
 55:         MatLoad(B,viewer);
 56:         PetscViewerDestroy(&viewer);
 57:         MatHeaderReplace(A,&B);
 58:       }
 59:       PetscCalloc2(N*N,&S,N*N,&T);
 60:       for (i=0; i<N; i++) { /* really easy problem used for testing */
 61:         S[i*(N+1)] = 1e+6;
 62:         T[i*(N+1)] = 1e-2;
 63:       }
 64:       MatCreateKAIJ(A,N,N,S,T,&KA);
 65:     }
 66:   }
 67:   if (!flg) {
 68:     if (size > 4) {
 69:       Mat B;
 70:       MatGetSize(A,&M,NULL);
 71:       MatCreate(PETSC_COMM_WORLD,&B);
 72:       if (rank > 3) {
 73:         MatSetSizes(B,0,0,M,M);
 74:       } else {
 75:         MatSetSizes(B,!rank?M-3*(M/4):M/4,!rank?M-3*(M/4):M/4,M,M);
 76:       }
 77:       PetscViewerBinaryOpen(PETSC_COMM_WORLD,name,FILE_MODE_READ,&viewer);
 78:       MatLoad(B,viewer);
 79:       PetscViewerDestroy(&viewer);
 80:       MatHeaderReplace(A,&B);
 81:     }
 82:   }
 83:   MatGetLocalSize(A,&m,NULL);
 84:   MatCreateDense(PETSC_COMM_WORLD,m,PETSC_DECIDE,PETSC_DECIDE,N,NULL,&B);
 85:   MatCreateDense(PETSC_COMM_WORLD,m,PETSC_DECIDE,PETSC_DECIDE,N,NULL,&X);
 86:   if (!breakdown) {
 87:     MatSetRandom(B,NULL);
 88:   }
 89:   KSPSetFromOptions(ksp);
 90:   if (!flg) {
 91:     if (!breakdown) {
 92:       KSPMatSolve(ksp,B,X);
 93:       KSPGetMatSolveBlockSize(ksp,&M);
 94:       if (M != PETSC_DECIDE) {
 95:         KSPSetMatSolveBlockSize(ksp,PETSC_DECIDE);
 96:         MatZeroEntries(X);
 97:         KSPMatSolve(ksp,B,X);
 98:       }
 99:       KSPGetPC(ksp,&pc);
100:       PetscObjectTypeCompare((PetscObject)pc,PCLU,&flg);
101:       if (flg) {
102:         PCFactorGetMatrix(pc,&F);
103:         MatMatSolve(F,B,B);
104:         MatAYPX(B,-1.0,X,SAME_NONZERO_PATTERN);
105:         MatNorm(B,NORM_INFINITY,&norm);
106:         if (norm > 100*PETSC_MACHINE_EPSILON) SETERRQ1(PetscObjectComm((PetscObject)ksp),PETSC_ERR_PLIB,"KSPMatSolve() and MatMatSolve() difference has nonzero norm %g",(double)norm);
107:       }
108:     } else {
109:       MatZeroEntries(B);
110:       KSPMatSolve(ksp,B,X);
111:       KSPGetConvergedReason(ksp,&reason);
112:       if (reason != KSP_CONVERGED_HAPPY_BREAKDOWN) SETERRQ1(PetscObjectComm((PetscObject)ksp),PETSC_ERR_PLIB,"KSPConvergedReason() %s != KSP_CONVERGED_HAPPY_BREAKDOWN",KSPConvergedReasons[reason]);
113:       MatDenseGetArrayWrite(B,&x);
114:       for (i=0; i<m*N; ++i) x[i] = 1.0;
115:       MatDenseRestoreArrayWrite(B,&x);
116:       KSPMatSolve(ksp,B,X);
117:       KSPGetConvergedReason(ksp,&reason);
118:       if (reason != KSP_DIVERGED_BREAKDOWN) SETERRQ1(PetscObjectComm((PetscObject)ksp),PETSC_ERR_PLIB,"KSPConvergedReason() %s != KSP_DIVERGED_BREAKDOWN",KSPConvergedReasons[reason]);
119:     }
120:   } else {
121:     KSPSetOperators(ksp,KA,KA);
122:     MatGetSize(KA,&M,NULL);
123:     VecCreateMPI(PETSC_COMM_WORLD,m*N,M,&cb);
124:     VecCreateMPI(PETSC_COMM_WORLD,m*N,M,&cx);
125:     VecSetRandom(cb,NULL);
126:     /* solving with MatKAIJ is equivalent to block solving with row-major RHS and solutions */
127:     /* only applies if MatKAIJGetScaledIdentity() returns true                              */
128:     KSPSolve(ksp,cb,cx);
129:     VecDestroy(&cx);
130:     VecDestroy(&cb);
131:   }
132:   MatDestroy(&X);
133:   MatDestroy(&B);
134:   PetscFree2(S,T);
135:   MatDestroy(&KA);
136:   MatDestroy(&A);
137:   KSPDestroy(&ksp);
138:   PetscFinalize();
139:   return ierr;
140: }

142: /*TEST

144:    testset:
145:       nsize: 2
146:       requires: datafilespath double !complex !define(PETSC_USE_64BIT_INDICES)
147:       args: -ksp_converged_reason -ksp_max_it 500 -f ${DATAFILESPATH}/matrices/hpddm/GCRODR/A_400.dat -mat_type {{aij sbaij}shared output}
148:       test:
149:          suffix: 1
150:          args:
151:       test:
152:          suffix: 2
153:          requires: hpddm
154:          args: -ksp_type hpddm -pc_type asm -ksp_hpddm_type {{gmres bgmres}separate output}
155:       test:
156:          suffix: 3
157:          requires: hpddm
158:          args: -ksp_type hpddm -ksp_hpddm_recycle 5 -ksp_hpddm_type {{gcrodr bgcrodr}separate output}
159:       test:
160:          nsize: 4
161:          suffix: 4
162:          requires: hpddm
163:          args: -ksp_rtol 1e-4 -ksp_type hpddm -ksp_hpddm_recycle 5 -ksp_hpddm_type bgcrodr -ksp_view_final_residual -N 12 -ksp_matsolve_block_size 5

165:    test:
166:       nsize: 1
167:       requires: hpddm datafilespath double !complex !define(PETSC_USE_64BIT_INDICES)
168:       suffix: preonly
169:       args: -N 6 -f ${DATAFILESPATH}/matrices/hpddm/GCRODR/A_400.dat -pc_type lu -ksp_type hpddm -ksp_hpddm_type preonly

171:    # to avoid breakdown failures, use -ksp_hpddm_deflation_tol, cf. KSPHPDDM documentation
172:    test:
173:       nsize: 1
174:       requires: hpddm datafilespath double !complex !define(PETSC_USE_64BIT_INDICES)
175:       suffix: breakdown
176:       output_file: output/ex77_preonly.out
177:       args: -N 3 -f ${DATAFILESPATH}/matrices/hpddm/GCRODR/A_400.dat -pc_type none -ksp_type hpddm -ksp_hpddm_type {{bcg bgmres bgcrodr bfbcg}shared output} -breakdown

179:    test:
180:       nsize: 2
181:       requires: hpddm datafilespath double !complex !define(PETSC_USE_64BIT_INDICES)
182:       args: -N 12 -ksp_converged_reason -ksp_max_it 500 -f ${DATAFILESPATH}/matrices/hpddm/GCRODR/A_400.dat -mat_type kaij -pc_type pbjacobi -ksp_type hpddm -ksp_hpddm_type {{gmres bgmres}separate output}

184:    test:
185:       nsize: 3
186:       requires: hpddm datafilespath double !complex !define(PETSC_USE_64BIT_INDICES)
187:       suffix: kaij_zero
188:       output_file: output/ex77_ksp_hpddm_type-bgmres.out
189:       args: -N 12 -ksp_converged_reason -ksp_max_it 500 -f ${DATAFILESPATH}/matrices/hpddm/GCRODR/A_400.dat -mat_type kaij -pc_type pbjacobi -ksp_type hpddm -ksp_hpddm_type bgmres

191:    test:
192:       nsize: 4
193:       requires: hpddm datafilespath double !complex !define(PETSC_USE_64BIT_INDICES) slepc
194:       suffix: 4_slepc
195:       output_file: output/ex77_4.out
196:       filter: sed "/^ksp_hpddm_recycle_ Linear eigensolve converged/d"
197:       args: -ksp_converged_reason -ksp_max_it 500 -f ${DATAFILESPATH}/matrices/hpddm/GCRODR/A_400.dat -ksp_rtol 1e-4 -ksp_type hpddm -ksp_hpddm_recycle 5 -ksp_hpddm_type bgcrodr -ksp_view_final_residual -N 12 -ksp_matsolve_block_size 5 -ksp_hpddm_recycle_redistribute 2 -ksp_hpddm_recycle_mat_type {{aij dense}shared output} -ksp_hpddm_recycle_eps_converged_reason -ksp_hpddm_recycle_st_pc_type redundant

199:    testset:
200:       nsize: 4
201:       requires: hpddm datafilespath double !complex !define(PETSC_USE_64BIT_INDICES) slepc elemental
202:       filter: sed "/^ksp_hpddm_recycle_ Linear eigensolve converged/d"
203:       args: -ksp_converged_reason -ksp_max_it 500 -f ${DATAFILESPATH}/matrices/hpddm/GCRODR/A_400.dat -ksp_rtol 1e-4 -ksp_type hpddm -ksp_hpddm_recycle 5 -ksp_hpddm_type bgcrodr -ksp_view_final_residual -N 12 -ksp_matsolve_block_size 5 -ksp_hpddm_recycle_redistribute 2 -ksp_hpddm_recycle_mat_type elemental -ksp_hpddm_recycle_eps_converged_reason
204:       test:
205:          suffix: 4_elemental
206:          output_file: output/ex77_4.out
207:       test:
208:          suffix: 4_elemental_matmat
209:          output_file: output/ex77_4.out
210:          args: -ksp_hpddm_recycle_eps_type subspace -ksp_hpddm_recycle_eps_target 0 -ksp_hpddm_recycle_eps_target_magnitude -ksp_hpddm_recycle_st_type sinvert -ksp_hpddm_recycle_bv_type mat -ksp_hpddm_recycle_bv_orthog_block svqb

212:    test:
213:       nsize: 5
214:       requires: hpddm datafilespath double !complex !define(PETSC_USE_64BIT_INDICES)
215:       suffix: 4_zero
216:       output_file: output/ex77_4.out
217:       args: -ksp_converged_reason -ksp_max_it 500 -f ${DATAFILESPATH}/matrices/hpddm/GCRODR/A_400.dat -ksp_rtol 1e-4 -ksp_type hpddm -ksp_hpddm_recycle 5 -ksp_hpddm_type bgcrodr -ksp_view_final_residual -N 12 -ksp_matsolve_block_size 5

219: TEST*/