Actual source code: ex58.c
petsc-3.8.4 2018-03-24
2: static char help[] = "Tests MatTranspose() and MatEqual() for MPIAIJ matrices.\n\n";
4: #include <petscmat.h>
7: int main(int argc,char **argv)
8: {
9: Mat A,B;
10: PetscInt m = 7,n,i,rstart,rend,cols[3];
12: PetscScalar v[3];
13: PetscBool equal;
14: const char *eq[2];
16: PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
17: PetscOptionsGetInt(NULL,NULL,"-m",&m,NULL);
18: n = m;
20: /* ------- Assemble matrix, --------- */
22: MatCreateAIJ(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,m,n,0,0,0,0,&A);
23: MatSetOption(A,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_FALSE);
24: MatGetOwnershipRange(A,&rstart,&rend);
25: if (!rstart) {
26: cols[0] = 0;
27: cols[1] = 1;
28: v[0] = 2.0; v[1] = -1.0;
29: MatSetValues(A,1,&rstart,2,cols,v,INSERT_VALUES);
30: rstart++;
31: }
32: if (rend == m) {
33: rend--;
34: cols[0] = rend-1;
35: cols[1] = rend;
36: v[0] = -1.0; v[1] = 2.0;
37: MatSetValues(A,1,&rend,2,cols,v,INSERT_VALUES);
38: }
39: v[0] = -1.0; v[1] = 2.0; v[2] = -1.0;
40: for (i=rstart; i<rend; i++) {
41: cols[0] = i-1;
42: cols[1] = i;
43: cols[2] = i+1;
44: MatSetValues(A,1,&i,3,cols,v,INSERT_VALUES);
45: }
46: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
47: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
49: MatTranspose(A,MAT_INITIAL_MATRIX,&B);
51: MatEqual(A,B,&equal);
53: eq[0] = "not equal";
54: eq[1] = "equal";
55: PetscPrintf(PETSC_COMM_WORLD,"Matrices are %s\n",eq[equal]);
57: /* Free data structures */
58: MatDestroy(&A);
59: MatDestroy(&B);
62: PetscFinalize();
63: return ierr;
64: }