xref: /petsc/src/mat/tests/ex25.c (revision 9371c9d470a9602b6d10a8bf50c9b2280a79e45a)
1c4762a1bSJed Brown 
2c4762a1bSJed Brown static char help[] = "Tests MatTranspose()\n\n";
3c4762a1bSJed Brown 
4c4762a1bSJed Brown #include <petscmat.h>
5c4762a1bSJed Brown 
6*9371c9d4SSatish Balay int main(int argc, char **args) {
7c4762a1bSJed Brown   Mat         C, A;
8c4762a1bSJed Brown   PetscScalar v;
9c4762a1bSJed Brown   PetscInt    i, j, m = 4, n = 4, Ii, J, Istart, Iend;
10c4762a1bSJed Brown   PetscMPIInt rank, size;
11c4762a1bSJed Brown   PetscBool   equal = PETSC_FALSE;
12c4762a1bSJed Brown 
13327415f7SBarry Smith   PetscFunctionBeginUser;
149566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &args, (char *)0, help));
159566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
169566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
17c4762a1bSJed Brown 
189566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-m", &m, NULL));
19c4762a1bSJed Brown   n = m;
20c4762a1bSJed Brown 
219566063dSJacob Faibussowitsch   PetscCall(MatCreateAIJ(PETSC_COMM_WORLD, PETSC_DECIDE, PETSC_DECIDE, m * n, m * n, 5, NULL, 5, NULL, &C));
22c4762a1bSJed Brown 
23c4762a1bSJed Brown   /* create the symmetric matrix for the five point stencil */
249566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRange(C, &Istart, &Iend));
25c4762a1bSJed Brown   for (Ii = Istart; Ii < Iend; Ii++) {
26*9371c9d4SSatish Balay     v = -1.0;
27*9371c9d4SSatish Balay     i = Ii / n;
28*9371c9d4SSatish Balay     j = Ii - i * n;
29*9371c9d4SSatish Balay     if (i > 0) {
30*9371c9d4SSatish Balay       J = Ii - n;
31*9371c9d4SSatish Balay       PetscCall(MatSetValues(C, 1, &Ii, 1, &J, &v, INSERT_VALUES));
32*9371c9d4SSatish Balay     }
33*9371c9d4SSatish Balay     if (i < m - 1) {
34*9371c9d4SSatish Balay       J = Ii + n;
35*9371c9d4SSatish Balay       PetscCall(MatSetValues(C, 1, &Ii, 1, &J, &v, INSERT_VALUES));
36*9371c9d4SSatish Balay     }
37*9371c9d4SSatish Balay     if (j > 0) {
38*9371c9d4SSatish Balay       J = Ii - 1;
39*9371c9d4SSatish Balay       PetscCall(MatSetValues(C, 1, &Ii, 1, &J, &v, INSERT_VALUES));
40*9371c9d4SSatish Balay     }
41*9371c9d4SSatish Balay     if (j < n - 1) {
42*9371c9d4SSatish Balay       J = Ii + 1;
43*9371c9d4SSatish Balay       PetscCall(MatSetValues(C, 1, &Ii, 1, &J, &v, INSERT_VALUES));
44*9371c9d4SSatish Balay     }
45*9371c9d4SSatish Balay     v = 4.0;
46*9371c9d4SSatish Balay     PetscCall(MatSetValues(C, 1, &Ii, 1, &Ii, &v, INSERT_VALUES));
47c4762a1bSJed Brown   }
489566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(C, MAT_FINAL_ASSEMBLY));
499566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(C, MAT_FINAL_ASSEMBLY));
50c4762a1bSJed Brown 
519566063dSJacob Faibussowitsch   PetscCall(MatTranspose(C, MAT_INITIAL_MATRIX, &A));
52c4762a1bSJed Brown 
539566063dSJacob Faibussowitsch   PetscCall(MatEqual(C, A, &equal));
5428b400f6SJacob Faibussowitsch   PetscCheck(equal, PETSC_COMM_SELF, PETSC_ERR_SUP, "C != C^T");
55c4762a1bSJed Brown 
569566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&C));
579566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&A));
589566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
59b122ec5aSJacob Faibussowitsch   return 0;
60c4762a1bSJed Brown }
61c4762a1bSJed Brown 
62c4762a1bSJed Brown /*TEST
63c4762a1bSJed Brown 
64c4762a1bSJed Brown    test:
65c4762a1bSJed Brown 
66c4762a1bSJed Brown TEST*/
67