Actual source code: ex79.c

  1: #include <petsc.h>

  3: static char help[] = "Solves a linear system with a block of right-hand sides, apply a preconditioner to the same block.\n\n";

  5: PetscErrorCode MatApply(PC pc, Mat X, Mat Y)
  6: {

 10:   MatCopy(X,Y,SAME_NONZERO_PATTERN);
 11:   return(0);
 12: }

 14: int main(int argc,char **args)
 15: {
 16:   Mat                X,B;         /* computed solutions and RHS */
 17:   Mat                A;           /* linear system matrix */
 18:   KSP                ksp;         /* linear solver context */
 19:   PC                 pc;          /* preconditioner context */
 20:   PetscInt           m = 10;
 21: #if defined(PETSC_USE_LOG)
 22:   PetscLogEvent      event;
 23: #endif
 24:   PetscEventPerfInfo info;
 25:   PetscErrorCode     ierr;

 27:   PetscInitialize(&argc,&args,NULL,help);if (ierr) return ierr;
 28:   PetscLogDefaultBegin();
 29:   PetscOptionsGetInt(NULL,NULL,"-m",&m,NULL);
 30:   MatCreateAIJ(PETSC_COMM_WORLD,m,m,PETSC_DECIDE,PETSC_DECIDE,m,NULL,m,NULL,&A);
 31:   MatCreateDense(PETSC_COMM_WORLD,m,PETSC_DECIDE,PETSC_DECIDE,m,NULL,&B);
 32:   MatCreateDense(PETSC_COMM_WORLD,m,PETSC_DECIDE,PETSC_DECIDE,m,NULL,&X);
 33:   MatSetRandom(A,NULL);
 34:   MatSetOption(A,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_FALSE);
 35:   MatShift(A,10.0);
 36:   MatSetRandom(B,NULL);
 37:   KSPCreate(PETSC_COMM_WORLD,&ksp);
 38:   KSPSetOperators(ksp,A,A);
 39:   KSPSetFromOptions(ksp);
 40:   KSPGetPC(ksp,&pc);
 41:   PCShellSetMatApply(pc,MatApply);
 42:   KSPMatSolve(ksp,B,X);
 43:   PCMatApply(pc,B,X);
 44:   MatDestroy(&X);
 45:   MatDestroy(&B);
 46:   MatDestroy(&A);
 47:   KSPDestroy(&ksp);
 48:   PetscLogEventRegister("PCApply",PC_CLASSID,&event);
 49:   PetscLogEventGetPerfInfo(PETSC_DETERMINE,event,&info);
 50:   if (PetscDefined(USE_LOG) && m > 1 && info.count) SETERRQ1(PetscObjectComm((PetscObject)ksp),PETSC_ERR_PLIB,"PCApply() called %d times",info.count);
 51:   PetscLogEventRegister("PCMatApply",PC_CLASSID,&event);
 52:   PetscLogEventGetPerfInfo(PETSC_DETERMINE,event,&info);
 53:   if (PetscDefined(USE_LOG) && m > 1 && !info.count) SETERRQ(PetscObjectComm((PetscObject)ksp),PETSC_ERR_PLIB,"PCMatApply() never called");
 54:   PetscFinalize();
 55:   return ierr;
 56: }

 58: /*TEST
 59:    # KSPHPDDM does either pseudo-blocking or "true" blocking, all tests should succeed with other -ksp_hpddm_type
 60:    testset:
 61:       nsize: 1
 62:       args: -pc_type {{bjacobi lu ilu mat cholesky icc none shell asm gasm}shared output}
 63:       test:
 64:          suffix: 1
 65:          output_file: output/ex77_preonly.out
 66:          args: -ksp_type preonly
 67:       test:
 68:          suffix: 1_hpddm
 69:          output_file: output/ex77_preonly.out
 70:          requires: hpddm
 71:          args: -ksp_type hpddm -ksp_hpddm_type preonly

 73:    testset:
 74:       nsize: 1
 75:       args: -pc_type ksp
 76:       test:
 77:          suffix: 2
 78:          output_file: output/ex77_preonly.out
 79:          args: -ksp_ksp_type preonly -ksp_type preonly
 80:       test:
 81:          suffix: 2_hpddm
 82:          output_file: output/ex77_preonly.out
 83:          requires: hpddm
 84:          args: -ksp_ksp_type hpddm -ksp_type hpddm -ksp_hpddm_type preonly -ksp_ksp_hpddm_type preonly

 86:    testset:
 87:       nsize: 1
 88:       requires: h2opus !thrust
 89:       args: -pc_type h2opus
 90:       test:
 91:          suffix: 3
 92:          output_file: output/ex77_preonly.out
 93:          args: -ksp_type preonly
 94:       test:
 95:          suffix: 3_hpddm
 96:          output_file: output/ex77_preonly.out
 97:          requires: hpddm
 98:          args: -ksp_type hpddm -ksp_hpddm_type preonly

100:    testset:
101:       nsize: 1
102:       requires: spai
103:       args: -pc_type spai
104:       test:
105:          suffix: 4
106:          output_file: output/ex77_preonly.out
107:          args: -ksp_type preonly
108:       test:
109:          suffix: 4_hpddm
110:          output_file: output/ex77_preonly.out
111:          requires: hpddm
112:          args: -ksp_type hpddm -ksp_hpddm_type preonly
113:    # special code path in PCMatApply() for PCBJACOBI when a block is shared by multiple processes
114:    testset:
115:       nsize: 2
116:       args: -pc_type bjacobi -pc_bjacobi_blocks 1 -sub_pc_type none
117:       test:
118:          suffix: 5
119:          output_file: output/ex77_preonly.out
120:          args: -ksp_type preonly -sub_ksp_type preonly
121:       test:
122:          suffix: 5_hpddm
123:          output_file: output/ex77_preonly.out
124:          requires: hpddm
125:          args: -ksp_type hpddm -ksp_hpddm_type preonly -sub_ksp_type hpddm
126:    # special code path in PCMatApply() for PCGASM when a block is shared by multiple processes
127:    testset:
128:       nsize: 2
129:       args: -pc_type gasm -pc_gasm_total_subdomains 1 -sub_pc_type none
130:       test:
131:          suffix: 6
132:          output_file: output/ex77_preonly.out
133:          args: -ksp_type preonly -sub_ksp_type preonly
134:       test:
135:          suffix: 6_hpddm
136:          output_file: output/ex77_preonly.out
137:          requires: hpddm
138:          args: -ksp_type hpddm -ksp_hpddm_type preonly -sub_ksp_type hpddm

140:    testset:
141:       nsize: 1
142:       requires: suitesparse
143:       args: -pc_type qr
144:       test:
145:          suffix: 7
146:          output_file: output/ex77_preonly.out
147:          args: -ksp_type preonly
148:       test:
149:          suffix: 7_hpddm
150:          output_file: output/ex77_preonly.out
151:          requires: hpddm
152:          args: -ksp_type hpddm -ksp_hpddm_type preonly

154: TEST*/