Actual source code: ex37.c

petsc-3.6.4 2016-04-12
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>

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

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

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

 42:   PetscObjectGetComm((PetscObject)A,&comm);
 43:   MPI_Comm_size(comm,&size);
 44:   MPI_Comm_rank(comm,&rank);

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

 53:   VecDuplicate(b,&x);
 54:   VecDuplicate(b,&u);
 55:   VecSet(x,0.0);

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

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

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

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

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

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

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

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

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

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

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

148:   PetscFinalize();
149:   return 0;
150: }