xref: /petsc/src/mat/tests/ex100.c (revision d71ae5a4db6382e7f06317b8d368875286fe9008)
1c4762a1bSJed Brown 
2c4762a1bSJed Brown static char help[] = "Tests various routines in MatMAIJ format.\n";
3c4762a1bSJed Brown 
4c4762a1bSJed Brown #include <petscmat.h>
5c4762a1bSJed Brown #define IMAX 15
6*d71ae5a4SJacob Faibussowitsch int main(int argc, char **args)
7*d71ae5a4SJacob Faibussowitsch {
8c4762a1bSJed Brown   Mat         A, B, MA;
9c4762a1bSJed Brown   PetscViewer fd;
10c4762a1bSJed Brown   char        file[PETSC_MAX_PATH_LEN];
11c4762a1bSJed Brown   PetscInt    m, n, M, N, dof = 1;
12c4762a1bSJed Brown   PetscMPIInt rank, size;
13c4762a1bSJed Brown   PetscBool   flg;
14c4762a1bSJed Brown 
15327415f7SBarry Smith   PetscFunctionBeginUser;
169566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &args, (char *)0, help));
179566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
189566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
19c4762a1bSJed Brown 
20c4762a1bSJed Brown   /* Load aij matrix A */
219566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetString(NULL, NULL, "-f", file, sizeof(file), &flg));
2228b400f6SJacob Faibussowitsch   PetscCheck(flg, PETSC_COMM_WORLD, PETSC_ERR_USER, "Must indicate binary file with the -f option");
239566063dSJacob Faibussowitsch   PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD, file, FILE_MODE_READ, &fd));
249566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD, &A));
259566063dSJacob Faibussowitsch   PetscCall(MatLoad(A, fd));
269566063dSJacob Faibussowitsch   PetscCall(PetscViewerDestroy(&fd));
27c4762a1bSJed Brown 
28c4762a1bSJed Brown   /* Get dof, then create maij matrix MA */
299566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-dof", &dof, NULL));
309566063dSJacob Faibussowitsch   PetscCall(MatCreateMAIJ(A, dof, &MA));
319566063dSJacob Faibussowitsch   PetscCall(MatGetLocalSize(MA, &m, &n));
329566063dSJacob Faibussowitsch   PetscCall(MatGetSize(MA, &M, &N));
33c4762a1bSJed Brown 
34c4762a1bSJed Brown   if (size == 1) {
359566063dSJacob Faibussowitsch     PetscCall(MatConvert(MA, MATSEQAIJ, MAT_INITIAL_MATRIX, &B));
36c4762a1bSJed Brown   } else {
379566063dSJacob Faibussowitsch     PetscCall(MatConvert(MA, MATMPIAIJ, MAT_INITIAL_MATRIX, &B));
38c4762a1bSJed Brown   }
39c4762a1bSJed Brown 
40c4762a1bSJed Brown   /* Test MatMult() */
419566063dSJacob Faibussowitsch   PetscCall(MatMultEqual(MA, B, 10, &flg));
4228b400f6SJacob Faibussowitsch   PetscCheck(flg, PETSC_COMM_WORLD, PETSC_ERR_CONV_FAILED, "Error: MatMul() for MAIJ matrix");
43c4762a1bSJed Brown   /* Test MatMultAdd() */
449566063dSJacob Faibussowitsch   PetscCall(MatMultAddEqual(MA, B, 10, &flg));
4528b400f6SJacob Faibussowitsch   PetscCheck(flg, PETSC_COMM_WORLD, PETSC_ERR_CONV_FAILED, "Error: MatMulAdd() for MAIJ matrix");
46c4762a1bSJed Brown 
47c4762a1bSJed Brown   /* Test MatMultTranspose() */
489566063dSJacob Faibussowitsch   PetscCall(MatMultTransposeEqual(MA, B, 10, &flg));
4928b400f6SJacob Faibussowitsch   PetscCheck(flg, PETSC_COMM_WORLD, PETSC_ERR_CONV_FAILED, "Error: MatMulAdd() for MAIJ matrix");
50c4762a1bSJed Brown 
51c4762a1bSJed Brown   /* Test MatMultTransposeAdd() */
529566063dSJacob Faibussowitsch   PetscCall(MatMultTransposeAddEqual(MA, B, 10, &flg));
5328b400f6SJacob Faibussowitsch   PetscCheck(flg, PETSC_COMM_WORLD, PETSC_ERR_CONV_FAILED, "Error: MatMulTransposeAdd() for MAIJ matrix");
54c4762a1bSJed Brown 
559566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&MA));
569566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&A));
579566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&B));
589566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
59b122ec5aSJacob Faibussowitsch   return 0;
60c4762a1bSJed Brown }
61c4762a1bSJed Brown 
62c4762a1bSJed Brown /*TEST
63c4762a1bSJed Brown 
64c4762a1bSJed Brown    build:
65c4762a1bSJed Brown       requires: !complex
66c4762a1bSJed Brown 
67c4762a1bSJed Brown    test:
68c4762a1bSJed Brown       nsize: {{1 3}}
69dfd57a17SPierre Jolivet       requires: datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES)
70c4762a1bSJed Brown       args: -f ${DATAFILESPATH}/matrices/arco1 -dof {{1 2 3 4 5 6 8 9 16}} -viewer_binary_skip_info
71c4762a1bSJed Brown 
72c4762a1bSJed Brown TEST*/
73