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