xref: /petsc/src/mat/tests/ex208.c (revision 5f80ce2ab25dff0f4601e710601cbbcecf323266)
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>
5c4762a1bSJed Brown int main(int argc,char **args)
6c4762a1bSJed Brown {
7c4762a1bSJed Brown   Mat               A,B;
8c4762a1bSJed Brown   PetscErrorCode    ierr;
9c4762a1bSJed Brown   PetscInt          m=3,n=4,i,nsubcomm;
10c4762a1bSJed Brown   PetscMPIInt       size,rank;
11c4762a1bSJed Brown 
12c4762a1bSJed Brown   ierr = PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr;
13*5f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
14*5f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
15c4762a1bSJed Brown 
16c4762a1bSJed Brown   nsubcomm = size;
17*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-nsubcomm",&nsubcomm,NULL));
18c4762a1bSJed Brown 
19*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreate(PETSC_COMM_WORLD, &A));
20*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetSizes(A, m, n, PETSC_DETERMINE, PETSC_DETERMINE));
21*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetType(A, MATAIJ));
22*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetFromOptions(A));
23*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetUp(A));
24c4762a1bSJed Brown 
25dd400576SPatrick Sanan   if (rank == 0) {
26c4762a1bSJed Brown     for (i=0;i<size*PetscMin(m,n);i++) {
27*5f80ce2aSJacob Faibussowitsch       CHKERRQ(MatSetValue(A, i, i, 1.0, INSERT_VALUES));
28c4762a1bSJed Brown     }
29c4762a1bSJed Brown   }
30*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
31*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
32*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatView(A,PETSC_VIEWER_STDOUT_WORLD));
33c4762a1bSJed Brown 
34*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreateRedundantMatrix(A, nsubcomm, MPI_COMM_NULL, MAT_INITIAL_MATRIX, &B));
35c4762a1bSJed Brown   if (nsubcomm==size) { /* B is a sequential matrix */
36dd400576SPatrick Sanan     if (rank == 0) {
37*5f80ce2aSJacob Faibussowitsch       CHKERRQ(MatView(B,PETSC_VIEWER_STDOUT_SELF));
38c4762a1bSJed Brown     }
39c4762a1bSJed Brown   } else {
40c4762a1bSJed Brown     MPI_Comm comm;
41*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscObjectGetComm((PetscObject)B,&comm));
42*5f80ce2aSJacob Faibussowitsch     CHKERRQ(MatView(B,PETSC_VIEWER_STDOUT_(comm)));
43c4762a1bSJed Brown   }
44c4762a1bSJed Brown 
45*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&A));
46*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&B));
47c4762a1bSJed Brown   ierr = PetscFinalize();
48c4762a1bSJed Brown   return ierr;
49c4762a1bSJed Brown }
50c4762a1bSJed Brown 
51c4762a1bSJed Brown /*TEST
52c4762a1bSJed Brown 
53c4762a1bSJed Brown    test:
54c4762a1bSJed Brown 
55c4762a1bSJed Brown    test:
56c4762a1bSJed Brown       suffix: 2
57c4762a1bSJed Brown       nsize: 3
58c4762a1bSJed Brown 
59c4762a1bSJed Brown    test:
60c4762a1bSJed Brown       suffix: baij
61c4762a1bSJed Brown       args: -mat_type baij
62c4762a1bSJed Brown 
63c4762a1bSJed Brown    test:
64c4762a1bSJed Brown       suffix: baij_2
65c4762a1bSJed Brown       nsize: 3
66c4762a1bSJed Brown       args: -mat_type baij
67c4762a1bSJed Brown 
68c4762a1bSJed Brown    test:
69c4762a1bSJed Brown       suffix: dense
70c4762a1bSJed Brown       args: -mat_type dense
71c4762a1bSJed Brown 
72c4762a1bSJed Brown    test:
73c4762a1bSJed Brown       suffix: dense_2
74c4762a1bSJed Brown       nsize: 3
75c4762a1bSJed Brown       args: -mat_type dense
76c4762a1bSJed Brown 
77c4762a1bSJed Brown TEST*/
78