xref: /petsc/src/mat/tests/ex39.c (revision 609caa7c8c030312b00807b4f015fd827bb80932)
1c4762a1bSJed Brown static char help[] = "Tests Elemental interface.\n\n";
2c4762a1bSJed Brown 
3c4762a1bSJed Brown #include <petscmat.h>
4c4762a1bSJed Brown 
main(int argc,char ** args)5d71ae5a4SJacob Faibussowitsch int main(int argc, char **args)
6d71ae5a4SJacob Faibussowitsch {
7c4762a1bSJed Brown   Mat             Cdense, B, C, Ct;
8c4762a1bSJed Brown   Vec             d;
9c4762a1bSJed Brown   PetscInt        i, j, m = 5, n, nrows, ncols;
10c4762a1bSJed Brown   const PetscInt *rows, *cols;
11c4762a1bSJed Brown   IS              isrows, iscols;
12c4762a1bSJed Brown   PetscScalar    *v;
13c4762a1bSJed Brown   PetscMPIInt     rank, size;
14c4762a1bSJed Brown   PetscReal       Cnorm;
15c4762a1bSJed Brown   PetscBool       flg, mats_view = PETSC_FALSE;
16c4762a1bSJed Brown 
17327415f7SBarry Smith   PetscFunctionBeginUser;
18c8025a54SPierre Jolivet   PetscCall(PetscInitialize(&argc, &args, NULL, help));
199566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
209566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
219566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-m", &m, NULL));
22c4762a1bSJed Brown   n = m;
239566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-n", &n, NULL));
249566063dSJacob Faibussowitsch   PetscCall(PetscOptionsHasName(NULL, NULL, "-mats_view", &mats_view));
25c4762a1bSJed Brown 
269566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD, &C));
279566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(C, m, n, PETSC_DECIDE, PETSC_DECIDE));
289566063dSJacob Faibussowitsch   PetscCall(MatSetType(C, MATELEMENTAL));
299566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(C));
309566063dSJacob Faibussowitsch   PetscCall(MatSetUp(C));
319566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipIS(C, &isrows, &iscols));
329566063dSJacob Faibussowitsch   PetscCall(ISGetLocalSize(isrows, &nrows));
339566063dSJacob Faibussowitsch   PetscCall(ISGetIndices(isrows, &rows));
349566063dSJacob Faibussowitsch   PetscCall(ISGetLocalSize(iscols, &ncols));
359566063dSJacob Faibussowitsch   PetscCall(ISGetIndices(iscols, &cols));
369566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nrows * ncols, &v));
37c4762a1bSJed Brown #if defined(PETSC_USE_COMPLEX)
38c4762a1bSJed Brown   PetscRandom rand;
39c4762a1bSJed Brown   PetscScalar rval;
409566063dSJacob Faibussowitsch   PetscCall(PetscRandomCreate(PETSC_COMM_WORLD, &rand));
419566063dSJacob Faibussowitsch   PetscCall(PetscRandomSetFromOptions(rand));
42c4762a1bSJed Brown   for (i = 0; i < nrows; i++) {
43c4762a1bSJed Brown     for (j = 0; j < ncols; j++) {
449566063dSJacob Faibussowitsch       PetscCall(PetscRandomGetValue(rand, &rval));
45c4762a1bSJed Brown       v[i * ncols + j] = rval;
46c4762a1bSJed Brown     }
47c4762a1bSJed Brown   }
489566063dSJacob Faibussowitsch   PetscCall(PetscRandomDestroy(&rand));
49c4762a1bSJed Brown #else
50c4762a1bSJed Brown   for (i = 0; i < nrows; i++) {
51ad540459SPierre Jolivet     for (j = 0; j < ncols; j++) v[i * ncols + j] = (PetscReal)(10000 * rank + 100 * rows[i] + cols[j]);
52c4762a1bSJed Brown   }
53c4762a1bSJed Brown #endif
549566063dSJacob Faibussowitsch   PetscCall(MatSetValues(C, nrows, rows, ncols, cols, v, INSERT_VALUES));
559566063dSJacob Faibussowitsch   PetscCall(ISRestoreIndices(isrows, &rows));
569566063dSJacob Faibussowitsch   PetscCall(ISRestoreIndices(iscols, &cols));
579566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(C, MAT_FINAL_ASSEMBLY));
589566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(C, MAT_FINAL_ASSEMBLY));
599566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&isrows));
609566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&iscols));
61c4762a1bSJed Brown 
62c4762a1bSJed Brown   /* Test MatView(), MatDuplicate() and out-of-place MatConvert() */
639566063dSJacob Faibussowitsch   PetscCall(MatDuplicate(C, MAT_COPY_VALUES, &B));
64c4762a1bSJed Brown   if (mats_view) {
659566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Duplicated C:\n"));
669566063dSJacob Faibussowitsch     PetscCall(MatView(B, PETSC_VIEWER_STDOUT_WORLD));
67c4762a1bSJed Brown   }
689566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&B));
699566063dSJacob Faibussowitsch   PetscCall(MatConvert(C, MATMPIDENSE, MAT_INITIAL_MATRIX, &Cdense));
709566063dSJacob Faibussowitsch   PetscCall(MatMultEqual(C, Cdense, 5, &flg));
7128b400f6SJacob Faibussowitsch   PetscCheck(flg, PETSC_COMM_SELF, PETSC_ERR_ARG_NOTSAMETYPE, "Cdense != C. MatConvert() fails");
72c4762a1bSJed Brown 
73c4762a1bSJed Brown   /* Test MatNorm() */
749566063dSJacob Faibussowitsch   PetscCall(MatNorm(C, NORM_1, &Cnorm));
75c4762a1bSJed Brown 
76c4762a1bSJed Brown   /* Test MatTranspose(), MatZeroEntries() and MatGetDiagonal() */
779566063dSJacob Faibussowitsch   PetscCall(MatTranspose(C, MAT_INITIAL_MATRIX, &Ct));
789566063dSJacob Faibussowitsch   PetscCall(MatConjugate(Ct));
79c4762a1bSJed Brown   if (mats_view) {
809566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "C's Transpose Conjugate:\n"));
819566063dSJacob Faibussowitsch     PetscCall(MatView(Ct, PETSC_VIEWER_STDOUT_WORLD));
82c4762a1bSJed Brown   }
839566063dSJacob Faibussowitsch   PetscCall(MatZeroEntries(Ct));
849566063dSJacob Faibussowitsch   PetscCall(VecCreate(PETSC_COMM_WORLD, &d));
859566063dSJacob Faibussowitsch   PetscCall(VecSetSizes(d, m > n ? n : m, PETSC_DECIDE));
869566063dSJacob Faibussowitsch   PetscCall(VecSetFromOptions(d));
879566063dSJacob Faibussowitsch   PetscCall(MatGetDiagonal(C, d));
88c4762a1bSJed Brown   if (mats_view) {
899566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Diagonal of C:\n"));
909566063dSJacob Faibussowitsch     PetscCall(VecView(d, PETSC_VIEWER_STDOUT_WORLD));
91c4762a1bSJed Brown   }
92c4762a1bSJed Brown   if (m > n) {
939566063dSJacob Faibussowitsch     PetscCall(MatDiagonalScale(C, NULL, d));
94c4762a1bSJed Brown   } else {
959566063dSJacob Faibussowitsch     PetscCall(MatDiagonalScale(C, d, NULL));
96c4762a1bSJed Brown   }
97c4762a1bSJed Brown   if (mats_view) {
989566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Diagonal Scaled C:\n"));
999566063dSJacob Faibussowitsch     PetscCall(MatView(C, PETSC_VIEWER_STDOUT_WORLD));
100c4762a1bSJed Brown   }
101c4762a1bSJed Brown 
102c4762a1bSJed Brown   /* Test MatAXPY(), MatAYPX() and in-place MatConvert() */
1039566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD, &B));
1049566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(B, m, n, PETSC_DECIDE, PETSC_DECIDE));
1059566063dSJacob Faibussowitsch   PetscCall(MatSetType(B, MATELEMENTAL));
1069566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(B));
1079566063dSJacob Faibussowitsch   PetscCall(MatSetUp(B));
1089566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipIS(B, &isrows, &iscols));
1099566063dSJacob Faibussowitsch   PetscCall(ISGetLocalSize(isrows, &nrows));
1109566063dSJacob Faibussowitsch   PetscCall(ISGetIndices(isrows, &rows));
1119566063dSJacob Faibussowitsch   PetscCall(ISGetLocalSize(iscols, &ncols));
1129566063dSJacob Faibussowitsch   PetscCall(ISGetIndices(iscols, &cols));
113c4762a1bSJed Brown   for (i = 0; i < nrows; i++) {
114ad540459SPierre Jolivet     for (j = 0; j < ncols; j++) v[i * ncols + j] = (PetscReal)(1000 * rows[i] + cols[j]);
115c4762a1bSJed Brown   }
1169566063dSJacob Faibussowitsch   PetscCall(MatSetValues(B, nrows, rows, ncols, cols, v, INSERT_VALUES));
1179566063dSJacob Faibussowitsch   PetscCall(PetscFree(v));
1189566063dSJacob Faibussowitsch   PetscCall(ISRestoreIndices(isrows, &rows));
1199566063dSJacob Faibussowitsch   PetscCall(ISRestoreIndices(iscols, &cols));
1209566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY));
1219566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY));
1229566063dSJacob Faibussowitsch   PetscCall(MatAXPY(B, 2.5, C, SAME_NONZERO_PATTERN));
1239566063dSJacob Faibussowitsch   PetscCall(MatAYPX(B, 3.75, C, SAME_NONZERO_PATTERN));
1249566063dSJacob Faibussowitsch   PetscCall(MatConvert(B, MATDENSE, MAT_INPLACE_MATRIX, &B));
125c4762a1bSJed Brown   if (mats_view) {
1269566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "B after MatAXPY and MatAYPX:\n"));
1279566063dSJacob Faibussowitsch     PetscCall(MatView(B, PETSC_VIEWER_STDOUT_WORLD));
128c4762a1bSJed Brown   }
1299566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&isrows));
1309566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&iscols));
1319566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&B));
132c4762a1bSJed Brown 
133c4762a1bSJed Brown   /* Test MatMatTransposeMult(): B = C*C^T */
134fb842aefSJose E. Roman   PetscCall(MatMatTransposeMult(C, C, MAT_INITIAL_MATRIX, PETSC_DETERMINE, &B));
1359566063dSJacob Faibussowitsch   PetscCall(MatScale(C, 2.0));
136fb842aefSJose E. Roman   PetscCall(MatMatTransposeMult(C, C, MAT_REUSE_MATRIX, PETSC_DETERMINE, &B));
1379566063dSJacob Faibussowitsch   PetscCall(MatMatTransposeMultEqual(C, C, B, 10, &flg));
13828b400f6SJacob Faibussowitsch   PetscCheck(flg, PETSC_COMM_WORLD, PETSC_ERR_ARG_INCOMP, "MatMatTransposeMult: B != C*B^T");
139c20d7725SJed Brown 
140c4762a1bSJed Brown   if (mats_view) {
1419566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "C MatMatTransposeMult C:\n"));
1429566063dSJacob Faibussowitsch     PetscCall(MatView(B, PETSC_VIEWER_STDOUT_WORLD));
143c4762a1bSJed Brown   }
144c4762a1bSJed Brown 
1459566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&Cdense));
1469566063dSJacob Faibussowitsch   PetscCall(PetscFree(v));
1479566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&B));
1489566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&C));
1499566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&Ct));
1509566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&d));
1519566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
152b122ec5aSJacob Faibussowitsch   return 0;
153c4762a1bSJed Brown }
154c4762a1bSJed Brown 
155c4762a1bSJed Brown /*TEST
156c4762a1bSJed Brown 
157c4762a1bSJed Brown    test:
158c4762a1bSJed Brown       nsize: 2
159c4762a1bSJed Brown       args: -m 3 -n 2
160c4762a1bSJed Brown       requires: elemental
161*3886731fSPierre Jolivet       output_file: output/empty.out
162c4762a1bSJed Brown 
163c4762a1bSJed Brown    test:
164c4762a1bSJed Brown       suffix: 2
165c4762a1bSJed Brown       nsize: 6
166c4762a1bSJed Brown       args: -m 2 -n 3
167c4762a1bSJed Brown       requires: elemental
168*3886731fSPierre Jolivet       output_file: output/empty.out
169c4762a1bSJed Brown 
170c20d7725SJed Brown    test:
171c20d7725SJed Brown       suffix: 3
172c20d7725SJed Brown       nsize: 1
173c20d7725SJed Brown       args: -m 2 -n 3
174c20d7725SJed Brown       requires: elemental
175*3886731fSPierre Jolivet       output_file: output/empty.out
176c20d7725SJed Brown 
177c4762a1bSJed Brown TEST*/
178