Actual source code: ex208.c

petsc-3.13.6 2020-09-29
Report Typos and Errors
  1: static char help[] = "Test MatCreateRedundantMatrix for rectangular matrix.\n\
  2:                       Contributed by Jose E. Roman, July 2017\n\n";

  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;

 12:   PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr;
 13:   MPI_Comm_size(PETSC_COMM_WORLD, &size);
 14:   MPI_Comm_rank(PETSC_COMM_WORLD, &rank);

 16:   nsubcomm = size;
 17:   PetscOptionsGetInt(NULL,NULL,"-nsubcomm",&nsubcomm,NULL);

 19:   MatCreate(PETSC_COMM_WORLD, &A); 
 20:   MatSetSizes(A, m, n, PETSC_DETERMINE, PETSC_DETERMINE);
 21:   MatSetType(A, MATAIJ);
 22:   MatSetFromOptions(A);
 23:   MatSetUp(A);

 25:   if (!rank) {
 26:     for (i=0;i<size*PetscMin(m,n);i++) {
 27:       MatSetValue(A, i, i, 1.0, INSERT_VALUES); 
 28:     }
 29:   }
 30:   MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY);
 31:   MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY);
 32:   MatView(A,PETSC_VIEWER_STDOUT_WORLD);

 34:   MatCreateRedundantMatrix(A, nsubcomm, MPI_COMM_NULL, MAT_INITIAL_MATRIX, &B);
 35:   if (nsubcomm==size) { /* B is a sequential matrix */
 36:     if (!rank) {
 37:       MatView(B,PETSC_VIEWER_STDOUT_SELF);
 38:     }
 39:   } else {
 40:     MPI_Comm comm;
 41:     PetscObjectGetComm((PetscObject)B,&comm);
 42:     MatView(B,PETSC_VIEWER_STDOUT_(comm));
 43:   }

 45:   MatDestroy(&A);
 46:   MatDestroy(&B);
 47:   PetscFinalize();
 48:   return ierr;
 49: }


 52: /*TEST

 54:    test:

 56:    test:
 57:       suffix: 2
 58:       nsize: 3

 60:    test:
 61:       suffix: baij
 62:       args: -mat_type baij

 64:    test:
 65:       suffix: baij_2
 66:       nsize: 3
 67:       args: -mat_type baij

 69:    test:
 70:       suffix: dense
 71:       args: -mat_type dense

 73:    test:
 74:       suffix: dense_2
 75:       nsize: 3
 76:       args: -mat_type dense

 78: TEST*/