Actual source code: ex47.c
petsc-3.10.5 2019-03-28
2: static char help[] = "Tests the various routines in MatBAIJ format.\n\
3: Input arguments are:\n\
4: -f <input_file> : file to load. For example see $PETSC_DIR/share/petsc/datafiles/matrices\n\n";
6: #include <petscmat.h>
8: int main(int argc,char **args)
9: {
10: Mat A,B,C;
11: PetscViewer va,vb,vc;
12: Vec x,y;
13: PetscErrorCode ierr;
14: PetscInt i,j,row,m,n,ncols1,ncols2,ct,m2,n2;
15: const PetscInt *cols1,*cols2;
16: char file[PETSC_MAX_PATH_LEN];
17: PetscBool tflg;
18: PetscScalar rval;
19: const PetscScalar *vals1,*vals2;
20: PetscReal norm1,norm2,rnorm;
21: PetscRandom r;
24: PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr;
25: #if defined(PETSC_USE_COMPLEX)
26: SETERRQ(PETSC_COMM_WORLD,1,"This example does not work with complex numbers");
27: #else
28: PetscOptionsGetString(NULL,NULL,"-f",file,PETSC_MAX_PATH_LEN,NULL);
30: /* Load the matrix as AIJ format */
31: PetscViewerBinaryOpen(PETSC_COMM_WORLD,file,FILE_MODE_READ,&va);
32: MatCreate(PETSC_COMM_WORLD,&A);
33: MatSetType(A,MATSEQAIJ);
34: MatLoad(A,va);
35: PetscViewerDestroy(&va);
37: /* Load the matrix as BAIJ format */
38: PetscViewerBinaryOpen(PETSC_COMM_WORLD,file,FILE_MODE_READ,&vb);
39: MatCreate(PETSC_COMM_WORLD,&B);
40: MatSetType(B,MATSEQBAIJ);
41: MatLoad(B,vb);
42: PetscViewerDestroy(&vb);
44: /* Load the matrix as BAIJ format */
45: PetscViewerBinaryOpen(PETSC_COMM_WORLD,file,FILE_MODE_READ,&vc);
46: MatCreate(PETSC_COMM_WORLD,&C);
47: MatSetType(C,MATSEQBAIJ);
48: MatSetFromOptions(C);
49: MatLoad(C,vc);
50: PetscViewerDestroy(&vc);
52: MatGetSize(A,&m,&n);
53: MatGetSize(B,&m2,&n2);
54: if (m!=m2) SETERRQ(PETSC_COMM_SELF,1,"Matrices are of different size. Cannot run this example");
56: /* Test MatEqual() */
57: MatEqual(B,C,&tflg);
58: if (!tflg) SETERRQ(PETSC_COMM_SELF,1,"MatEqual() failed");
60: /* Test MatGetDiagonal() */
61: VecCreateSeq(PETSC_COMM_SELF,m,&x);
62: VecCreateSeq(PETSC_COMM_SELF,m,&y);
64: MatGetDiagonal(A,x);
65: MatGetDiagonal(B,y);
67: VecEqual(x,y,&tflg);
68: if (!tflg) SETERRQ(PETSC_COMM_SELF,1,"MatGetDiagonal() failed");
70: /* Test MatDiagonalScale() */
71: PetscRandomCreate(PETSC_COMM_SELF,&r);
72: PetscRandomSetFromOptions(r);
73: VecSetRandom(x,r);
74: VecSetRandom(y,r);
76: MatDiagonalScale(A,x,y);
77: MatDiagonalScale(B,x,y);
78: MatMult(A,x,y);
79: VecNorm(y,NORM_2,&norm1);
80: MatMult(B,x,y);
81: VecNorm(y,NORM_2,&norm2);
82: rnorm = ((norm1-norm2)*100)/norm1;
83: if (rnorm<-0.1 || rnorm>0.01) {
84: SETERRQ2(PETSC_COMM_SELF,1,"MatDiagonalScale() failed Norm1 %g Norm2 %g",(double)norm1,(double)norm2);
85: }
87: /* Test MatGetRow()/ MatRestoreRow() */
88: for (ct=0; ct<100; ct++) {
89: PetscRandomGetValue(r,&rval);
90: row = (int)(rval*m);
91: MatGetRow(A,row,&ncols1,&cols1,&vals1);
92: MatGetRow(B,row,&ncols2,&cols2,&vals2);
94: for (i=0,j=0; i<ncols1 && j<ncols2; i++) {
95: while (cols2[j] != cols1[i]) j++;
96: if (vals1[i] != vals2[j]) SETERRQ(PETSC_COMM_SELF,1,"MatGetRow() failed - vals incorrect.");
97: }
98: if (i<ncols1) SETERRQ(PETSC_COMM_SELF,1,"MatGetRow() failed - cols incorrect");
100: MatRestoreRow(A,row,&ncols1,&cols1,&vals1);
101: MatRestoreRow(B,row,&ncols2,&cols2,&vals2);
102: }
104: MatDestroy(&A);
105: MatDestroy(&B);
106: MatDestroy(&C);
107: VecDestroy(&x);
108: VecDestroy(&y);
109: PetscRandomDestroy(&r);
110: #endif
111: PetscFinalize();
112: return ierr;
113: }
115: /*TEST
117: build:
118: requires: !complex
120: test:
121: args: -f ${DATAFILESPATH}/matrices/cfd.1.10 -mat_block_size 5
122: requires: !complex double datafilespath !define(PETSC_USE_64BIT_INDICES)
124: TEST*/