Actual source code: ex108.c

  1: static char help[] = "Testing MatCreateSeqBAIJWithArrays() and MatCreateSeqSBAIJWithArrays().\n\n";

 3:  #include petscmat.h

  7: int main(int argc,char **argv) {
  8:   Mat            A,B,As;
  9:   PetscInt       i,j,k,nz,n,*ai,*aj,asi[]={0,2,3,4,6,7};
 10:   PetscInt       asj[]={0,4,1,2,3,4,4};
 11:   PetscScalar    asa[7],*aa;
 12:   PetscRandom    rctx;
 14:   PetscMPIInt    size;
 15:   PetscTruth     flg;

 17:   PetscInitialize(&argc,&argv,(char *)0,help);
 18:   MPI_Comm_size(PETSC_COMM_WORLD,&size);
 19:   if (size != 1) SETERRQ(PETSC_ERR_SUP,"This is a uniprocessor example only!");

 21:   /* Create a aij matrix for checking */
 22:   MatCreateSeqAIJ(PETSC_COMM_SELF,5,5,2,PETSC_NULL,&A);
 23:   PetscRandomCreate(PETSC_COMM_WORLD,&rctx);
 24:   PetscRandomSetFromOptions(rctx);
 25: 
 26:   k = 0;
 27:   for (i=0; i<5; i++) {
 28:     nz = asi[i+1] - asi[i];  /* length of i_th row of A */
 29:     for (j=0; j<nz; j++){
 30:       PetscRandomGetValue(rctx,&asa[k]);
 31:       MatSetValues(A,1,&i,1,&asj[k],&asa[k],INSERT_VALUES);
 32:       MatSetValues(A,1,&i,1,&asj[k],&asa[k],INSERT_VALUES);
 33:       if (i != asj[k]){ /* insert symmetric entry */
 34:         MatSetValues(A,1,&asj[k],1,&i,&asa[k],INSERT_VALUES);
 35:       }
 36:       k++;
 37:     }
 38:   }
 39:   MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
 40:   MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);

 42:   /* Create a baij matrix using MatCreateSeqBAIJWithArrays() */
 43:   MatGetRowIJ(A,0,PETSC_FALSE,&n,&ai,&aj,&flg);
 44:   MatGetArray(A,&aa);
 45:   MatCreateSeqBAIJWithArrays(PETSC_COMM_SELF,1,5,5,ai,aj,aa,&B);
 46:   MatRestoreArray(A,&aa);
 47:   MatMultEqual(A,B,10,&flg);
 48:   if (!flg) SETERRQ(PETSC_ERR_ARG_NOTSAMETYPE,"MatMult(A,B) are NOT equal");

 50:   /* Create a sbaij matrix using MatCreateSeqSBAIJWithArrays() */
 51:   MatCreateSeqSBAIJWithArrays(PETSC_COMM_SELF,1,5,5,asi,asj,asa,&As);
 52:   MatMultEqual(A,As,10,&flg);
 53:   if (!flg) SETERRQ(PETSC_ERR_ARG_NOTSAMETYPE,"MatMult(A,As) are NOT equal");
 54: 
 55:   /* Free spaces */
 56:   PetscRandomDestroy(rctx);
 57:   MatDestroy(A);
 58:   MatDestroy(B);
 59:   MatDestroy(As);
 60:   PetscFinalize();
 61:   return(0);
 62: }