1 static char help[] = "Testing MatCreateSeqBAIJWithArrays() and MatCreateSeqSBAIJWithArrays().\n\n"; 2 3 #include <petscmat.h> 4 5 int main(int argc,char **argv) 6 { 7 Mat A,B,As; 8 const PetscInt *ai,*aj; 9 PetscInt i,j,k,nz,n,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; 13 PetscErrorCode ierr; 14 PetscMPIInt size; 15 PetscBool flg; 16 17 ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr; 18 CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size)); 19 PetscCheckFalse(size != 1,PETSC_COMM_WORLD,PETSC_ERR_SUP,"This is a uniprocessor example only!"); 20 21 /* Create a aij matrix for checking */ 22 CHKERRQ(MatCreateSeqAIJ(PETSC_COMM_SELF,5,5,2,NULL,&A)); 23 CHKERRQ(PetscRandomCreate(PETSC_COMM_WORLD,&rctx)); 24 CHKERRQ(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 CHKERRQ(PetscRandomGetValue(rctx,&asa[k])); 31 CHKERRQ(MatSetValues(A,1,&i,1,&asj[k],&asa[k],INSERT_VALUES)); 32 CHKERRQ(MatSetValues(A,1,&i,1,&asj[k],&asa[k],INSERT_VALUES)); 33 if (i != asj[k]) { /* insert symmetric entry */ 34 CHKERRQ(MatSetValues(A,1,&asj[k],1,&i,&asa[k],INSERT_VALUES)); 35 } 36 k++; 37 } 38 } 39 CHKERRQ(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY)); 40 CHKERRQ(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY)); 41 42 /* Create a baij matrix using MatCreateSeqBAIJWithArrays() */ 43 CHKERRQ(MatGetRowIJ(A,0,PETSC_FALSE,PETSC_FALSE,&n,&ai,&aj,&flg)); 44 CHKERRQ(MatSeqAIJGetArray(A,&aa)); 45 /* WARNING: This sharing is dangerous if either A or B is later assembled */ 46 CHKERRQ(MatCreateSeqBAIJWithArrays(PETSC_COMM_SELF,1,5,5,(PetscInt*)ai,(PetscInt*)aj,aa,&B)); 47 CHKERRQ(MatSeqAIJRestoreArray(A,&aa)); 48 CHKERRQ(MatRestoreRowIJ(A,0,PETSC_FALSE,PETSC_FALSE,&n,&ai,&aj,&flg)); 49 CHKERRQ(MatMultEqual(A,B,10,&flg)); 50 PetscCheck(flg,PETSC_COMM_SELF,PETSC_ERR_ARG_NOTSAMETYPE,"MatMult(A,B) are NOT equal"); 51 52 /* Create a sbaij matrix using MatCreateSeqSBAIJWithArrays() */ 53 CHKERRQ(MatCreateSeqSBAIJWithArrays(PETSC_COMM_SELF,1,5,5,asi,asj,asa,&As)); 54 CHKERRQ(MatMultEqual(A,As,10,&flg)); 55 PetscCheck(flg,PETSC_COMM_SELF,PETSC_ERR_ARG_NOTSAMETYPE,"MatMult(A,As) are NOT equal"); 56 57 /* Free spaces */ 58 CHKERRQ(PetscRandomDestroy(&rctx)); 59 CHKERRQ(MatDestroy(&A)); 60 CHKERRQ(MatDestroy(&B)); 61 CHKERRQ(MatDestroy(&As)); 62 ierr = PetscFinalize(); 63 return ierr; 64 } 65