xref: /petsc/src/mat/tests/ex58.c (revision d71ae5a4db6382e7f06317b8d368875286fe9008)
1c4762a1bSJed Brown 
2c4762a1bSJed Brown static char help[] = "Tests MatTranspose() and MatEqual() for MPIAIJ matrices.\n\n";
3c4762a1bSJed Brown 
4c4762a1bSJed Brown #include <petscmat.h>
5c4762a1bSJed Brown 
6*d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
7*d71ae5a4SJacob Faibussowitsch {
8c4762a1bSJed Brown   Mat         A, B;
9c4762a1bSJed Brown   PetscInt    m = 7, n, i, rstart, rend, cols[3];
10c4762a1bSJed Brown   PetscScalar v[3];
11c4762a1bSJed Brown   PetscBool   equal;
12c4762a1bSJed Brown   const char *eq[2];
13c4762a1bSJed Brown 
14327415f7SBarry Smith   PetscFunctionBeginUser;
159566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, (char *)0, help));
169566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-m", &m, NULL));
17c4762a1bSJed Brown   n = m;
18c4762a1bSJed Brown 
19c4762a1bSJed Brown   /* ------- Assemble matrix, --------- */
20c4762a1bSJed Brown 
219566063dSJacob Faibussowitsch   PetscCall(MatCreateAIJ(PETSC_COMM_WORLD, PETSC_DECIDE, PETSC_DECIDE, m, n, 0, 0, 0, 0, &A));
229566063dSJacob Faibussowitsch   PetscCall(MatSetOption(A, MAT_NEW_NONZERO_LOCATION_ERR, PETSC_FALSE));
239566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRange(A, &rstart, &rend));
24c4762a1bSJed Brown   if (!rstart) {
25c4762a1bSJed Brown     cols[0] = 0;
26c4762a1bSJed Brown     cols[1] = 1;
279371c9d4SSatish Balay     v[0]    = 2.0;
289371c9d4SSatish Balay     v[1]    = -1.0;
299566063dSJacob Faibussowitsch     PetscCall(MatSetValues(A, 1, &rstart, 2, cols, v, INSERT_VALUES));
30c4762a1bSJed Brown     rstart++;
31c4762a1bSJed Brown   }
32c4762a1bSJed Brown   if (rend == m) {
33c4762a1bSJed Brown     rend--;
34c4762a1bSJed Brown     cols[0] = rend - 1;
35c4762a1bSJed Brown     cols[1] = rend;
369371c9d4SSatish Balay     v[0]    = -1.0;
379371c9d4SSatish Balay     v[1]    = 2.0;
389566063dSJacob Faibussowitsch     PetscCall(MatSetValues(A, 1, &rend, 2, cols, v, INSERT_VALUES));
39c4762a1bSJed Brown   }
409371c9d4SSatish Balay   v[0] = -1.0;
419371c9d4SSatish Balay   v[1] = 2.0;
429371c9d4SSatish Balay   v[2] = -1.0;
43c4762a1bSJed Brown   for (i = rstart; i < rend; i++) {
44c4762a1bSJed Brown     cols[0] = i - 1;
45c4762a1bSJed Brown     cols[1] = i;
46c4762a1bSJed Brown     cols[2] = i + 1;
479566063dSJacob Faibussowitsch     PetscCall(MatSetValues(A, 1, &i, 3, cols, v, INSERT_VALUES));
48c4762a1bSJed Brown   }
499566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
509566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
51c4762a1bSJed Brown 
529566063dSJacob Faibussowitsch   PetscCall(MatTranspose(A, MAT_INITIAL_MATRIX, &B));
53c4762a1bSJed Brown 
549566063dSJacob Faibussowitsch   PetscCall(MatEqual(A, B, &equal));
55c4762a1bSJed Brown 
56c4762a1bSJed Brown   eq[0] = "not equal";
57c4762a1bSJed Brown   eq[1] = "equal";
589566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Matrices are %s\n", eq[equal]));
59c4762a1bSJed Brown 
609566063dSJacob Faibussowitsch   PetscCall(MatTranspose(A, MAT_REUSE_MATRIX, &B));
619566063dSJacob Faibussowitsch   PetscCall(MatEqual(A, B, &equal));
629566063dSJacob Faibussowitsch   if (!equal) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "MatTranspose with MAT_REUSE_MATRIX failed"));
63c4762a1bSJed Brown 
64c4762a1bSJed Brown   /* Free data structures */
659566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&A));
669566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&B));
67c4762a1bSJed Brown 
689566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
69b122ec5aSJacob Faibussowitsch   return 0;
70c4762a1bSJed Brown }
71c4762a1bSJed Brown 
72c4762a1bSJed Brown /*TEST
73c4762a1bSJed Brown 
74c4762a1bSJed Brown     test:
75c4762a1bSJed Brown 
76c4762a1bSJed Brown     test:
77c4762a1bSJed Brown       suffix: 2
78c4762a1bSJed Brown       nsize: 2
79c4762a1bSJed Brown       output_file: output/ex58_1.out
80c4762a1bSJed Brown 
81c4762a1bSJed Brown TEST*/
82