xref: /petsc/src/mat/tests/ex233.c (revision 16ce926e8e1f977214c72b888607cdae2c3b4104)
1c4762a1bSJed Brown static char help[] = "Tests MatMPI{AIJ,BAIJ,SBAIJ}SetPreallocationCSR\n\n";
2c4762a1bSJed Brown 
3c4762a1bSJed Brown #include <petscmat.h>
4c4762a1bSJed Brown 
main(int argc,char ** argv)5d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
6d71ae5a4SJacob Faibussowitsch {
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 
17327415f7SBarry Smith   PetscFunctionBeginUser;
189566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
199566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
209566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
21*467446fbSPierre Jolivet   PetscCheck(size > 1, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "This is an example with more than one process");
22c4762a1bSJed Brown   if (rank) {
23c4762a1bSJed Brown     PetscInt i;
24c4762a1bSJed Brown     for (i = 0; i < 3; i++) ia[i] = 0;
25c4762a1bSJed Brown     for (i = 0; i < 5; i++) ia2[i] = 0;
26c4762a1bSJed Brown   }
279566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-rect", &rect, NULL));
289566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD, &ssbaij));
299566063dSJacob Faibussowitsch   PetscCall(MatSetBlockSize(ssbaij, 2));
30c4762a1bSJed Brown   if (rect) {
319566063dSJacob Faibussowitsch     PetscCall(MatSetType(ssbaij, MATMPIBAIJ));
329566063dSJacob Faibussowitsch     PetscCall(MatSetSizes(ssbaij, 4, 6, PETSC_DECIDE, PETSC_DECIDE));
33c4762a1bSJed Brown   } else {
349566063dSJacob Faibussowitsch     PetscCall(MatSetType(ssbaij, MATMPISBAIJ));
359566063dSJacob Faibussowitsch     PetscCall(MatSetSizes(ssbaij, 4, 4, PETSC_DECIDE, PETSC_DECIDE));
36c4762a1bSJed Brown   }
379566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(ssbaij));
389566063dSJacob Faibussowitsch   PetscCall(MatMPIAIJSetPreallocationCSR(ssbaij, ia2, ja2, c2));
399566063dSJacob Faibussowitsch   PetscCall(MatMPIBAIJSetPreallocationCSR(ssbaij, 2, ia, ja, c));
409566063dSJacob Faibussowitsch   PetscCall(MatMPISBAIJSetPreallocationCSR(ssbaij, 2, ia, ja, c));
419566063dSJacob Faibussowitsch   PetscCall(MatViewFromOptions(ssbaij, NULL, "-view"));
429566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&ssbaij));
439566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
44b122ec5aSJacob Faibussowitsch   return 0;
45c4762a1bSJed Brown }
46c4762a1bSJed Brown 
47c4762a1bSJed Brown /*TEST
48c4762a1bSJed Brown 
49c4762a1bSJed Brown   test:
50c4762a1bSJed Brown     filter: grep -v type | sed -e "s/\.//g"
51c4762a1bSJed Brown     suffix: aijbaij_csr
52c4762a1bSJed Brown     nsize: 2
53c4762a1bSJed Brown     args: -mat_type {{aij baij}} -view -rect {{0 1}}
54c4762a1bSJed Brown 
55c4762a1bSJed Brown   test:
56c4762a1bSJed Brown     filter: sed -e "s/\.//g"
57c4762a1bSJed Brown     suffix: sbaij_csr
58c4762a1bSJed Brown     nsize: 2
59c4762a1bSJed Brown     args: -mat_type sbaij -view -rect {{0 1}}
60c4762a1bSJed Brown 
61c4762a1bSJed Brown TEST*/
62