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