xref: /petsc/src/mat/tests/ex58.c (revision b122ec5aa1bd4469eb4e0673542fb7de3f411254)
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 
6c4762a1bSJed Brown int main(int argc,char **argv)
7c4762a1bSJed Brown {
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 
14*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscInitialize(&argc,&argv,(char*)0,help));
155f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-m",&m,NULL));
16c4762a1bSJed Brown   n    = m;
17c4762a1bSJed Brown 
18c4762a1bSJed Brown   /* ------- Assemble matrix, --------- */
19c4762a1bSJed Brown 
205f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreateAIJ(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,m,n,0,0,0,0,&A));
215f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetOption(A,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_FALSE));
225f80ce2aSJacob Faibussowitsch   CHKERRQ(MatGetOwnershipRange(A,&rstart,&rend));
23c4762a1bSJed Brown   if (!rstart) {
24c4762a1bSJed Brown     cols[0] = 0;
25c4762a1bSJed Brown     cols[1] = 1;
26c4762a1bSJed Brown     v[0]    = 2.0; v[1] = -1.0;
275f80ce2aSJacob Faibussowitsch     CHKERRQ(MatSetValues(A,1,&rstart,2,cols,v,INSERT_VALUES));
28c4762a1bSJed Brown     rstart++;
29c4762a1bSJed Brown   }
30c4762a1bSJed Brown   if (rend == m) {
31c4762a1bSJed Brown     rend--;
32c4762a1bSJed Brown     cols[0] = rend-1;
33c4762a1bSJed Brown     cols[1] = rend;
34c4762a1bSJed Brown     v[0]    = -1.0; v[1] = 2.0;
355f80ce2aSJacob Faibussowitsch     CHKERRQ(MatSetValues(A,1,&rend,2,cols,v,INSERT_VALUES));
36c4762a1bSJed Brown   }
37c4762a1bSJed Brown   v[0] = -1.0; v[1] = 2.0; v[2] = -1.0;
38c4762a1bSJed Brown   for (i=rstart; i<rend; i++) {
39c4762a1bSJed Brown     cols[0] = i-1;
40c4762a1bSJed Brown     cols[1] = i;
41c4762a1bSJed Brown     cols[2] = i+1;
425f80ce2aSJacob Faibussowitsch     CHKERRQ(MatSetValues(A,1,&i,3,cols,v,INSERT_VALUES));
43c4762a1bSJed Brown   }
445f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
455f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
46c4762a1bSJed Brown 
475f80ce2aSJacob Faibussowitsch   CHKERRQ(MatTranspose(A,MAT_INITIAL_MATRIX,&B));
48c4762a1bSJed Brown 
495f80ce2aSJacob Faibussowitsch   CHKERRQ(MatEqual(A,B,&equal));
50c4762a1bSJed Brown 
51c4762a1bSJed Brown   eq[0] = "not equal";
52c4762a1bSJed Brown   eq[1] = "equal";
535f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"Matrices are %s\n",eq[equal]));
54c4762a1bSJed Brown 
555f80ce2aSJacob Faibussowitsch   CHKERRQ(MatTranspose(A,MAT_REUSE_MATRIX,&B));
565f80ce2aSJacob Faibussowitsch   CHKERRQ(MatEqual(A,B,&equal));
575f80ce2aSJacob Faibussowitsch   if (!equal) CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"MatTranspose with MAT_REUSE_MATRIX failed"));
58c4762a1bSJed Brown 
59c4762a1bSJed Brown   /* Free data structures */
605f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&A));
615f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&B));
62c4762a1bSJed Brown 
63*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscFinalize());
64*b122ec5aSJacob Faibussowitsch   return 0;
65c4762a1bSJed Brown }
66c4762a1bSJed Brown 
67c4762a1bSJed Brown /*TEST
68c4762a1bSJed Brown 
69c4762a1bSJed Brown     test:
70c4762a1bSJed Brown 
71c4762a1bSJed Brown     test:
72c4762a1bSJed Brown       suffix: 2
73c4762a1bSJed Brown       nsize: 2
74c4762a1bSJed Brown       output_file: output/ex58_1.out
75c4762a1bSJed Brown 
76c4762a1bSJed Brown TEST*/
77