Actual source code: ex122.c

petsc-3.6.1 2015-08-06
Report Typos and Errors
  1: static char help[] = "Test MatMatMult() for AIJ and Dense matrices.\n\n";

  3: #include <petscmat.h>

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

 17:   PetscInitialize(&argc,&argv,(char*)0,help);
 18:   PetscOptionsGetInt(NULL,"-M",&M,NULL);
 19:   PetscOptionsGetInt(NULL,"-N",&N,NULL);
 20:   PetscOptionsGetReal(NULL,"-fill",&fill,NULL);

 22:   PetscRandomCreate(PETSC_COMM_WORLD,&r);
 23:   PetscRandomSetFromOptions(r);

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

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


 47:   /* Test MatMatMult() */
 48:   MatMatMult(B,A,MAT_INITIAL_MATRIX,fill,&C);

 50:   /*
 51:   PetscViewerSetFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_MATLAB);
 52:   MatView(C,0);
 53:   MatView(B,0);
 54:   MatView(A,0);
 55:   */


 58:   MatMatMultSymbolic(B,A,fill,&D);
 59:   for (i=0; i<2; i++) {
 60:     MatMatMultNumeric(B,A,D);
 61:   }
 62:   MatEqual(C,D,&equal);
 63:   if (!equal) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"C != D");

 65:   MatDestroy(&D);
 66:   MatDestroy(&C);
 67:   MatDestroy(&B);
 68:   MatDestroy(&A);
 69:   PetscFinalize();
 70:   return(0);
 71: }