xref: /petsc/src/mat/tests/ex131.c (revision 9566063d113dddea24716c546802770db7481bc0)
1c4762a1bSJed Brown 
2c4762a1bSJed Brown static char help[] = "Tests MatMult() on MatLoad() matrix \n\n";
3c4762a1bSJed Brown 
4c4762a1bSJed Brown #include <petscmat.h>
5c4762a1bSJed Brown 
6c4762a1bSJed Brown int main(int argc,char **args)
7c4762a1bSJed Brown {
8c4762a1bSJed Brown   Mat            A;
9c4762a1bSJed Brown   Vec            x,b;
10c4762a1bSJed Brown   PetscViewer    fd;              /* viewer */
11c4762a1bSJed Brown   char           file[PETSC_MAX_PATH_LEN]; /* input file name */
12c4762a1bSJed Brown   PetscBool      flg;
13c4762a1bSJed Brown 
14*9566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc,&args,(char*)0,help));
15c4762a1bSJed Brown   /* Determine file from which we read the matrix A */
16*9566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetString(NULL,NULL,"-f",file,sizeof(file),&flg));
1728b400f6SJacob Faibussowitsch   PetscCheck(flg,PETSC_COMM_WORLD,PETSC_ERR_USER,"Must indicate binary file with the -f option");
18c4762a1bSJed Brown 
19c4762a1bSJed Brown   /* Load matrix A */
20*9566063dSJacob Faibussowitsch   PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD,file,FILE_MODE_READ,&fd));
21*9566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD,&A));
22*9566063dSJacob Faibussowitsch   PetscCall(MatLoad(A,fd));
23c4762a1bSJed Brown   flg  = PETSC_FALSE;
24*9566063dSJacob Faibussowitsch   PetscCall(VecCreate(PETSC_COMM_WORLD,&x));
25*9566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetString(NULL,NULL,"-vec",file,sizeof(file),&flg));
26c4762a1bSJed Brown   if (flg) {
27c4762a1bSJed Brown     if (file[0] == '0') {
28c4762a1bSJed Brown       PetscInt    m;
29c4762a1bSJed Brown       PetscScalar one = 1.0;
30*9566063dSJacob Faibussowitsch       PetscCall(PetscInfo(0,"Using vector of ones for RHS\n"));
31*9566063dSJacob Faibussowitsch       PetscCall(MatGetLocalSize(A,&m,NULL));
32*9566063dSJacob Faibussowitsch       PetscCall(VecSetSizes(x,m,PETSC_DECIDE));
33*9566063dSJacob Faibussowitsch       PetscCall(VecSetFromOptions(x));
34*9566063dSJacob Faibussowitsch       PetscCall(VecSet(x,one));
35c4762a1bSJed Brown     }
36c4762a1bSJed Brown   } else {
37*9566063dSJacob Faibussowitsch     PetscCall(VecLoad(x,fd));
38*9566063dSJacob Faibussowitsch     PetscCall(PetscViewerDestroy(&fd));
39c4762a1bSJed Brown   }
40*9566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(x,&b));
41*9566063dSJacob Faibussowitsch   PetscCall(MatMult(A,x,b));
42c4762a1bSJed Brown 
43c4762a1bSJed Brown   /* Print (for testing only) */
44*9566063dSJacob Faibussowitsch   PetscCall(MatView(A,0));
45*9566063dSJacob Faibussowitsch   PetscCall(VecView(b,0));
46c4762a1bSJed Brown   /* Free data structures */
47*9566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&A));
48*9566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&x));
49*9566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&b));
50*9566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
51b122ec5aSJacob Faibussowitsch   return 0;
52c4762a1bSJed Brown }
53