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 ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRQ(ierr); 19 if (size != 1) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"This is a uniprocessor example only!"); 20 21 /* Create a aij matrix for checking */ 22 ierr = MatCreateSeqAIJ(PETSC_COMM_SELF,5,5,2,NULL,&A);CHKERRQ(ierr); 23 ierr = PetscRandomCreate(PETSC_COMM_WORLD,&rctx);CHKERRQ(ierr); 24 ierr = PetscRandomSetFromOptions(rctx);CHKERRQ(ierr); 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 ierr = PetscRandomGetValue(rctx,&asa[k]);CHKERRQ(ierr); 31 ierr = MatSetValues(A,1,&i,1,&asj[k],&asa[k],INSERT_VALUES);CHKERRQ(ierr); 32 ierr = MatSetValues(A,1,&i,1,&asj[k],&asa[k],INSERT_VALUES);CHKERRQ(ierr); 33 if (i != asj[k]) { /* insert symmetric entry */ 34 ierr = MatSetValues(A,1,&asj[k],1,&i,&asa[k],INSERT_VALUES);CHKERRQ(ierr); 35 } 36 k++; 37 } 38 } 39 ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 40 ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 41 42 /* Create a baij matrix using MatCreateSeqBAIJWithArrays() */ 43 ierr = MatGetRowIJ(A,0,PETSC_FALSE,PETSC_FALSE,&n,&ai,&aj,&flg);CHKERRQ(ierr); 44 ierr = MatSeqAIJGetArray(A,&aa);CHKERRQ(ierr); 45 /* WARNING: This sharing is dangerous if either A or B is later assembled */ 46 ierr = MatCreateSeqBAIJWithArrays(PETSC_COMM_SELF,1,5,5,(PetscInt*)ai,(PetscInt*)aj,aa,&B);CHKERRQ(ierr); 47 ierr = MatSeqAIJRestoreArray(A,&aa);CHKERRQ(ierr); 48 ierr = MatRestoreRowIJ(A,0,PETSC_FALSE,PETSC_FALSE,&n,&ai,&aj,&flg);CHKERRQ(ierr); 49 ierr = MatMultEqual(A,B,10,&flg);CHKERRQ(ierr); 50 if (!flg) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NOTSAMETYPE,"MatMult(A,B) are NOT equal"); 51 52 /* Create a sbaij matrix using MatCreateSeqSBAIJWithArrays() */ 53 ierr = MatCreateSeqSBAIJWithArrays(PETSC_COMM_SELF,1,5,5,asi,asj,asa,&As);CHKERRQ(ierr); 54 ierr = MatMultEqual(A,As,10,&flg);CHKERRQ(ierr); 55 if (!flg) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NOTSAMETYPE,"MatMult(A,As) are NOT equal"); 56 57 /* Free spaces */ 58 ierr = PetscRandomDestroy(&rctx);CHKERRQ(ierr); 59 ierr = MatDestroy(&A);CHKERRQ(ierr); 60 ierr = MatDestroy(&B);CHKERRQ(ierr); 61 ierr = MatDestroy(&As);CHKERRQ(ierr); 62 ierr = PetscFinalize(); 63 return ierr; 64 } 65