xref: /petsc/src/mat/tests/ex233.c (revision b122ec5aa1bd4469eb4e0673542fb7de3f411254)
1 static char help[] = "Tests MatMPI{AIJ,BAIJ,SBAIJ}SetPreallocationCSR\n\n";
2 
3 #include <petscmat.h>
4 
5 int main(int argc,char **argv)
6 {
7   PetscInt       ia[3]={0,2,4};
8   PetscInt       ja[4]={0,1,0,1};
9   PetscScalar    c[16]={0,1,2,3,4,5,6,7,8,9,10,11,12,13,14,15};
10   PetscInt       ia2[5]={0,4,8,12,16};
11   PetscInt       ja2[16]={0,1,2,3,0,1,2,3,0,1,2,3,0,1,2,3};
12   PetscScalar    c2[16]={0,1,4,5,2,3,6,7,8,9,12,13,10,11,14,15};
13   PetscMPIInt    size,rank;
14   Mat            ssbaij;
15   PetscBool      rect = PETSC_FALSE;
16 
17   CHKERRQ(PetscInitialize(&argc,&argv,NULL,help));
18   CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size));
19   CHKERRMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank));
20   PetscCheckFalse(size < 2,PETSC_COMM_WORLD,PETSC_ERR_SUP,"This is an example with more then one processors");
21   if (rank) {
22     PetscInt i;
23     for (i = 0; i < 3; i++) ia[i] = 0;
24     for (i = 0; i < 5; i++) ia2[i] = 0;
25   }
26   CHKERRQ(PetscOptionsGetBool(NULL,NULL,"-rect",&rect,NULL));
27   CHKERRQ(MatCreate(PETSC_COMM_WORLD,&ssbaij));
28   CHKERRQ(MatSetBlockSize(ssbaij,2));
29   if (rect) {
30     CHKERRQ(MatSetType(ssbaij,MATMPIBAIJ));
31     CHKERRQ(MatSetSizes(ssbaij,4,6,PETSC_DECIDE,PETSC_DECIDE));
32   } else {
33     CHKERRQ(MatSetType(ssbaij,MATMPISBAIJ));
34     CHKERRQ(MatSetSizes(ssbaij,4,4,PETSC_DECIDE,PETSC_DECIDE));
35   }
36   CHKERRQ(MatSetFromOptions(ssbaij));
37   CHKERRQ(MatMPIAIJSetPreallocationCSR(ssbaij,ia2,ja2,c2));
38   CHKERRQ(MatMPIBAIJSetPreallocationCSR(ssbaij,2,ia,ja,c));
39   CHKERRQ(MatMPISBAIJSetPreallocationCSR(ssbaij,2,ia,ja,c));
40   CHKERRQ(MatViewFromOptions(ssbaij,NULL,"-view"));
41   CHKERRQ(MatDestroy(&ssbaij));
42   CHKERRQ(PetscFinalize());
43   return 0;
44 }
45 
46 /*TEST
47 
48   test:
49     filter: grep -v type | sed -e "s/\.//g"
50     suffix: aijbaij_csr
51     nsize: 2
52     args: -mat_type {{aij baij}} -view -rect {{0 1}}
53 
54   test:
55     filter: sed -e "s/\.//g"
56     suffix: sbaij_csr
57     nsize: 2
58     args: -mat_type sbaij -view -rect {{0 1}}
59 
60 TEST*/
61