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> 59371c9d4SSatish Balay int main(int argc, char **args) { 6c4762a1bSJed Brown Mat A, B; 7c4762a1bSJed Brown PetscInt m = 3, n = 4, i, nsubcomm; 8c4762a1bSJed Brown PetscMPIInt size, rank; 9c4762a1bSJed Brown 10327415f7SBarry Smith PetscFunctionBeginUser; 119566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &args, (char *)0, help)); 129566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size)); 139566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank)); 14c4762a1bSJed Brown 15c4762a1bSJed Brown nsubcomm = size; 169566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-nsubcomm", &nsubcomm, NULL)); 17c4762a1bSJed Brown 189566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &A)); 199566063dSJacob Faibussowitsch PetscCall(MatSetSizes(A, m, n, PETSC_DETERMINE, PETSC_DETERMINE)); 209566063dSJacob Faibussowitsch PetscCall(MatSetType(A, MATAIJ)); 219566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(A)); 229566063dSJacob Faibussowitsch PetscCall(MatSetUp(A)); 23c4762a1bSJed Brown 24dd400576SPatrick Sanan if (rank == 0) { 25*48a46eb9SPierre Jolivet for (i = 0; i < size * PetscMin(m, n); i++) PetscCall(MatSetValue(A, i, i, 1.0, INSERT_VALUES)); 26c4762a1bSJed Brown } 279566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY)); 289566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY)); 299566063dSJacob Faibussowitsch PetscCall(MatView(A, PETSC_VIEWER_STDOUT_WORLD)); 30c4762a1bSJed Brown 319566063dSJacob Faibussowitsch PetscCall(MatCreateRedundantMatrix(A, nsubcomm, MPI_COMM_NULL, MAT_INITIAL_MATRIX, &B)); 32c4762a1bSJed Brown if (nsubcomm == size) { /* B is a sequential matrix */ 33*48a46eb9SPierre Jolivet if (rank == 0) PetscCall(MatView(B, PETSC_VIEWER_STDOUT_SELF)); 34c4762a1bSJed Brown } else { 35c4762a1bSJed Brown MPI_Comm comm; 369566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)B, &comm)); 379566063dSJacob Faibussowitsch PetscCall(MatView(B, PETSC_VIEWER_STDOUT_(comm))); 38c4762a1bSJed Brown } 39c4762a1bSJed Brown 409566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 419566063dSJacob Faibussowitsch PetscCall(MatDestroy(&B)); 429566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 43b122ec5aSJacob Faibussowitsch return 0; 44c4762a1bSJed Brown } 45c4762a1bSJed Brown 46c4762a1bSJed Brown /*TEST 47c4762a1bSJed Brown 48c4762a1bSJed Brown test: 49c4762a1bSJed Brown 50c4762a1bSJed Brown test: 51c4762a1bSJed Brown suffix: 2 52c4762a1bSJed Brown nsize: 3 53c4762a1bSJed Brown 54c4762a1bSJed Brown test: 55c4762a1bSJed Brown suffix: baij 56c4762a1bSJed Brown args: -mat_type baij 57c4762a1bSJed Brown 58c4762a1bSJed Brown test: 59c4762a1bSJed Brown suffix: baij_2 60c4762a1bSJed Brown nsize: 3 61c4762a1bSJed Brown args: -mat_type baij 62c4762a1bSJed Brown 63c4762a1bSJed Brown test: 64c4762a1bSJed Brown suffix: dense 65c4762a1bSJed Brown args: -mat_type dense 66c4762a1bSJed Brown 67c4762a1bSJed Brown test: 68c4762a1bSJed Brown suffix: dense_2 69c4762a1bSJed Brown nsize: 3 70c4762a1bSJed Brown args: -mat_type dense 71c4762a1bSJed Brown 72c4762a1bSJed Brown TEST*/ 73