xref: /petsc/src/mat/tests/ex108.c (revision 2c71b3e237ead271e4f3aa1505f92bf476e3413d)
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   PetscErrorCode ierr;
14c4762a1bSJed Brown   PetscMPIInt    size;
15c4762a1bSJed Brown   PetscBool      flg;
16c4762a1bSJed Brown 
17c4762a1bSJed Brown   ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
18ffc4695bSBarry Smith   ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRMPI(ierr);
19*2c71b3e2SJacob Faibussowitsch   PetscCheckFalse(size != 1,PETSC_COMM_WORLD,PETSC_ERR_SUP,"This is a uniprocessor example only!");
20c4762a1bSJed Brown 
21c4762a1bSJed Brown   /* Create a aij matrix for checking */
22c4762a1bSJed Brown   ierr = MatCreateSeqAIJ(PETSC_COMM_SELF,5,5,2,NULL,&A);CHKERRQ(ierr);
23c4762a1bSJed Brown   ierr = PetscRandomCreate(PETSC_COMM_WORLD,&rctx);CHKERRQ(ierr);
24c4762a1bSJed Brown   ierr = PetscRandomSetFromOptions(rctx);CHKERRQ(ierr);
25c4762a1bSJed Brown 
26c4762a1bSJed Brown   k = 0;
27c4762a1bSJed Brown   for (i=0; i<5; i++) {
28c4762a1bSJed Brown     nz = asi[i+1] - asi[i];  /* length of i_th row of A */
29c4762a1bSJed Brown     for (j=0; j<nz; j++) {
30c4762a1bSJed Brown       ierr = PetscRandomGetValue(rctx,&asa[k]);CHKERRQ(ierr);
31c4762a1bSJed Brown       ierr = MatSetValues(A,1,&i,1,&asj[k],&asa[k],INSERT_VALUES);CHKERRQ(ierr);
32c4762a1bSJed Brown       ierr = MatSetValues(A,1,&i,1,&asj[k],&asa[k],INSERT_VALUES);CHKERRQ(ierr);
33c4762a1bSJed Brown       if (i != asj[k]) { /* insert symmetric entry */
34c4762a1bSJed Brown         ierr = MatSetValues(A,1,&asj[k],1,&i,&asa[k],INSERT_VALUES);CHKERRQ(ierr);
35c4762a1bSJed Brown       }
36c4762a1bSJed Brown       k++;
37c4762a1bSJed Brown     }
38c4762a1bSJed Brown   }
39c4762a1bSJed Brown   ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
40c4762a1bSJed Brown   ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
41c4762a1bSJed Brown 
42c4762a1bSJed Brown   /* Create a baij matrix using MatCreateSeqBAIJWithArrays() */
43c4762a1bSJed Brown   ierr = MatGetRowIJ(A,0,PETSC_FALSE,PETSC_FALSE,&n,&ai,&aj,&flg);CHKERRQ(ierr);
44c4762a1bSJed Brown   ierr = MatSeqAIJGetArray(A,&aa);CHKERRQ(ierr);
45c4762a1bSJed Brown   /* WARNING: This sharing is dangerous if either A or B is later assembled */
46c4762a1bSJed Brown   ierr = MatCreateSeqBAIJWithArrays(PETSC_COMM_SELF,1,5,5,(PetscInt*)ai,(PetscInt*)aj,aa,&B);CHKERRQ(ierr);
47c4762a1bSJed Brown   ierr = MatSeqAIJRestoreArray(A,&aa);CHKERRQ(ierr);
48c4762a1bSJed Brown   ierr = MatRestoreRowIJ(A,0,PETSC_FALSE,PETSC_FALSE,&n,&ai,&aj,&flg);CHKERRQ(ierr);
49c4762a1bSJed Brown   ierr = MatMultEqual(A,B,10,&flg);CHKERRQ(ierr);
50*2c71b3e2SJacob Faibussowitsch   PetscCheckFalse(!flg,PETSC_COMM_SELF,PETSC_ERR_ARG_NOTSAMETYPE,"MatMult(A,B) are NOT equal");
51c4762a1bSJed Brown 
52c4762a1bSJed Brown   /* Create a sbaij matrix using MatCreateSeqSBAIJWithArrays() */
53c4762a1bSJed Brown   ierr = MatCreateSeqSBAIJWithArrays(PETSC_COMM_SELF,1,5,5,asi,asj,asa,&As);CHKERRQ(ierr);
54c4762a1bSJed Brown   ierr = MatMultEqual(A,As,10,&flg);CHKERRQ(ierr);
55*2c71b3e2SJacob Faibussowitsch   PetscCheckFalse(!flg,PETSC_COMM_SELF,PETSC_ERR_ARG_NOTSAMETYPE,"MatMult(A,As) are NOT equal");
56c4762a1bSJed Brown 
57c4762a1bSJed Brown   /* Free spaces */
58c4762a1bSJed Brown   ierr = PetscRandomDestroy(&rctx);CHKERRQ(ierr);
59c4762a1bSJed Brown   ierr = MatDestroy(&A);CHKERRQ(ierr);
60c4762a1bSJed Brown   ierr = MatDestroy(&B);CHKERRQ(ierr);
61c4762a1bSJed Brown   ierr = MatDestroy(&As);CHKERRQ(ierr);
62c4762a1bSJed Brown   ierr = PetscFinalize();
63c4762a1bSJed Brown   return ierr;
64c4762a1bSJed Brown }
65