xref: /petsc/src/mat/tests/ex195.c (revision f3fa974cd8ac3307c23c9d7e81b13b5f402f9e6e)
1c4762a1bSJed Brown /*
2c4762a1bSJed Brown  * ex195.c
3c4762a1bSJed Brown  *
4c4762a1bSJed Brown  *  Created on: Aug 24, 2015
5c4762a1bSJed Brown  *      Author: Fande Kong <fdkong.jd@gmail.com>
6c4762a1bSJed Brown  */
7c4762a1bSJed Brown 
8c4762a1bSJed Brown static char help[] = " Demonstrate the use of MatConvert_Nest_AIJ\n";
9c4762a1bSJed Brown 
10c4762a1bSJed Brown #include <petscmat.h>
11c4762a1bSJed Brown 
12d71ae5a4SJacob Faibussowitsch int main(int argc, char **args)
13d71ae5a4SJacob Faibussowitsch {
14c20d7725SJed Brown   Mat         A1, A2, A3, A4, A5, B, C, C1, nest;
15c4762a1bSJed Brown   Mat         mata[4];
16c4762a1bSJed Brown   Mat         aij;
17c4762a1bSJed Brown   MPI_Comm    comm;
18c4762a1bSJed Brown   PetscInt    m, M, n, istart, iend, ii, i, J, j, K = 10;
19c4762a1bSJed Brown   PetscScalar v;
20c4762a1bSJed Brown   PetscMPIInt size;
21c20d7725SJed Brown   PetscBool   equal;
22c4762a1bSJed Brown 
23327415f7SBarry Smith   PetscFunctionBeginUser;
249566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &args, (char *)0, help));
25c4762a1bSJed Brown   comm = PETSC_COMM_WORLD;
269566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(comm, &size));
27c20d7725SJed Brown 
28c4762a1bSJed Brown   /*
29c4762a1bSJed Brown      Assemble the matrix for the five point stencil, YET AGAIN
30c4762a1bSJed Brown   */
319566063dSJacob Faibussowitsch   PetscCall(MatCreate(comm, &A1));
32c4762a1bSJed Brown   m = 2, n = 2;
339566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(A1, PETSC_DECIDE, PETSC_DECIDE, m * n, m * n));
349566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(A1));
359566063dSJacob Faibussowitsch   PetscCall(MatSetUp(A1));
369566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRange(A1, &istart, &iend));
37c4762a1bSJed Brown   for (ii = istart; ii < iend; ii++) {
389371c9d4SSatish Balay     v = -1.0;
399371c9d4SSatish Balay     i = ii / n;
409371c9d4SSatish Balay     j = ii - i * n;
419371c9d4SSatish Balay     if (i > 0) {
429371c9d4SSatish Balay       J = ii - n;
439371c9d4SSatish Balay       PetscCall(MatSetValues(A1, 1, &ii, 1, &J, &v, INSERT_VALUES));
449371c9d4SSatish Balay     }
459371c9d4SSatish Balay     if (i < m - 1) {
469371c9d4SSatish Balay       J = ii + n;
479371c9d4SSatish Balay       PetscCall(MatSetValues(A1, 1, &ii, 1, &J, &v, INSERT_VALUES));
489371c9d4SSatish Balay     }
499371c9d4SSatish Balay     if (j > 0) {
509371c9d4SSatish Balay       J = ii - 1;
519371c9d4SSatish Balay       PetscCall(MatSetValues(A1, 1, &ii, 1, &J, &v, INSERT_VALUES));
529371c9d4SSatish Balay     }
539371c9d4SSatish Balay     if (j < n - 1) {
549371c9d4SSatish Balay       J = ii + 1;
559371c9d4SSatish Balay       PetscCall(MatSetValues(A1, 1, &ii, 1, &J, &v, INSERT_VALUES));
569371c9d4SSatish Balay     }
579371c9d4SSatish Balay     v = 4.0;
589371c9d4SSatish Balay     PetscCall(MatSetValues(A1, 1, &ii, 1, &ii, &v, INSERT_VALUES));
59c4762a1bSJed Brown   }
609566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A1, MAT_FINAL_ASSEMBLY));
619566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A1, MAT_FINAL_ASSEMBLY));
629566063dSJacob Faibussowitsch   PetscCall(MatView(A1, PETSC_VIEWER_STDOUT_WORLD));
63c20d7725SJed Brown 
649566063dSJacob Faibussowitsch   PetscCall(MatDuplicate(A1, MAT_COPY_VALUES, &A2));
659566063dSJacob Faibussowitsch   PetscCall(MatDuplicate(A1, MAT_COPY_VALUES, &A3));
669566063dSJacob Faibussowitsch   PetscCall(MatDuplicate(A1, MAT_COPY_VALUES, &A4));
67c20d7725SJed Brown 
68c4762a1bSJed Brown   /*create a nest matrix */
699566063dSJacob Faibussowitsch   PetscCall(MatCreate(comm, &nest));
709566063dSJacob Faibussowitsch   PetscCall(MatSetType(nest, MATNEST));
71c4762a1bSJed Brown   mata[0] = A1, mata[1] = A2, mata[2] = A3, mata[3] = A4;
729566063dSJacob Faibussowitsch   PetscCall(MatNestSetSubMats(nest, 2, NULL, 2, NULL, mata));
739566063dSJacob Faibussowitsch   PetscCall(MatSetUp(nest));
749566063dSJacob Faibussowitsch   PetscCall(MatConvert(nest, MATAIJ, MAT_INITIAL_MATRIX, &aij));
759566063dSJacob Faibussowitsch   PetscCall(MatView(aij, PETSC_VIEWER_STDOUT_WORLD));
76c20d7725SJed Brown 
77c20d7725SJed Brown   /* create a dense matrix */
789566063dSJacob Faibussowitsch   PetscCall(MatGetSize(nest, &M, NULL));
799566063dSJacob Faibussowitsch   PetscCall(MatGetLocalSize(nest, &m, NULL));
809566063dSJacob Faibussowitsch   PetscCall(MatCreateDense(comm, m, PETSC_DECIDE, M, K, NULL, &B));
81*f3fa974cSJacob Faibussowitsch   PetscCall(MatSetRandom(B, NULL));
82c20d7725SJed Brown 
83c20d7725SJed Brown   /* C = nest*B_dense */
849566063dSJacob Faibussowitsch   PetscCall(MatMatMult(nest, B, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &C));
859566063dSJacob Faibussowitsch   PetscCall(MatMatMult(nest, B, MAT_REUSE_MATRIX, PETSC_DEFAULT, &C));
869566063dSJacob Faibussowitsch   PetscCall(MatMatMultEqual(nest, B, C, 10, &equal));
8728b400f6SJacob Faibussowitsch   PetscCheck(equal, PETSC_COMM_WORLD, PETSC_ERR_PLIB, "Error in C != nest*B_dense");
883fbca975SStefano Zampini 
893fbca975SStefano Zampini   /* Test B = nest*C, reuse C and B with MatProductCreateWithMat() */
903fbca975SStefano Zampini   /* C has been obtained from nest*B. Clear internal data structures related to factors to prevent circular references */
919566063dSJacob Faibussowitsch   PetscCall(MatProductClear(C));
929566063dSJacob Faibussowitsch   PetscCall(MatProductCreateWithMat(nest, C, NULL, B));
939566063dSJacob Faibussowitsch   PetscCall(MatProductSetType(B, MATPRODUCT_AB));
949566063dSJacob Faibussowitsch   PetscCall(MatProductSetFromOptions(B));
959566063dSJacob Faibussowitsch   PetscCall(MatProductSymbolic(B));
969566063dSJacob Faibussowitsch   PetscCall(MatProductNumeric(B));
979566063dSJacob Faibussowitsch   PetscCall(MatMatMultEqual(nest, C, B, 10, &equal));
9828b400f6SJacob Faibussowitsch   PetscCheck(equal, PETSC_COMM_WORLD, PETSC_ERR_PLIB, "Error in B != nest*C_dense");
999566063dSJacob Faibussowitsch   PetscCall(MatConvert(nest, MATAIJ, MAT_INPLACE_MATRIX, &nest));
1009566063dSJacob Faibussowitsch   PetscCall(MatEqual(nest, aij, &equal));
10128b400f6SJacob Faibussowitsch   PetscCheck(equal, PETSC_COMM_WORLD, PETSC_ERR_PLIB, "Error in aij != nest");
1029566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&nest));
103c20d7725SJed Brown 
104c20d7725SJed Brown   if (size > 1) {                           /* Do not know why this test fails for size = 1 */
1059566063dSJacob Faibussowitsch     PetscCall(MatCreateTranspose(A1, &A5)); /* A1 is symmetric */
106c20d7725SJed Brown     mata[0] = A5;
1079566063dSJacob Faibussowitsch     PetscCall(MatCreate(comm, &nest));
1089566063dSJacob Faibussowitsch     PetscCall(MatSetType(nest, MATNEST));
1099566063dSJacob Faibussowitsch     PetscCall(MatNestSetSubMats(nest, 2, NULL, 2, NULL, mata));
1109566063dSJacob Faibussowitsch     PetscCall(MatSetUp(nest));
1119566063dSJacob Faibussowitsch     PetscCall(MatMatMult(nest, B, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &C1));
1129566063dSJacob Faibussowitsch     PetscCall(MatMatMult(nest, B, MAT_REUSE_MATRIX, PETSC_DEFAULT, &C1));
113c20d7725SJed Brown 
1149566063dSJacob Faibussowitsch     PetscCall(MatMatMultEqual(nest, B, C1, 10, &equal));
11528b400f6SJacob Faibussowitsch     PetscCheck(equal, PETSC_COMM_WORLD, PETSC_ERR_PLIB, "Error in C1 != C");
1169566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&C1));
1179566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&A5));
1189566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&nest));
119c20d7725SJed Brown   }
120c20d7725SJed Brown 
1219566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&C));
1229566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&B));
1239566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&aij));
1249566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&A1));
1259566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&A2));
1269566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&A3));
1279566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&A4));
1289566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
129b122ec5aSJacob Faibussowitsch   return 0;
130c4762a1bSJed Brown }
131c4762a1bSJed Brown 
132c4762a1bSJed Brown /*TEST
133c4762a1bSJed Brown 
134c4762a1bSJed Brown    test:
135c4762a1bSJed Brown       nsize: 2
136c4762a1bSJed Brown 
1373fbca975SStefano Zampini    test:
1383fbca975SStefano Zampini       suffix: 2
1393fbca975SStefano Zampini 
140c4762a1bSJed Brown TEST*/
141