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