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