xref: /petsc/src/mat/tests/ex58.c (revision c4762a1b19cd2af06abeed90e8f9d34fb975dd94)
1*c4762a1bSJed Brown 
2*c4762a1bSJed Brown static char help[] = "Tests MatTranspose() and MatEqual() for MPIAIJ matrices.\n\n";
3*c4762a1bSJed Brown 
4*c4762a1bSJed Brown #include <petscmat.h>
5*c4762a1bSJed Brown 
6*c4762a1bSJed Brown int main(int argc,char **argv)
7*c4762a1bSJed Brown {
8*c4762a1bSJed Brown   Mat            A,B;
9*c4762a1bSJed Brown   PetscInt       m = 7,n,i,rstart,rend,cols[3];
10*c4762a1bSJed Brown   PetscErrorCode ierr;
11*c4762a1bSJed Brown   PetscScalar    v[3];
12*c4762a1bSJed Brown   PetscBool      equal;
13*c4762a1bSJed Brown   const char     *eq[2];
14*c4762a1bSJed Brown 
15*c4762a1bSJed Brown   ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
16*c4762a1bSJed Brown   ierr = PetscOptionsGetInt(NULL,NULL,"-m",&m,NULL);CHKERRQ(ierr);
17*c4762a1bSJed Brown   n    = m;
18*c4762a1bSJed Brown 
19*c4762a1bSJed Brown   /* ------- Assemble matrix, --------- */
20*c4762a1bSJed Brown 
21*c4762a1bSJed Brown   ierr = MatCreateAIJ(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,m,n,0,0,0,0,&A);CHKERRQ(ierr);
22*c4762a1bSJed Brown   ierr = MatSetOption(A,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_FALSE);CHKERRQ(ierr);
23*c4762a1bSJed Brown   ierr = MatGetOwnershipRange(A,&rstart,&rend);CHKERRQ(ierr);
24*c4762a1bSJed Brown   if (!rstart) {
25*c4762a1bSJed Brown     cols[0] = 0;
26*c4762a1bSJed Brown     cols[1] = 1;
27*c4762a1bSJed Brown     v[0]    = 2.0; v[1] = -1.0;
28*c4762a1bSJed Brown     ierr    = MatSetValues(A,1,&rstart,2,cols,v,INSERT_VALUES);CHKERRQ(ierr);
29*c4762a1bSJed Brown     rstart++;
30*c4762a1bSJed Brown   }
31*c4762a1bSJed Brown   if (rend == m) {
32*c4762a1bSJed Brown     rend--;
33*c4762a1bSJed Brown     cols[0] = rend-1;
34*c4762a1bSJed Brown     cols[1] = rend;
35*c4762a1bSJed Brown     v[0]    = -1.0; v[1] = 2.0;
36*c4762a1bSJed Brown     ierr    = MatSetValues(A,1,&rend,2,cols,v,INSERT_VALUES);CHKERRQ(ierr);
37*c4762a1bSJed Brown   }
38*c4762a1bSJed Brown   v[0] = -1.0; v[1] = 2.0; v[2] = -1.0;
39*c4762a1bSJed Brown   for (i=rstart; i<rend; i++) {
40*c4762a1bSJed Brown     cols[0] = i-1;
41*c4762a1bSJed Brown     cols[1] = i;
42*c4762a1bSJed Brown     cols[2] = i+1;
43*c4762a1bSJed Brown     ierr    = MatSetValues(A,1,&i,3,cols,v,INSERT_VALUES);CHKERRQ(ierr);
44*c4762a1bSJed Brown   }
45*c4762a1bSJed Brown   ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
46*c4762a1bSJed Brown   ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
47*c4762a1bSJed Brown 
48*c4762a1bSJed Brown   ierr = MatTranspose(A,MAT_INITIAL_MATRIX,&B);CHKERRQ(ierr);
49*c4762a1bSJed Brown 
50*c4762a1bSJed Brown   ierr = MatEqual(A,B,&equal);CHKERRQ(ierr);
51*c4762a1bSJed Brown 
52*c4762a1bSJed Brown   eq[0] = "not equal";
53*c4762a1bSJed Brown   eq[1] = "equal";
54*c4762a1bSJed Brown   ierr  = PetscPrintf(PETSC_COMM_WORLD,"Matrices are %s\n",eq[equal]);CHKERRQ(ierr);
55*c4762a1bSJed Brown 
56*c4762a1bSJed Brown   ierr = MatTranspose(A,MAT_REUSE_MATRIX,&B);CHKERRQ(ierr);
57*c4762a1bSJed Brown   ierr = MatEqual(A,B,&equal);CHKERRQ(ierr);
58*c4762a1bSJed Brown   if (!equal) { ierr = PetscPrintf(PETSC_COMM_WORLD,"MatTranspose with MAT_REUSE_MATRIX failed");CHKERRQ(ierr); }
59*c4762a1bSJed Brown 
60*c4762a1bSJed Brown   /* Free data structures */
61*c4762a1bSJed Brown   ierr = MatDestroy(&A);CHKERRQ(ierr);
62*c4762a1bSJed Brown   ierr = MatDestroy(&B);CHKERRQ(ierr);
63*c4762a1bSJed Brown 
64*c4762a1bSJed Brown 
65*c4762a1bSJed Brown   ierr = PetscFinalize();
66*c4762a1bSJed Brown   return ierr;
67*c4762a1bSJed Brown }
68*c4762a1bSJed Brown 
69*c4762a1bSJed Brown /*TEST
70*c4762a1bSJed Brown 
71*c4762a1bSJed Brown     test:
72*c4762a1bSJed Brown 
73*c4762a1bSJed Brown     test:
74*c4762a1bSJed Brown       suffix: 2
75*c4762a1bSJed Brown       nsize: 2
76*c4762a1bSJed Brown       output_file: output/ex58_1.out
77*c4762a1bSJed Brown 
78*c4762a1bSJed Brown TEST*/
79