xref: /petsc/src/mat/tests/ex134.c (revision 910b42cbfb1b6524ef4594118dd5e013aee7da5d)
1c4762a1bSJed Brown static const char help[] = "Test parallel assembly of SBAIJ matrices\n\n";
2c4762a1bSJed Brown 
3c4762a1bSJed Brown #include <petscmat.h>
4c4762a1bSJed Brown 
5d71ae5a4SJacob Faibussowitsch PetscErrorCode Assemble(MPI_Comm comm, PetscInt bs, MatType mtype)
6d71ae5a4SJacob Faibussowitsch {
7c4762a1bSJed Brown   const PetscInt    rc[]   = {0, 1, 2, 3};
89371c9d4SSatish Balay   const PetscScalar vals[] = {100, 2,  3,  4,  5,  600, 7,  8,  9,  100, 11, 1200, 13, 14, 15, 1600, 17, 18, 19, 20, 21, 22, 23, 24, 25, 26, 27, 2800, 29, 30, 31, 32,
99371c9d4SSatish Balay                               33,  34, 35, 36, 37, 38,  39, 40, 41, 42,  43, 44,   45, 46, 47, 48,   49, 50, 51, 52, 53, 54, 55, 56, 57, 58, 49, 60,   61, 62, 63, 64};
10c4762a1bSJed Brown   Mat               A;
11c4762a1bSJed Brown #if defined(PETSC_HAVE_MUMPS) || defined(PETSC_HAVE_MKL_CPARDISO)
12c4762a1bSJed Brown   Mat           F;
1383863137SStefano Zampini   MatSolverType stype = MATSOLVERPETSC;
14c4762a1bSJed Brown   PetscRandom   rdm;
15c4762a1bSJed Brown   Vec           b, x, y;
16c4762a1bSJed Brown   PetscInt      i, j;
17*910b42cbSStefano Zampini   PetscReal     norm2, tol = 100 * PETSC_SQRT_MACHINE_EPSILON;
18c4762a1bSJed Brown   PetscBool     issbaij;
19c4762a1bSJed Brown #endif
20c4762a1bSJed Brown   PetscViewer viewer;
21c4762a1bSJed Brown 
22c4762a1bSJed Brown   PetscFunctionBegin;
239566063dSJacob Faibussowitsch   PetscCall(MatCreate(comm, &A));
249566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(A, PETSC_DECIDE, PETSC_DECIDE, 4 * bs, 4 * bs));
259566063dSJacob Faibussowitsch   PetscCall(MatSetType(A, mtype));
269566063dSJacob Faibussowitsch   PetscCall(MatMPIBAIJSetPreallocation(A, bs, 2, NULL, 2, NULL));
279566063dSJacob Faibussowitsch   PetscCall(MatMPISBAIJSetPreallocation(A, bs, 2, NULL, 2, NULL));
289566063dSJacob Faibussowitsch   PetscCall(MatSetOption(A, MAT_NEW_NONZERO_ALLOCATION_ERR, PETSC_TRUE));
29c4762a1bSJed Brown   /* All processes contribute a global matrix */
309566063dSJacob Faibussowitsch   PetscCall(MatSetValuesBlocked(A, 4, rc, 4, rc, vals, ADD_VALUES));
319566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
329566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
339566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(comm, "Matrix %s(%" PetscInt_FMT ")\n", mtype, bs));
349566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIGetStdout(comm, &viewer));
359566063dSJacob Faibussowitsch   PetscCall(PetscViewerPushFormat(viewer, PETSC_VIEWER_ASCII_INFO_DETAIL));
369566063dSJacob Faibussowitsch   PetscCall(MatView(A, viewer));
379566063dSJacob Faibussowitsch   PetscCall(PetscViewerPopFormat(viewer));
389566063dSJacob Faibussowitsch   PetscCall(MatView(A, viewer));
39c4762a1bSJed Brown #if defined(PETSC_HAVE_MUMPS) || defined(PETSC_HAVE_MKL_CPARDISO)
409566063dSJacob Faibussowitsch   PetscCall(PetscStrcmp(mtype, MATMPISBAIJ, &issbaij));
4148a46eb9SPierre Jolivet   if (!issbaij) PetscCall(MatShift(A, 10));
429566063dSJacob Faibussowitsch   PetscCall(PetscRandomCreate(PETSC_COMM_WORLD, &rdm));
439566063dSJacob Faibussowitsch   PetscCall(PetscRandomSetFromOptions(rdm));
449566063dSJacob Faibussowitsch   PetscCall(MatCreateVecs(A, &x, &y));
459566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(x, &b));
46c4762a1bSJed Brown   for (j = 0; j < 2; j++) {
47c4762a1bSJed Brown   #if defined(PETSC_HAVE_MUMPS)
48c4762a1bSJed Brown     if (j == 0) stype = MATSOLVERMUMPS;
49c4762a1bSJed Brown   #else
50c4762a1bSJed Brown     if (j == 0) continue;
51c4762a1bSJed Brown   #endif
52c4762a1bSJed Brown   #if defined(PETSC_HAVE_MKL_CPARDISO)
539131cd0dSPierre Jolivet     if (j == 1) {
549131cd0dSPierre Jolivet       if (issbaij && PetscDefined(USE_COMPLEX)) continue;
559131cd0dSPierre Jolivet       stype = MATSOLVERMKL_CPARDISO;
569131cd0dSPierre Jolivet     }
57c4762a1bSJed Brown   #else
58c4762a1bSJed Brown     if (j == 1) continue;
59c4762a1bSJed Brown   #endif
60c4762a1bSJed Brown     if (issbaij) {
619566063dSJacob Faibussowitsch       PetscCall(MatGetFactor(A, stype, MAT_FACTOR_CHOLESKY, &F));
629566063dSJacob Faibussowitsch       PetscCall(MatCholeskyFactorSymbolic(F, A, NULL, NULL));
639566063dSJacob Faibussowitsch       PetscCall(MatCholeskyFactorNumeric(F, A, NULL));
64c4762a1bSJed Brown     } else {
659566063dSJacob Faibussowitsch       PetscCall(MatGetFactor(A, stype, MAT_FACTOR_LU, &F));
669566063dSJacob Faibussowitsch       PetscCall(MatLUFactorSymbolic(F, A, NULL, NULL, NULL));
679566063dSJacob Faibussowitsch       PetscCall(MatLUFactorNumeric(F, A, NULL));
68c4762a1bSJed Brown     }
69c4762a1bSJed Brown     for (i = 0; i < 10; i++) {
709566063dSJacob Faibussowitsch       PetscCall(VecSetRandom(b, rdm));
719566063dSJacob Faibussowitsch       PetscCall(MatSolve(F, b, y));
72c4762a1bSJed Brown       /* Check the error */
739566063dSJacob Faibussowitsch       PetscCall(MatMult(A, y, x));
749566063dSJacob Faibussowitsch       PetscCall(VecAXPY(x, -1.0, b));
759566063dSJacob Faibussowitsch       PetscCall(VecNorm(x, NORM_2, &norm2));
7648a46eb9SPierre Jolivet       if (norm2 > tol) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Error:MatSolve(), norm2: %g\n", (double)norm2));
77c4762a1bSJed Brown     }
789566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&F));
79c4762a1bSJed Brown   }
809566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&x));
819566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&y));
829566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&b));
839566063dSJacob Faibussowitsch   PetscCall(PetscRandomDestroy(&rdm));
84c4762a1bSJed Brown #endif
859566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&A));
863ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
87c4762a1bSJed Brown }
88c4762a1bSJed Brown 
89d71ae5a4SJacob Faibussowitsch int main(int argc, char *argv[])
90d71ae5a4SJacob Faibussowitsch {
91c4762a1bSJed Brown   MPI_Comm    comm;
92c4762a1bSJed Brown   PetscMPIInt size;
93c4762a1bSJed Brown 
94327415f7SBarry Smith   PetscFunctionBeginUser;
959566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
96c4762a1bSJed Brown   comm = PETSC_COMM_WORLD;
979566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(comm, &size));
98be096a46SBarry Smith   PetscCheck(size == 2, comm, PETSC_ERR_USER, "This example must be run with exactly two processes");
999566063dSJacob Faibussowitsch   PetscCall(Assemble(comm, 2, MATMPIBAIJ));
1009566063dSJacob Faibussowitsch   PetscCall(Assemble(comm, 2, MATMPISBAIJ));
1019566063dSJacob Faibussowitsch   PetscCall(Assemble(comm, 1, MATMPIBAIJ));
1029566063dSJacob Faibussowitsch   PetscCall(Assemble(comm, 1, MATMPISBAIJ));
1039566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
104b122ec5aSJacob Faibussowitsch   return 0;
105c4762a1bSJed Brown }
106c4762a1bSJed Brown 
107c4762a1bSJed Brown /*TEST
108c4762a1bSJed Brown 
109c4762a1bSJed Brown    test:
110c4762a1bSJed Brown       nsize: 2
11197929ea7SJunchao Zhang       args: -mat_ignore_lower_triangular
112c4762a1bSJed Brown       filter: sed -e "s~mem [0-9]*~mem~g"
113c4762a1bSJed Brown 
114c4762a1bSJed Brown TEST*/
115