Actual source code: ex122.c

petsc-3.8.4 2018-03-24
Report Typos and Errors
  1: static char help[] = "Test MatMatMult() for AIJ and Dense matrices.\n\n";

  3:  #include <petscmat.h>

  5: int main(int argc,char **argv)
  6: {
  7:   Mat            A,B,C,D;
  8:   PetscInt       i,M=10,N=5;
 10:   PetscRandom    r;
 11:   PetscBool      equal=PETSC_FALSE;
 12:   PetscReal      fill = 1.0;
 13:   PetscInt       nza,am,an;

 15:   PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
 16:   PetscOptionsGetInt(NULL,NULL,"-M",&M,NULL);
 17:   PetscOptionsGetInt(NULL,NULL,"-N",&N,NULL);
 18:   PetscOptionsGetReal(NULL,NULL,"-fill",&fill,NULL);

 20:   PetscRandomCreate(PETSC_COMM_WORLD,&r);
 21:   PetscRandomSetFromOptions(r);

 23:   /* create a aij matrix A */
 24:   MatCreate(PETSC_COMM_WORLD,&A);
 25:   MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,N,M);
 26:   MatSetType(A,MATAIJ);
 27:   nza  = (PetscInt)(.3*M); /* num of nozeros in each row of A */
 28:   MatSeqAIJSetPreallocation(A,nza,NULL);
 29:   MatMPIAIJSetPreallocation(A,nza,NULL,nza,NULL);
 30:   MatSetRandom(A,r);

 32:   /* create a dense matrix B */
 33:   MatGetLocalSize(A,&am,&an);
 34:   MatCreate(PETSC_COMM_WORLD,&B);
 35:   MatSetSizes(B,PETSC_DECIDE,am,N,PETSC_DECIDE);
 36:   MatSetType(B,MATDENSE);
 37:   MatSeqDenseSetPreallocation(B,NULL);
 38:   MatMPIDenseSetPreallocation(B,NULL);
 39:   MatSetRandom(B,r);
 40:   PetscRandomDestroy(&r);
 41:   MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
 42:   MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);


 45:   /* Test MatMatMult() */
 46:   MatMatMult(B,A,MAT_INITIAL_MATRIX,fill,&C);

 48:   MatMatMultSymbolic(B,A,fill,&D);
 49:   for (i=0; i<2; i++) {
 50:     MatMatMultNumeric(B,A,D);
 51:   }
 52:   MatEqual(C,D,&equal);
 53:   if (!equal) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"C != D");

 55:   MatDestroy(&D);
 56:   MatDestroy(&C);
 57:   MatDestroy(&B);
 58:   MatDestroy(&A);
 59:   PetscFinalize();
 60:   return ierr;
 61: }