xref: /petsc/src/mat/tests/ex108.c (revision b122ec5aa1bd4469eb4e0673542fb7de3f411254)
1c4762a1bSJed Brown static char help[] = "Testing MatCreateSeqBAIJWithArrays() and MatCreateSeqSBAIJWithArrays().\n\n";
2c4762a1bSJed Brown 
3c4762a1bSJed Brown #include <petscmat.h>
4c4762a1bSJed Brown 
5c4762a1bSJed Brown int main(int argc,char **argv)
6c4762a1bSJed Brown {
7c4762a1bSJed Brown   Mat            A,B,As;
8c4762a1bSJed Brown   const PetscInt *ai,*aj;
9c4762a1bSJed Brown   PetscInt       i,j,k,nz,n,asi[]={0,2,3,4,6,7};
10c4762a1bSJed Brown   PetscInt       asj[]={0,4,1,2,3,4,4};
11c4762a1bSJed Brown   PetscScalar    asa[7],*aa;
12c4762a1bSJed Brown   PetscRandom    rctx;
13c4762a1bSJed Brown   PetscMPIInt    size;
14c4762a1bSJed Brown   PetscBool      flg;
15c4762a1bSJed Brown 
16*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscInitialize(&argc,&argv,(char*)0,help));
175f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size));
182c71b3e2SJacob Faibussowitsch   PetscCheckFalse(size != 1,PETSC_COMM_WORLD,PETSC_ERR_SUP,"This is a uniprocessor example only!");
19c4762a1bSJed Brown 
20c4762a1bSJed Brown   /* Create a aij matrix for checking */
215f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreateSeqAIJ(PETSC_COMM_SELF,5,5,2,NULL,&A));
225f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscRandomCreate(PETSC_COMM_WORLD,&rctx));
235f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscRandomSetFromOptions(rctx));
24c4762a1bSJed Brown 
25c4762a1bSJed Brown   k = 0;
26c4762a1bSJed Brown   for (i=0; i<5; i++) {
27c4762a1bSJed Brown     nz = asi[i+1] - asi[i];  /* length of i_th row of A */
28c4762a1bSJed Brown     for (j=0; j<nz; j++) {
295f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscRandomGetValue(rctx,&asa[k]));
305f80ce2aSJacob Faibussowitsch       CHKERRQ(MatSetValues(A,1,&i,1,&asj[k],&asa[k],INSERT_VALUES));
315f80ce2aSJacob Faibussowitsch       CHKERRQ(MatSetValues(A,1,&i,1,&asj[k],&asa[k],INSERT_VALUES));
32c4762a1bSJed Brown       if (i != asj[k]) { /* insert symmetric entry */
335f80ce2aSJacob Faibussowitsch         CHKERRQ(MatSetValues(A,1,&asj[k],1,&i,&asa[k],INSERT_VALUES));
34c4762a1bSJed Brown       }
35c4762a1bSJed Brown       k++;
36c4762a1bSJed Brown     }
37c4762a1bSJed Brown   }
385f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
395f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
40c4762a1bSJed Brown 
41c4762a1bSJed Brown   /* Create a baij matrix using MatCreateSeqBAIJWithArrays() */
425f80ce2aSJacob Faibussowitsch   CHKERRQ(MatGetRowIJ(A,0,PETSC_FALSE,PETSC_FALSE,&n,&ai,&aj,&flg));
435f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSeqAIJGetArray(A,&aa));
44c4762a1bSJed Brown   /* WARNING: This sharing is dangerous if either A or B is later assembled */
455f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreateSeqBAIJWithArrays(PETSC_COMM_SELF,1,5,5,(PetscInt*)ai,(PetscInt*)aj,aa,&B));
465f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSeqAIJRestoreArray(A,&aa));
475f80ce2aSJacob Faibussowitsch   CHKERRQ(MatRestoreRowIJ(A,0,PETSC_FALSE,PETSC_FALSE,&n,&ai,&aj,&flg));
485f80ce2aSJacob Faibussowitsch   CHKERRQ(MatMultEqual(A,B,10,&flg));
4928b400f6SJacob Faibussowitsch   PetscCheck(flg,PETSC_COMM_SELF,PETSC_ERR_ARG_NOTSAMETYPE,"MatMult(A,B) are NOT equal");
50c4762a1bSJed Brown 
51c4762a1bSJed Brown   /* Create a sbaij matrix using MatCreateSeqSBAIJWithArrays() */
525f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreateSeqSBAIJWithArrays(PETSC_COMM_SELF,1,5,5,asi,asj,asa,&As));
535f80ce2aSJacob Faibussowitsch   CHKERRQ(MatMultEqual(A,As,10,&flg));
5428b400f6SJacob Faibussowitsch   PetscCheck(flg,PETSC_COMM_SELF,PETSC_ERR_ARG_NOTSAMETYPE,"MatMult(A,As) are NOT equal");
55c4762a1bSJed Brown 
56c4762a1bSJed Brown   /* Free spaces */
575f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscRandomDestroy(&rctx));
585f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&A));
595f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&B));
605f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&As));
61*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscFinalize());
62*b122ec5aSJacob Faibussowitsch   return 0;
63c4762a1bSJed Brown }
64