xref: /petsc/src/mat/tests/ex25.c (revision b122ec5aa1bd4469eb4e0673542fb7de3f411254)
1c4762a1bSJed Brown 
2c4762a1bSJed Brown static char help[] = "Tests MatTranspose()\n\n";
3c4762a1bSJed Brown 
4c4762a1bSJed Brown #include <petscmat.h>
5c4762a1bSJed Brown 
6c4762a1bSJed Brown int main(int argc,char **args)
7c4762a1bSJed Brown {
8c4762a1bSJed Brown   Mat            C,A;
9c4762a1bSJed Brown   PetscScalar    v;
10c4762a1bSJed Brown   PetscInt       i,j,m = 4,n = 4,Ii,J,Istart,Iend;
11c4762a1bSJed Brown   PetscMPIInt    rank,size;
12c4762a1bSJed Brown   PetscBool      equal=PETSC_FALSE;
13c4762a1bSJed Brown 
14*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscInitialize(&argc,&args,(char*)0,help));
155f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank));
165f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size));
17c4762a1bSJed Brown 
185f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-m",&m,NULL));
19c4762a1bSJed Brown   n    = m;
20c4762a1bSJed Brown 
215f80ce2aSJacob Faibussowitsch   CHKERRQ(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 */
245f80ce2aSJacob Faibussowitsch   CHKERRQ(MatGetOwnershipRange(C,&Istart,&Iend));
25c4762a1bSJed Brown   for (Ii=Istart; Ii<Iend; Ii++) {
26c4762a1bSJed Brown     v = -1.0; i = Ii/n; j = Ii - i*n;
275f80ce2aSJacob Faibussowitsch     if (i>0)   {J = Ii - n; CHKERRQ(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));}
285f80ce2aSJacob Faibussowitsch     if (i<m-1) {J = Ii + n; CHKERRQ(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));}
295f80ce2aSJacob Faibussowitsch     if (j>0)   {J = Ii - 1; CHKERRQ(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));}
305f80ce2aSJacob Faibussowitsch     if (j<n-1) {J = Ii + 1; CHKERRQ(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));}
315f80ce2aSJacob Faibussowitsch     v = 4.0; CHKERRQ(MatSetValues(C,1,&Ii,1,&Ii,&v,INSERT_VALUES));
32c4762a1bSJed Brown   }
335f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY));
345f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY));
35c4762a1bSJed Brown 
365f80ce2aSJacob Faibussowitsch   CHKERRQ(MatTranspose(C,MAT_INITIAL_MATRIX,&A));
37c4762a1bSJed Brown 
385f80ce2aSJacob Faibussowitsch   CHKERRQ(MatEqual(C,A,&equal));
3928b400f6SJacob Faibussowitsch   PetscCheck(equal,PETSC_COMM_SELF,PETSC_ERR_SUP,"C != C^T");
40c4762a1bSJed Brown 
415f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&C));
425f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&A));
43*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscFinalize());
44*b122ec5aSJacob Faibussowitsch   return 0;
45c4762a1bSJed Brown }
46c4762a1bSJed Brown 
47c4762a1bSJed Brown /*TEST
48c4762a1bSJed Brown 
49c4762a1bSJed Brown    test:
50c4762a1bSJed Brown 
51c4762a1bSJed Brown TEST*/
52