Actual source code: ex109.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,N,Istart,Iend,n=7,j,J,Ii,m=8,am,an;
  9:   PetscScalar    v;
 11:   PetscRandom    r;
 12:   PetscBool      equal=PETSC_FALSE;
 13:   PetscReal      fill = 1.0;
 14:   PetscMPIInt    size;

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

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

 24:   /* Create a aij matrix A */
 25:   M    = N = m*n;
 26:   MatCreate(PETSC_COMM_WORLD,&A);
 27:   MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,M,N);
 28:   MatSetType(A,MATAIJ);
 29:   MatSetFromOptions(A);
 30:   MatMPIAIJSetPreallocation(A,5,NULL,5,NULL);
 31:   MatSeqAIJSetPreallocation(A,5,NULL);

 33:   MatGetOwnershipRange(A,&Istart,&Iend);
 34:   am   = Iend - Istart;
 35:   for (Ii=Istart; Ii<Iend; Ii++) {
 36:     v = -1.0; i = Ii/n; j = Ii - i*n;
 37:     if (i>0)   {J = Ii - n; MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES);}
 38:     if (i<m-1) {J = Ii + n; MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES);}
 39:     if (j>0)   {J = Ii - 1; MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES);}
 40:     if (j<n-1) {J = Ii + 1; MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES);}
 41:     v = 4.0; MatSetValues(A,1,&Ii,1,&Ii,&v,INSERT_VALUES);
 42:   }
 43:   MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
 44:   MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);

 46:   /* Create a dense matrix B */
 47:   MatGetLocalSize(A,&am,&an);
 48:   MatCreate(PETSC_COMM_WORLD,&B);
 49:   MatSetSizes(B,an,PETSC_DECIDE,PETSC_DECIDE,M);
 50:   MatSetType(B,MATDENSE);
 51:   MatSeqDenseSetPreallocation(B,NULL);
 52:   MatMPIDenseSetPreallocation(B,NULL);
 53:   MatSetFromOptions(B);
 54:   MatSetRandom(B,r);
 55:   PetscRandomDestroy(&r);
 56:   MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
 57:   MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);

 59:   /* Test C = A*B (aij*dense) */
 60:   MatMatMult(A,B,MAT_INITIAL_MATRIX,fill,&C);
 61:   MatMatMult(A,B,MAT_REUSE_MATRIX,fill,&C);

 63:   MatMatMultSymbolic(A,B,fill,&D);
 64:   for (i=0; i<2; i++) {
 65:     MatMatMultNumeric(A,B,D);
 66:   }
 67:   MatEqual(C,D,&equal);
 68:   if (!equal) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"C != D");
 69:   MatDestroy(&D);

 71:   /* Test D = C*A (dense*aij) */
 72:   MatMatMult(C,A,MAT_INITIAL_MATRIX,fill,&D);
 73:   MatMatMult(C,A,MAT_REUSE_MATRIX,fill,&D);
 74:   MatDestroy(&D);

 76:   /* Test D = A*C (aij*dense) */
 77:   MatMatMult(A,C,MAT_INITIAL_MATRIX,fill,&D);
 78:   MatMatMult(A,C,MAT_REUSE_MATRIX,fill,&D);
 79:   MatDestroy(&D);

 81:   /* Test D = B*C (dense*dense) */
 82:   MPI_Comm_size(PETSC_COMM_WORLD,&size);
 83:   if (size == 1) {
 84:     MatMatMult(B,C,MAT_INITIAL_MATRIX,fill,&D);
 85:     MatMatMult(B,C,MAT_REUSE_MATRIX,fill,&D);
 86:     MatDestroy(&D);
 87:   }

 89:   MatDestroy(&C);
 90:   MatDestroy(&B);
 91:   MatDestroy(&A);
 92:   PetscFinalize();
 93:   return ierr;
 94: }