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