xref: /petsc/src/mat/tests/ex234.c (revision 2f613bf53f46f9356e00a2ca2bd69453be72fc31)
1 static char help[] = "Basic test of various routines with SBAIJ matrices\n\n";
2 
3 #include <petscmat.h>
4 
5 int main(int argc,char **argv)
6 {
7   PetscErrorCode ierr;
8   PetscInt       ia[3]={0,2,4};
9   PetscInt       ja[4]={0,1,0,1};
10   PetscScalar    c[16]={0,1,2,3,4,5,6,7,8,9,10,11,12,13,14,15};
11   PetscMPIInt    size;
12   Mat            ssbaij,msbaij;
13   Vec            x,y;
14 
15   ierr = PetscInitialize(&argc,&argv,NULL,help);if (ierr) return ierr;
16   ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRMPI(ierr);
17   if (size != 2) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"This is an example with two processors only!");
18   ierr = MatCreate(PETSC_COMM_SELF,&ssbaij);CHKERRQ(ierr);
19   ierr = MatSetType(ssbaij,MATSEQSBAIJ);CHKERRQ(ierr);
20   ierr = MatSetBlockSize(ssbaij,2);CHKERRQ(ierr);
21   ierr = MatSetSizes(ssbaij,4,8,4,8);CHKERRQ(ierr);
22   ierr = MatSeqSBAIJSetPreallocationCSR(ssbaij,2,ia,ja,c);CHKERRQ(ierr);
23   ierr = MatCreateMPIMatConcatenateSeqMat(PETSC_COMM_WORLD,ssbaij,PETSC_DECIDE,MAT_INITIAL_MATRIX,&msbaij);CHKERRQ(ierr);
24   ierr = MatView(msbaij,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
25   ierr = MatDestroy(&msbaij);CHKERRQ(ierr);
26   ierr = MatCreateMPIMatConcatenateSeqMat(PETSC_COMM_WORLD,ssbaij,4,MAT_INITIAL_MATRIX,&msbaij);CHKERRQ(ierr);
27   ierr = MatView(msbaij,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
28   ierr = MatCreateVecs(msbaij,&x,&y);CHKERRQ(ierr);
29   ierr = VecSet(x,1);CHKERRQ(ierr);
30   ierr = MatMult(msbaij,x,y);CHKERRQ(ierr);
31   ierr = VecView(y,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
32   ierr = MatMultAdd(msbaij,x,x,y);CHKERRQ(ierr);
33   ierr = VecView(y,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
34   ierr = MatGetDiagonal(msbaij,y);CHKERRQ(ierr);
35   ierr = VecView(y,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
36   ierr = VecDestroy(&x);CHKERRQ(ierr);
37   ierr = VecDestroy(&y);CHKERRQ(ierr);
38   ierr = MatDestroy(&msbaij);CHKERRQ(ierr);
39   ierr = MatDestroy(&ssbaij);CHKERRQ(ierr);
40   ierr = PetscFinalize();
41   return ierr;
42 }
43 
44 /*TEST
45 
46    test:
47      nsize: 2
48      filter: sed "s?\.??g"
49 
50 TEST*/
51