xref: /petsc/src/mat/tests/ex100.c (revision c4762a1b19cd2af06abeed90e8f9d34fb975dd94)
1*c4762a1bSJed Brown 
2*c4762a1bSJed Brown static char help[] = "Tests various routines in MatMAIJ format.\n";
3*c4762a1bSJed Brown 
4*c4762a1bSJed Brown #include <petscmat.h>
5*c4762a1bSJed Brown #define IMAX 15
6*c4762a1bSJed Brown int main(int argc,char **args)
7*c4762a1bSJed Brown {
8*c4762a1bSJed Brown   Mat            A,B,MA;
9*c4762a1bSJed Brown   PetscViewer    fd;
10*c4762a1bSJed Brown   char           file[PETSC_MAX_PATH_LEN];
11*c4762a1bSJed Brown   PetscInt       m,n,M,N,dof=1;
12*c4762a1bSJed Brown   PetscMPIInt    rank,size;
13*c4762a1bSJed Brown   PetscErrorCode ierr;
14*c4762a1bSJed Brown   PetscBool      flg;
15*c4762a1bSJed Brown 
16*c4762a1bSJed Brown   ierr = PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr;
17*c4762a1bSJed Brown   ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank);CHKERRQ(ierr);
18*c4762a1bSJed Brown   ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRQ(ierr);
19*c4762a1bSJed Brown 
20*c4762a1bSJed Brown   /* Load aij matrix A */
21*c4762a1bSJed Brown   ierr = PetscOptionsGetString(NULL,NULL,"-f",file,PETSC_MAX_PATH_LEN,&flg);CHKERRQ(ierr);
22*c4762a1bSJed Brown   if (!flg) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_USER,"Must indicate binary file with the -f option");
23*c4762a1bSJed Brown   ierr = PetscViewerBinaryOpen(PETSC_COMM_WORLD,file,FILE_MODE_READ,&fd);CHKERRQ(ierr);
24*c4762a1bSJed Brown   ierr = MatCreate(PETSC_COMM_WORLD,&A);CHKERRQ(ierr);
25*c4762a1bSJed Brown   ierr = MatLoad(A,fd);CHKERRQ(ierr);
26*c4762a1bSJed Brown   ierr = PetscViewerDestroy(&fd);CHKERRQ(ierr);
27*c4762a1bSJed Brown 
28*c4762a1bSJed Brown   /* Get dof, then create maij matrix MA */
29*c4762a1bSJed Brown   ierr = PetscOptionsGetInt(NULL,NULL,"-dof",&dof,NULL);CHKERRQ(ierr);
30*c4762a1bSJed Brown   ierr = MatCreateMAIJ(A,dof,&MA);CHKERRQ(ierr);
31*c4762a1bSJed Brown   ierr = MatGetLocalSize(MA,&m,&n);CHKERRQ(ierr);
32*c4762a1bSJed Brown   ierr = MatGetSize(MA,&M,&N);CHKERRQ(ierr);
33*c4762a1bSJed Brown 
34*c4762a1bSJed Brown   if (size == 1) {
35*c4762a1bSJed Brown     ierr = MatConvert(MA,MATSEQAIJ,MAT_INITIAL_MATRIX,&B);CHKERRQ(ierr);
36*c4762a1bSJed Brown   } else {
37*c4762a1bSJed Brown     ierr = MatConvert(MA,MATMPIAIJ,MAT_INITIAL_MATRIX,&B);CHKERRQ(ierr);
38*c4762a1bSJed Brown   }
39*c4762a1bSJed Brown 
40*c4762a1bSJed Brown   /* Test MatMult() */
41*c4762a1bSJed Brown   ierr = MatMultEqual(MA,B,10,&flg);CHKERRQ(ierr);
42*c4762a1bSJed Brown   if (!flg) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_CONV_FAILED,"Error: MatMul() for MAIJ matrix");
43*c4762a1bSJed Brown   /* Test MatMultAdd() */
44*c4762a1bSJed Brown   ierr = MatMultAddEqual(MA,B,10,&flg);CHKERRQ(ierr);
45*c4762a1bSJed Brown   if (!flg) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_CONV_FAILED,"Error: MatMulAdd() for MAIJ matrix");
46*c4762a1bSJed Brown 
47*c4762a1bSJed Brown   /* Test MatMultTranspose() */
48*c4762a1bSJed Brown   ierr = MatMultTransposeEqual(MA,B,10,&flg);CHKERRQ(ierr);
49*c4762a1bSJed Brown   if (!flg) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_CONV_FAILED,"Error: MatMulAdd() for MAIJ matrix");
50*c4762a1bSJed Brown 
51*c4762a1bSJed Brown   /* Test MatMultTransposeAdd() */
52*c4762a1bSJed Brown   ierr = MatMultTransposeAddEqual(MA,B,10,&flg);CHKERRQ(ierr);
53*c4762a1bSJed Brown   if (!flg) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_CONV_FAILED,"Error: MatMulTransposeAdd() for MAIJ matrix");
54*c4762a1bSJed Brown 
55*c4762a1bSJed Brown   ierr = MatDestroy(&MA);CHKERRQ(ierr);
56*c4762a1bSJed Brown   ierr = MatDestroy(&A);CHKERRQ(ierr);
57*c4762a1bSJed Brown   ierr = MatDestroy(&B);CHKERRQ(ierr);
58*c4762a1bSJed Brown   ierr = PetscFinalize();
59*c4762a1bSJed Brown   return ierr;
60*c4762a1bSJed Brown }
61*c4762a1bSJed Brown 
62*c4762a1bSJed Brown 
63*c4762a1bSJed Brown /*TEST
64*c4762a1bSJed Brown 
65*c4762a1bSJed Brown    build:
66*c4762a1bSJed Brown       requires: !complex
67*c4762a1bSJed Brown 
68*c4762a1bSJed Brown    test:
69*c4762a1bSJed Brown       nsize: {{1 3}}
70*c4762a1bSJed Brown       requires: datafilespath !complex double !define(PETSC_USE_64BIT_INDICES)
71*c4762a1bSJed Brown       args: -f ${DATAFILESPATH}/matrices/arco1 -dof {{1 2 3 4 5 6 8 9 16}} -viewer_binary_skip_info
72*c4762a1bSJed Brown 
73*c4762a1bSJed Brown TEST*/
74