Actual source code: ex37.c

petsc-3.8.3 2017-12-09
Report Typos and Errors

  2: static char help[] = "Test MatGetMultiProcBlock() and MatCreateRedundantMatrix() \n\
  3: Reads a PETSc matrix and vector from a file and solves a linear system.\n\n";
  4: /*
  5:   Example:
  6:   mpiexec -n 4 ./ex37 -f <input_file> -nsubcomm 2 -subcomm_view -subcomm_type <1 or 2>
  7: */

  9:  #include <petscksp.h>

 11: int main(int argc,char **args)
 12: {
 13:   KSP            subksp;
 14:   Mat            A,subA;
 15:   Vec            x,b,u,subb,subx,subu;
 16:   PetscViewer    fd;
 17:   char           file[PETSC_MAX_PATH_LEN];
 18:   PetscBool      flg;
 20:   PetscInt       i,m,n,its;
 21:   PetscReal      norm;
 22:   PetscMPIInt    rank,size;
 23:   MPI_Comm       comm,subcomm;
 24:   PetscSubcomm   psubcomm;
 25:   PetscInt       nsubcomm=1,id;
 26:   PetscScalar    *barray,*xarray,*uarray,*array,one=1.0;
 27:   PetscInt       type=1;

 29:   PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr;
 30:   /* Load the matrix */
 31:   PetscOptionsGetString(NULL,NULL,"-f",file,PETSC_MAX_PATH_LEN,&flg);
 32:   if (!flg) SETERRQ(PETSC_COMM_WORLD,1,"Must indicate binary file with the -f option");
 33:   PetscViewerBinaryOpen(PETSC_COMM_WORLD,file,FILE_MODE_READ,&fd);

 35:   /* Load the matrix; then destroy the viewer.*/
 36:   MatCreate(PETSC_COMM_WORLD,&A);
 37:   MatLoad(A,fd);
 38:   PetscViewerDestroy(&fd);

 40:   PetscObjectGetComm((PetscObject)A,&comm);
 41:   MPI_Comm_size(comm,&size);
 42:   MPI_Comm_rank(comm,&rank);

 44:   /* Create rhs vector b */
 45:   MatGetLocalSize(A,&m,NULL);
 46:   VecCreate(PETSC_COMM_WORLD,&b);
 47:   VecSetSizes(b,m,PETSC_DECIDE);
 48:   VecSetFromOptions(b);
 49:   VecSet(b,one);

 51:   VecDuplicate(b,&x);
 52:   VecDuplicate(b,&u);
 53:   VecSet(x,0.0);

 55:   /* Test MatGetMultiProcBlock() */
 56:   PetscOptionsGetInt(NULL,NULL,"-nsubcomm",&nsubcomm,NULL);
 57:   PetscOptionsGetInt(NULL,NULL,"-subcomm_type",&type,NULL);

 59:   PetscSubcommCreate(comm,&psubcomm);
 60:   PetscSubcommSetNumber(psubcomm,nsubcomm);
 61:   if (type == PETSC_SUBCOMM_GENERAL) { /* user provides color, subrank and duprank */
 62:     PetscMPIInt color,subrank,duprank,subsize;
 63:     duprank = size-1 - rank;
 64:     subsize = size/nsubcomm;
 65:     if (subsize*nsubcomm != size) SETERRQ2(comm,PETSC_ERR_SUP,"This example requires nsubcomm %D divides size %D",nsubcomm,size);
 66:     color   = duprank/subsize;
 67:     subrank = duprank - color*subsize;
 68:     PetscSubcommSetTypeGeneral(psubcomm,color,subrank);
 69:   } else if (type == PETSC_SUBCOMM_CONTIGUOUS) {
 70:     PetscSubcommSetType(psubcomm,PETSC_SUBCOMM_CONTIGUOUS);
 71:   } else if (type == PETSC_SUBCOMM_INTERLACED) {
 72:     PetscSubcommSetType(psubcomm,PETSC_SUBCOMM_INTERLACED);
 73:   } else SETERRQ1(psubcomm->parent,PETSC_ERR_SUP,"PetscSubcommType %D is not supported yet",type);
 74:   PetscSubcommSetFromOptions(psubcomm);
 75:   subcomm = PetscSubcommChild(psubcomm);

 77:   /* Test MatCreateRedundantMatrix() */
 78:   if (size > 1) {
 79:     MatCreateRedundantMatrix(A,nsubcomm,subcomm,MAT_INITIAL_MATRIX,&subA);
 80:     MatCreateRedundantMatrix(A,nsubcomm,subcomm,MAT_REUSE_MATRIX,&subA);
 81:     MatDestroy(&subA);
 82:   }

 84:   /* Create subA */
 85:   MatGetMultiProcBlock(A,subcomm,MAT_INITIAL_MATRIX,&subA);
 86:   MatGetMultiProcBlock(A,subcomm,MAT_REUSE_MATRIX,&subA);

 88:   /* Create sub vectors without arrays. Place b's and x's local arrays into subb and subx */
 89:   MatGetLocalSize(subA,&m,&n);
 90:   VecCreateMPIWithArray(subcomm,1,m,PETSC_DECIDE,NULL,&subb);
 91:   VecCreateMPIWithArray(subcomm,1,n,PETSC_DECIDE,NULL,&subx);
 92:   VecCreateMPIWithArray(subcomm,1,n,PETSC_DECIDE,NULL,&subu);

 94:   VecGetArray(b,&barray);
 95:   VecGetArray(x,&xarray);
 96:   VecGetArray(u,&uarray);
 97:   VecPlaceArray(subb,barray);
 98:   VecPlaceArray(subx,xarray);
 99:   VecPlaceArray(subu,uarray);

101:   /* Create linear solvers associated with subA */
102:   KSPCreate(subcomm,&subksp);
103:   KSPSetOperators(subksp,subA,subA);
104:   KSPSetFromOptions(subksp);

106:   /* Solve sub systems */
107:   KSPSolve(subksp,subb,subx);
108:   KSPGetIterationNumber(subksp,&its);

110:   /* check residual */
111:   MatMult(subA,subx,subu);
112:   VecAXPY(subu,-1.0,subb);
113:   VecNorm(u,NORM_2,&norm);
114:   if (norm > 1.e-4 && !rank) {
115:     PetscPrintf(PETSC_COMM_WORLD,"[%D]  Number of iterations = %3D\n",rank,its);
116:     PetscPrintf(PETSC_COMM_WORLD,"Error: Residual norm of each block |subb - subA*subx |= %g\n",(double)norm);
117:   }
118:   VecResetArray(subb);
119:   VecResetArray(subx);
120:   VecResetArray(subu);

122:   PetscOptionsGetInt(NULL,NULL,"-subvec_view",&id,&flg);
123:   if (flg && rank == id) {
124:     PetscPrintf(PETSC_COMM_SELF,"[%D] subb:\n", rank);
125:     VecGetArray(subb,&array);
126:     for (i=0; i<m; i++) {PetscPrintf(PETSC_COMM_SELF,"%g\n",(double)PetscRealPart(array[i]));}
127:     VecRestoreArray(subb,&array);
128:     PetscPrintf(PETSC_COMM_SELF,"[%D] subx:\n", rank);
129:     VecGetArray(subx,&array);
130:     for (i=0; i<m; i++) {PetscPrintf(PETSC_COMM_SELF,"%g\n",(double)PetscRealPart(array[i]));}
131:     VecRestoreArray(subx,&array);
132:   }

134:   VecRestoreArray(x,&xarray);
135:   VecRestoreArray(b,&barray);
136:   VecRestoreArray(u,&uarray);
137:   MatDestroy(&subA);
138:   VecDestroy(&subb);
139:   VecDestroy(&subx);
140:   VecDestroy(&subu);
141:   KSPDestroy(&subksp);
142:   PetscSubcommDestroy(&psubcomm);
143:   MatDestroy(&A); VecDestroy(&b);
144:   VecDestroy(&u); VecDestroy(&x);

146:   PetscFinalize();
147:   return ierr;
148: }