xref: /petsc/src/mat/tests/ex208.c (revision 732aec7a18f2199fb53bb9a2f3aef439a834ce31)
1c4762a1bSJed Brown static char help[] = "Test MatCreateRedundantMatrix for rectangular matrix.\n\
2c4762a1bSJed Brown                       Contributed by Jose E. Roman, July 2017\n\n";
3c4762a1bSJed Brown 
4c4762a1bSJed Brown #include <petscmat.h>
main(int argc,char ** args)5d71ae5a4SJacob Faibussowitsch int main(int argc, char **args)
6d71ae5a4SJacob Faibussowitsch {
7c4762a1bSJed Brown   Mat         A, B;
8c4762a1bSJed Brown   PetscInt    m = 3, n = 4, i, nsubcomm;
9c4762a1bSJed Brown   PetscMPIInt size, rank;
10c4762a1bSJed Brown 
11327415f7SBarry Smith   PetscFunctionBeginUser;
12*c8025a54SPierre Jolivet   PetscCall(PetscInitialize(&argc, &args, NULL, help));
139566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
149566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
15c4762a1bSJed Brown 
16c4762a1bSJed Brown   nsubcomm = size;
179566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-nsubcomm", &nsubcomm, NULL));
18c4762a1bSJed Brown 
199566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD, &A));
209566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(A, m, n, PETSC_DETERMINE, PETSC_DETERMINE));
219566063dSJacob Faibussowitsch   PetscCall(MatSetType(A, MATAIJ));
229566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(A));
239566063dSJacob Faibussowitsch   PetscCall(MatSetUp(A));
24c4762a1bSJed Brown 
25dd400576SPatrick Sanan   if (rank == 0) {
2648a46eb9SPierre Jolivet     for (i = 0; i < size * PetscMin(m, n); i++) PetscCall(MatSetValue(A, i, i, 1.0, INSERT_VALUES));
27c4762a1bSJed Brown   }
289566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
299566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
309566063dSJacob Faibussowitsch   PetscCall(MatView(A, PETSC_VIEWER_STDOUT_WORLD));
31c4762a1bSJed Brown 
329566063dSJacob Faibussowitsch   PetscCall(MatCreateRedundantMatrix(A, nsubcomm, MPI_COMM_NULL, MAT_INITIAL_MATRIX, &B));
33c4762a1bSJed Brown   if (nsubcomm == size) { /* B is a sequential matrix */
3448a46eb9SPierre Jolivet     if (rank == 0) PetscCall(MatView(B, PETSC_VIEWER_STDOUT_SELF));
35c4762a1bSJed Brown   } else {
36c4762a1bSJed Brown     MPI_Comm comm;
379566063dSJacob Faibussowitsch     PetscCall(PetscObjectGetComm((PetscObject)B, &comm));
389566063dSJacob Faibussowitsch     PetscCall(MatView(B, PETSC_VIEWER_STDOUT_(comm)));
39c4762a1bSJed Brown   }
40c4762a1bSJed Brown 
419566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&A));
429566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&B));
439566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
44b122ec5aSJacob Faibussowitsch   return 0;
45c4762a1bSJed Brown }
46c4762a1bSJed Brown 
47c4762a1bSJed Brown /*TEST
48c4762a1bSJed Brown 
49c4762a1bSJed Brown    test:
50c4762a1bSJed Brown 
51c4762a1bSJed Brown    test:
52c4762a1bSJed Brown       suffix: 2
53c4762a1bSJed Brown       nsize: 3
54c4762a1bSJed Brown 
55c4762a1bSJed Brown    test:
56c4762a1bSJed Brown       suffix: baij
57c4762a1bSJed Brown       args: -mat_type baij
58c4762a1bSJed Brown 
59c4762a1bSJed Brown    test:
60c4762a1bSJed Brown       suffix: baij_2
61c4762a1bSJed Brown       nsize: 3
62c4762a1bSJed Brown       args: -mat_type baij
63c4762a1bSJed Brown 
64c4762a1bSJed Brown    test:
65c4762a1bSJed Brown       suffix: dense
66c4762a1bSJed Brown       args: -mat_type dense
67c4762a1bSJed Brown 
68c4762a1bSJed Brown    test:
69c4762a1bSJed Brown       suffix: dense_2
70c4762a1bSJed Brown       nsize: 3
71c4762a1bSJed Brown       args: -mat_type dense
72c4762a1bSJed Brown 
73c4762a1bSJed Brown TEST*/
74