Actual source code: ex49.c


  2: static char help[] = "Tests SeqSBAIJ factorizations for different block sizes\n\n";

  4: #include <petscksp.h>

  6: int main(int argc, char **args)
  7: {
  8:   Vec         x, b, u;
  9:   Mat         A, A2;
 10:   KSP         ksp;
 11:   PetscRandom rctx;
 12:   PetscReal   norm;
 13:   PetscInt    i, j, k, l, n = 27, its, bs = 2, Ii, J;
 14:   PetscBool   test_hermitian = PETSC_FALSE, convert = PETSC_FALSE;
 15:   PetscScalar v;

 17:   PetscFunctionBeginUser;
 18:   PetscCall(PetscInitialize(&argc, &args, (char *)0, help));
 19:   PetscCall(PetscOptionsGetInt(NULL, NULL, "-bs", &bs, NULL));
 20:   PetscCall(PetscOptionsGetInt(NULL, NULL, "-n", &n, NULL));
 21:   PetscCall(PetscOptionsGetBool(NULL, NULL, "-herm", &test_hermitian, NULL));
 22:   PetscCall(PetscOptionsGetBool(NULL, NULL, "-conv", &convert, NULL));

 24:   PetscCall(MatCreate(PETSC_COMM_SELF, &A));
 25:   PetscCall(MatSetSizes(A, n * bs, n * bs, PETSC_DETERMINE, PETSC_DETERMINE));
 26:   PetscCall(MatSetBlockSize(A, bs));
 27:   PetscCall(MatSetType(A, MATSEQSBAIJ));
 28:   PetscCall(MatSetFromOptions(A));
 29:   PetscCall(MatSeqSBAIJSetPreallocation(A, bs, n, NULL));
 30:   PetscCall(MatSeqBAIJSetPreallocation(A, bs, n, NULL));
 31:   PetscCall(MatSeqAIJSetPreallocation(A, n * bs, NULL));
 32:   PetscCall(MatMPIAIJSetPreallocation(A, n * bs, NULL, n * bs, NULL));

 34:   PetscCall(PetscRandomCreate(PETSC_COMM_SELF, &rctx));
 35:   for (i = 0; i < n; i++) {
 36:     for (j = i; j < n; j++) {
 37:       PetscCall(PetscRandomGetValue(rctx, &v));
 38:       if (PetscRealPart(v) < .1 || i == j) {
 39:         for (k = 0; k < bs; k++) {
 40:           for (l = 0; l < bs; l++) {
 41:             Ii = i * bs + k;
 42:             J  = j * bs + l;
 43:             PetscCall(PetscRandomGetValue(rctx, &v));
 44:             if (Ii == J) v = PetscRealPart(v + 3 * n * bs);
 45:             PetscCall(MatSetValue(A, Ii, J, v, INSERT_VALUES));
 46:             if (test_hermitian) {
 47:               PetscCall(MatSetValue(A, J, Ii, PetscConj(v), INSERT_VALUES));
 48:             } else {
 49:               PetscCall(MatSetValue(A, J, Ii, v, INSERT_VALUES));
 50:             }
 51:           }
 52:         }
 53:       }
 54:     }
 55:   }
 56:   PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
 57:   PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));

 59:   /* With complex numbers:
 60:      - PETSc cholesky does not support hermitian matrices
 61:      - CHOLMOD only supports hermitian matrices
 62:      - SUPERLU_DIST seems supporting both
 63:   */
 64:   if (test_hermitian) PetscCall(MatSetOption(A, MAT_HERMITIAN, PETSC_TRUE));

 66:   {
 67:     Mat M;
 68:     PetscCall(MatComputeOperator(A, MATAIJ, &M));
 69:     PetscCall(MatViewFromOptions(M, NULL, "-expl_view"));
 70:     PetscCall(MatDestroy(&M));
 71:   }

 73:   A2 = NULL;
 74:   if (convert) PetscCall(MatConvert(A, MATAIJ, MAT_INITIAL_MATRIX, &A2));

 76:   PetscCall(VecCreate(PETSC_COMM_SELF, &u));
 77:   PetscCall(VecSetSizes(u, PETSC_DECIDE, n * bs));
 78:   PetscCall(VecSetFromOptions(u));
 79:   PetscCall(VecDuplicate(u, &b));
 80:   PetscCall(VecDuplicate(b, &x));

 82:   PetscCall(VecSet(u, 1.0));
 83:   PetscCall(MatMult(A, u, b));

 85:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 86:                 Create the linear solver and set various options
 87:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 89:   /*
 90:      Create linear solver context
 91:   */
 92:   PetscCall(KSPCreate(PETSC_COMM_SELF, &ksp));

 94:   /*
 95:      Set operators.
 96:   */
 97:   PetscCall(KSPSetOperators(ksp, A2 ? A2 : A, A));

 99:   PetscCall(KSPSetFromOptions(ksp));

101:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
102:                       Solve the linear system
103:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

105:   PetscCall(KSPSolve(ksp, b, x));

107:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
108:                       Check solution and clean up
109:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

111:   /*
112:      Check the error
113:   */
114:   PetscCall(VecAXPY(x, -1.0, u));
115:   PetscCall(VecNorm(x, NORM_2, &norm));
116:   PetscCall(KSPGetIterationNumber(ksp, &its));

118:   /*
119:      Print convergence information.  PetscPrintf() produces a single
120:      print statement from all processes that share a communicator.
121:      An alternative is PetscFPrintf(), which prints to a file.
122:   */
123:   if (norm > 100 * PETSC_SMALL) PetscCall(PetscPrintf(PETSC_COMM_SELF, "Norm of residual %g iterations %" PetscInt_FMT " bs %" PetscInt_FMT "\n", (double)norm, its, bs));

125:   /*
126:      Free work space.  All PETSc objects should be destroyed when they
127:      are no longer needed.
128:   */
129:   PetscCall(KSPDestroy(&ksp));
130:   PetscCall(VecDestroy(&u));
131:   PetscCall(VecDestroy(&x));
132:   PetscCall(VecDestroy(&b));
133:   PetscCall(MatDestroy(&A));
134:   PetscCall(MatDestroy(&A2));
135:   PetscCall(PetscRandomDestroy(&rctx));

137:   /*
138:      Always call PetscFinalize() before exiting a program.  This routine
139:        - finalizes the PETSc libraries as well as MPI
140:        - provides summary and diagnostic information if certain runtime
141:          options are chosen (e.g., -log_view).
142:   */
143:   PetscCall(PetscFinalize());
144:   return 0;
145: }

147: /*TEST

149:    test:
150:       args: -mat_type {{aij baij sbaij}} -bs {{1 2 3 4 5 6 7 8 9 10 11 12}} -pc_type cholesky -herm 0 -conv {{0 1}}

152:    test:
153:       nsize: {{1 4}}
154:       suffix: cholmod
155:       requires: suitesparse
156:       args: -mat_type {{aij sbaij}} -bs 1 -pc_type cholesky -pc_factor_mat_solver_type cholmod -herm -conv {{0 1}}

158:    test:
159:       nsize: {{1 4}}
160:       suffix: superlu_dist
161:       requires: superlu_dist
162:       output_file: output/ex49_cholmod.out
163:       args: -mat_type mpiaij -bs 3 -pc_type cholesky -pc_factor_mat_solver_type superlu_dist -herm {{0 1}} -conv {{0 1}}

165:    test:
166:       suffix: mkl_pardiso
167:       requires: mkl_pardiso
168:       output_file: output/ex49_1.out
169:       args: -bs {{1 3}} -pc_type cholesky -pc_factor_mat_solver_type mkl_pardiso

171:    test:
172:       suffix: cg
173:       requires: complex
174:       output_file: output/ex49_cg.out
175:       args: -herm -ksp_cg_type hermitian -mat_type aij -ksp_type cg -pc_type jacobi -ksp_rtol 4e-07

177:    test:
178:       suffix: pipecg2
179:       requires: complex
180:       output_file: output/ex49_pipecg2.out
181:       args: -herm -mat_type aij -ksp_type pipecg2 -pc_type jacobi -ksp_rtol 4e-07 -ksp_norm_type {{preconditioned unpreconditioned natural}}

183: TEST*/