1c4762a1bSJed Brown static char help[] = "Tests the various routines in MatBAIJ format.\n\
2c4762a1bSJed Brown Input arguments are:\n\
3c4762a1bSJed Brown -f <input_file> : file to load. For example see $PETSC_DIR/share/petsc/datafiles/matrices\n\n";
4c4762a1bSJed Brown
5c4762a1bSJed Brown #include <petscmat.h>
6c4762a1bSJed Brown
main(int argc,char ** args)7d71ae5a4SJacob Faibussowitsch int main(int argc, char **args)
8d71ae5a4SJacob Faibussowitsch {
9c4762a1bSJed Brown Mat A, B, C;
10c4762a1bSJed Brown PetscViewer va, vb, vc;
11c4762a1bSJed Brown Vec x, y;
12c4762a1bSJed Brown PetscInt i, j, row, m, n, ncols1, ncols2, ct, m2, n2;
13c4762a1bSJed Brown const PetscInt *cols1, *cols2;
14c4762a1bSJed Brown char file[PETSC_MAX_PATH_LEN];
15c4762a1bSJed Brown PetscBool tflg;
16c4762a1bSJed Brown PetscScalar rval;
17c4762a1bSJed Brown const PetscScalar *vals1, *vals2;
18c4762a1bSJed Brown PetscReal norm1, norm2, rnorm;
19c4762a1bSJed Brown PetscRandom r;
20c4762a1bSJed Brown
21327415f7SBarry Smith PetscFunctionBeginUser;
22c8025a54SPierre Jolivet PetscCall(PetscInitialize(&argc, &args, NULL, help));
239566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetString(NULL, NULL, "-f", file, sizeof(file), NULL));
24c4762a1bSJed Brown
25c4762a1bSJed Brown /* Load the matrix as AIJ format */
269566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD, file, FILE_MODE_READ, &va));
279566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &A));
289566063dSJacob Faibussowitsch PetscCall(MatSetType(A, MATSEQAIJ));
299566063dSJacob Faibussowitsch PetscCall(MatLoad(A, va));
309566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&va));
31c4762a1bSJed Brown
32c4762a1bSJed Brown /* Load the matrix as BAIJ format */
339566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD, file, FILE_MODE_READ, &vb));
349566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &B));
359566063dSJacob Faibussowitsch PetscCall(MatSetType(B, MATSEQBAIJ));
369566063dSJacob Faibussowitsch PetscCall(MatLoad(B, vb));
379566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&vb));
38c4762a1bSJed Brown
39c4762a1bSJed Brown /* Load the matrix as BAIJ format */
409566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD, file, FILE_MODE_READ, &vc));
419566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &C));
429566063dSJacob Faibussowitsch PetscCall(MatSetType(C, MATSEQBAIJ));
439566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(C));
449566063dSJacob Faibussowitsch PetscCall(MatLoad(C, vc));
459566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&vc));
46c4762a1bSJed Brown
479566063dSJacob Faibussowitsch PetscCall(MatGetSize(A, &m, &n));
489566063dSJacob Faibussowitsch PetscCall(MatGetSize(B, &m2, &n2));
4908401ef6SPierre Jolivet PetscCheck(m == m2, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Matrices are of different size. Cannot run this example");
50c4762a1bSJed Brown
51c4762a1bSJed Brown /* Test MatEqual() */
529566063dSJacob Faibussowitsch PetscCall(MatEqual(B, C, &tflg));
5328b400f6SJacob Faibussowitsch PetscCheck(tflg, PETSC_COMM_SELF, PETSC_ERR_PLIB, "MatEqual() failed");
54c4762a1bSJed Brown
55c4762a1bSJed Brown /* Test MatGetDiagonal() */
569566063dSJacob Faibussowitsch PetscCall(VecCreateSeq(PETSC_COMM_SELF, m, &x));
579566063dSJacob Faibussowitsch PetscCall(VecCreateSeq(PETSC_COMM_SELF, m, &y));
58c4762a1bSJed Brown
599566063dSJacob Faibussowitsch PetscCall(MatGetDiagonal(A, x));
609566063dSJacob Faibussowitsch PetscCall(MatGetDiagonal(B, y));
61c4762a1bSJed Brown
629566063dSJacob Faibussowitsch PetscCall(VecEqual(x, y, &tflg));
6328b400f6SJacob Faibussowitsch PetscCheck(tflg, PETSC_COMM_SELF, PETSC_ERR_PLIB, "MatGetDiagonal() failed");
64c4762a1bSJed Brown
65c4762a1bSJed Brown /* Test MatDiagonalScale() */
669566063dSJacob Faibussowitsch PetscCall(PetscRandomCreate(PETSC_COMM_SELF, &r));
679566063dSJacob Faibussowitsch PetscCall(PetscRandomSetFromOptions(r));
689566063dSJacob Faibussowitsch PetscCall(VecSetRandom(x, r));
699566063dSJacob Faibussowitsch PetscCall(VecSetRandom(y, r));
70c4762a1bSJed Brown
719566063dSJacob Faibussowitsch PetscCall(MatDiagonalScale(A, x, y));
729566063dSJacob Faibussowitsch PetscCall(MatDiagonalScale(B, x, y));
739566063dSJacob Faibussowitsch PetscCall(MatMult(A, x, y));
749566063dSJacob Faibussowitsch PetscCall(VecNorm(y, NORM_2, &norm1));
759566063dSJacob Faibussowitsch PetscCall(MatMult(B, x, y));
769566063dSJacob Faibussowitsch PetscCall(VecNorm(y, NORM_2, &norm2));
77c4762a1bSJed Brown rnorm = ((norm1 - norm2) * 100) / norm1;
78bc30f867SBarry Smith PetscCheck(rnorm >= -0.1 && rnorm <= 0.01, PETSC_COMM_SELF, PETSC_ERR_PLIB, "MatDiagonalScale() failed Norm1 %g Norm2 %g", (double)norm1, (double)norm2);
79c4762a1bSJed Brown
80c4762a1bSJed Brown /* Test MatGetRow()/ MatRestoreRow() */
81c4762a1bSJed Brown for (ct = 0; ct < 100; ct++) {
829566063dSJacob Faibussowitsch PetscCall(PetscRandomGetValue(r, &rval));
83c4762a1bSJed Brown row = (int)(rval * m);
849566063dSJacob Faibussowitsch PetscCall(MatGetRow(A, row, &ncols1, &cols1, &vals1));
859566063dSJacob Faibussowitsch PetscCall(MatGetRow(B, row, &ncols2, &cols2, &vals2));
86c4762a1bSJed Brown
87c4762a1bSJed Brown for (i = 0, j = 0; i < ncols1 && j < ncols2; i++) {
88c4762a1bSJed Brown while (cols2[j] != cols1[i]) j++;
8908401ef6SPierre Jolivet PetscCheck(vals1[i] == vals2[j], PETSC_COMM_SELF, PETSC_ERR_PLIB, "MatGetRow() failed - vals incorrect.");
90c4762a1bSJed Brown }
9108401ef6SPierre Jolivet PetscCheck(i >= ncols1, PETSC_COMM_SELF, PETSC_ERR_PLIB, "MatGetRow() failed - cols incorrect");
92c4762a1bSJed Brown
939566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(A, row, &ncols1, &cols1, &vals1));
949566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(B, row, &ncols2, &cols2, &vals2));
95c4762a1bSJed Brown }
96c4762a1bSJed Brown
979566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A));
989566063dSJacob Faibussowitsch PetscCall(MatDestroy(&B));
999566063dSJacob Faibussowitsch PetscCall(MatDestroy(&C));
1009566063dSJacob Faibussowitsch PetscCall(VecDestroy(&x));
1019566063dSJacob Faibussowitsch PetscCall(VecDestroy(&y));
1029566063dSJacob Faibussowitsch PetscCall(PetscRandomDestroy(&r));
1039566063dSJacob Faibussowitsch PetscCall(PetscFinalize());
104b122ec5aSJacob Faibussowitsch return 0;
105c4762a1bSJed Brown }
106c4762a1bSJed Brown
107c4762a1bSJed Brown /*TEST
108c4762a1bSJed Brown
109c4762a1bSJed Brown build:
110c4762a1bSJed Brown requires: !complex
111c4762a1bSJed Brown
112c4762a1bSJed Brown test:
113c4762a1bSJed Brown args: -f ${DATAFILESPATH}/matrices/cfd.1.10 -mat_block_size 5
114dfd57a17SPierre Jolivet requires: !complex double datafilespath !defined(PETSC_USE_64BIT_INDICES)
115*3886731fSPierre Jolivet output_file: output/empty.out
116c4762a1bSJed Brown
117c4762a1bSJed Brown TEST*/
118