Actual source code: ex76.c

petsc-master 2019-11-13
Report Typos and Errors
  1:  #include <petscksp.h>

  3: static char help[] = "Solves a linear system using PCHPDDM.\n\n";

  5: int main(int argc,char **args)
  6: {
  7:   Vec            x,b;        /* computed solution and RHS */
  8:   Mat            A,aux;      /* linear system matrix */
  9:   KSP            ksp;        /* linear solver context */
 10:   PC             pc;
 11:   IS             is,sizes;
 12:   const PetscInt *idx;
 13:   PetscMPIInt    rank,size;
 14:   const char     *deft = MATAIJ;
 15:   PetscBool      flg;
 16:   char           dir[PETSC_MAX_PATH_LEN],name[PETSC_MAX_PATH_LEN],type[256];
 17:   PetscViewer    viewer;

 20:   PetscInitialize(&argc,&args,NULL,help);if (ierr) return ierr;
 21:   MPI_Comm_size(PETSC_COMM_WORLD,&size);
 22:   if (size != 4) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_USER,"This example requires 4 processes");
 23:   MPI_Comm_rank(PETSC_COMM_WORLD,&rank);
 24:   MatCreate(PETSC_COMM_WORLD,&A);
 25:   MatCreate(PETSC_COMM_SELF,&aux);
 26:   ISCreate(PETSC_COMM_SELF,&is);
 27:   PetscStrcpy(dir,".");
 28:   PetscOptionsGetString(NULL,NULL,"-load_dir",dir,sizeof(dir),NULL);
 29:   /* loading matrices */
 30:   PetscSNPrintf(name,sizeof(name),"%s/sizes_%d_%d.dat",dir,rank,size);
 31:   PetscViewerBinaryOpen(PETSC_COMM_SELF,name,FILE_MODE_READ,&viewer);
 32:   ISCreate(PETSC_COMM_SELF,&sizes);
 33:   ISLoad(sizes,viewer);
 34:   ISGetIndices(sizes,&idx);
 35:   MatSetSizes(A,idx[0],idx[1],idx[2],idx[3]);
 36:   ISRestoreIndices(sizes,&idx);
 37:   ISDestroy(&sizes);
 38:   PetscViewerDestroy(&viewer);
 39:   MatSetUp(A);
 40:   PetscSNPrintf(name,sizeof(name),"%s/A.dat",dir);
 41:   PetscViewerBinaryOpen(PETSC_COMM_WORLD,name,FILE_MODE_READ,&viewer);
 42:   MatLoad(A,viewer);
 43:   PetscViewerDestroy(&viewer);
 44:   PetscSNPrintf(name,sizeof(name),"%s/is_%d_%d.dat",dir,rank,size);
 45:   PetscViewerBinaryOpen(PETSC_COMM_SELF,name,FILE_MODE_READ,&viewer);
 46:   ISLoad(is,viewer);
 47:   PetscViewerDestroy(&viewer);
 48:   PetscSNPrintf(name,sizeof(name),"%s/Neumann_%d_%d.dat",dir,rank,size);
 49:   PetscViewerBinaryOpen(PETSC_COMM_SELF,name,FILE_MODE_READ,&viewer);
 50:   MatSetBlockSizesFromMats(aux,A,A);
 51:   MatLoad(aux,viewer);
 52:   PetscViewerDestroy(&viewer);
 53:   MatSetOption(A,MAT_SYMMETRIC,PETSC_TRUE);
 54:   MatSetOption(aux,MAT_SYMMETRIC,PETSC_TRUE);
 55:   /* ready for testing */
 56:   PetscOptionsBegin(PETSC_COMM_WORLD,"","","");
 57:   PetscOptionsFList("-mat_type","Matrix type","MatSetType",MatList,deft,type,256,&flg);
 58:   PetscOptionsEnd();
 59:   if (flg) {
 60:     MatConvert(A,type,MAT_INPLACE_MATRIX,&A);
 61:     MatConvert(aux,type,MAT_INPLACE_MATRIX,&aux);
 62:   }
 63:   KSPCreate(PETSC_COMM_WORLD,&ksp);
 64:   KSPSetOperators(ksp,A,A);
 65:   KSPGetPC(ksp,&pc);
 66:   PCSetType(pc,PCHPDDM);
 67: #if defined(PETSC_HAVE_HPDDM)
 68:   PCHPDDMSetAuxiliaryMat(pc,is,aux,NULL,NULL);
 69:   PCHPDDMHasNeumannMat(pc,PETSC_FALSE); /* PETSC_TRUE is fine as well, just testing */
 70: #endif
 71:   ISDestroy(&is);
 72:   MatDestroy(&aux);
 73:   KSPSetFromOptions(ksp);
 74:   MatCreateVecs(A,&x,&b);
 75:   VecSet(b,1.0);
 76:   KSPSolve(ksp,b,x);
 77:   VecDestroy(&x);
 78:   VecDestroy(&b);
 79:   KSPDestroy(&ksp);
 80:   MatDestroy(&A);
 81:   PetscFinalize();
 82:   return ierr;
 83: }

 85: /*TEST

 87:    test:
 88:       requires: hpddm slepc datafilespath double !complex !define(PETSC_USE_64BIT_INDICES)
 89:       nsize: 4
 90:       args: -ksp_rtol 1e-3 -ksp_converged_reason -pc_type {{bjacobi hpddm}shared output} -pc_hpddm_coarse_sub_pc_type lu -sub_pc_type lu -options_left no -load_dir ${DATAFILESPATH}/matrices/hpddm/GENEO

 92:    test:
 93:       requires: hpddm slepc datafilespath double !complex !define(PETSC_USE_64BIT_INDICES)
 94:       suffix: geneo
 95:       nsize: 4
 96:       args: -ksp_converged_reason -pc_type hpddm -pc_hpddm_levels_1_sub_pc_type cholesky -pc_hpddm_levels_1_eps_nev {{5 15}separate output} -pc_hpddm_levels_1_st_pc_type cholesky -pc_hpddm_coarse_p {{1 2}shared output} -pc_hpddm_coarse_pc_type redundant -mat_type {{aij baij sbaij}shared output} -load_dir ${DATAFILESPATH}/matrices/hpddm/GENEO

 98:    test:
 99:       requires: hpddm slepc datafilespath double !complex !define(PETSC_USE_64BIT_INDICES)
100:       suffix: fgmres_geneo_20_p_2
101:       nsize: 4
102:       args: -ksp_converged_reason -pc_type hpddm -pc_hpddm_levels_1_sub_pc_type lu -pc_hpddm_levels_1_eps_nev 20 -pc_hpddm_coarse_p 2 -pc_hpddm_coarse_pc_type redundant -ksp_type fgmres -pc_hpddm_coarse_mat_type {{baij sbaij}shared output} -load_dir ${DATAFILESPATH}/matrices/hpddm/GENEO

104:    test:
105:       requires: hpddm slepc datafilespath double !complex !define(PETSC_USE_64BIT_INDICES)
106:       suffix: fgmres_geneo_20_p_2_geneo
107:       output_file: output/ex76_fgmres_geneo_20_p_2.out
108:       nsize: 4
109:       args: -ksp_converged_reason -pc_type hpddm -pc_hpddm_levels_1_sub_pc_type cholesky -pc_hpddm_levels_1_eps_nev 20 -pc_hpddm_levels_2_p 2 -pc_hpddm_levels_2_mat_type {{baij sbaij}shared output} -pc_hpddm_levels_2_eps_nev {{5 20}separate output} -pc_hpddm_levels_2_sub_pc_type cholesky -pc_hpddm_levels_2_ksp_type gmres -ksp_type fgmres -pc_hpddm_coarse_mat_type {{baij sbaij}shared output} -mat_type {{aij sbaij}shared output} -load_dir ${DATAFILESPATH}/matrices/hpddm/GENEO

111: TEST*/