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 PetscInt m=3,n=4,i,nsubcomm; 9c4762a1bSJed Brown PetscMPIInt size,rank; 10c4762a1bSJed Brown 11*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscInitialize(&argc,&args,(char*)0,help)); 125f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size)); 135f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank)); 14c4762a1bSJed Brown 15c4762a1bSJed Brown nsubcomm = size; 165f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-nsubcomm",&nsubcomm,NULL)); 17c4762a1bSJed Brown 185f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreate(PETSC_COMM_WORLD, &A)); 195f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetSizes(A, m, n, PETSC_DETERMINE, PETSC_DETERMINE)); 205f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetType(A, MATAIJ)); 215f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetFromOptions(A)); 225f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetUp(A)); 23c4762a1bSJed Brown 24dd400576SPatrick Sanan if (rank == 0) { 25c4762a1bSJed Brown for (i=0;i<size*PetscMin(m,n);i++) { 265f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValue(A, i, i, 1.0, INSERT_VALUES)); 27c4762a1bSJed Brown } 28c4762a1bSJed Brown } 295f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY)); 305f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY)); 315f80ce2aSJacob Faibussowitsch CHKERRQ(MatView(A,PETSC_VIEWER_STDOUT_WORLD)); 32c4762a1bSJed Brown 335f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateRedundantMatrix(A, nsubcomm, MPI_COMM_NULL, MAT_INITIAL_MATRIX, &B)); 34c4762a1bSJed Brown if (nsubcomm==size) { /* B is a sequential matrix */ 35dd400576SPatrick Sanan if (rank == 0) { 365f80ce2aSJacob Faibussowitsch CHKERRQ(MatView(B,PETSC_VIEWER_STDOUT_SELF)); 37c4762a1bSJed Brown } 38c4762a1bSJed Brown } else { 39c4762a1bSJed Brown MPI_Comm comm; 405f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectGetComm((PetscObject)B,&comm)); 415f80ce2aSJacob Faibussowitsch CHKERRQ(MatView(B,PETSC_VIEWER_STDOUT_(comm))); 42c4762a1bSJed Brown } 43c4762a1bSJed Brown 445f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&A)); 455f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&B)); 46*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscFinalize()); 47*b122ec5aSJacob Faibussowitsch return 0; 48c4762a1bSJed Brown } 49c4762a1bSJed Brown 50c4762a1bSJed Brown /*TEST 51c4762a1bSJed Brown 52c4762a1bSJed Brown test: 53c4762a1bSJed Brown 54c4762a1bSJed Brown test: 55c4762a1bSJed Brown suffix: 2 56c4762a1bSJed Brown nsize: 3 57c4762a1bSJed Brown 58c4762a1bSJed Brown test: 59c4762a1bSJed Brown suffix: baij 60c4762a1bSJed Brown args: -mat_type baij 61c4762a1bSJed Brown 62c4762a1bSJed Brown test: 63c4762a1bSJed Brown suffix: baij_2 64c4762a1bSJed Brown nsize: 3 65c4762a1bSJed Brown args: -mat_type baij 66c4762a1bSJed Brown 67c4762a1bSJed Brown test: 68c4762a1bSJed Brown suffix: dense 69c4762a1bSJed Brown args: -mat_type dense 70c4762a1bSJed Brown 71c4762a1bSJed Brown test: 72c4762a1bSJed Brown suffix: dense_2 73c4762a1bSJed Brown nsize: 3 74c4762a1bSJed Brown args: -mat_type dense 75c4762a1bSJed Brown 76c4762a1bSJed Brown TEST*/ 77