Actual source code: ex75.c

petsc-3.13.6 2020-09-29
Report Typos and Errors
  1:  #include <petsc.h>

  3: static char help[] = "Solves a series of linear systems using KSPHPDDM.\n\n";

  5: int main(int argc,char **args)
  6: {
  7:   Vec            x,b;        /* computed solution and RHS */
  8:   Mat            A;          /* linear system matrix */
  9:   KSP            ksp;        /* linear solver context */
 10: #if defined(PETSC_HAVE_HPDDM)
 11:   Mat            U;          /* deflation space */
 12:   PetscBool      flg;
 13: #endif
 14:   PetscInt       i,j,nmat = 10;
 15:   PetscViewer    viewer;
 16:   char           dir[PETSC_MAX_PATH_LEN],name[256];
 17:   PetscBool      reset = PETSC_FALSE;

 20:   PetscInitialize(&argc,&args,NULL,help);if (ierr) return ierr;
 21:   PetscStrcpy(dir,".");
 22:   PetscOptionsGetString(NULL,NULL,"-load_dir",dir,sizeof(dir),NULL);
 23:   PetscOptionsGetInt(NULL,NULL,"-nmat",&nmat,NULL);
 24:   PetscOptionsGetBool(NULL,NULL,"-reset",&reset,NULL);
 25:   MatCreate(PETSC_COMM_WORLD,&A);
 26:   KSPCreate(PETSC_COMM_WORLD,&ksp);
 27:   KSPSetOperators(ksp,A,A);
 28:   for (i=0; i<nmat; i++) {
 29:     j = i+400;
 30:     PetscSNPrintf(name,sizeof(name),"%s/A_%d.dat",dir,j);
 31:     PetscViewerBinaryOpen(PETSC_COMM_WORLD,name,FILE_MODE_READ,&viewer);
 32:     MatLoad(A,viewer);
 33:     PetscViewerDestroy(&viewer);
 34:     if (i == 0) {
 35:       MatCreateVecs(A,&x,&b);
 36:     }
 37:     PetscSNPrintf(name,sizeof(name),"%s/rhs_%d.dat",dir,j);
 38:     PetscViewerBinaryOpen(PETSC_COMM_WORLD,name,FILE_MODE_READ,&viewer);
 39:     VecLoad(b,viewer);
 40:     PetscViewerDestroy(&viewer);
 41:     KSPSetFromOptions(ksp);
 42:     KSPSolve(ksp,b,x);
 43: #if defined(PETSC_HAVE_HPDDM)
 44:     PetscObjectTypeCompare((PetscObject)ksp,KSPHPDDM,&flg);
 45:     if (flg && reset) {
 46:       KSPHPDDMGetDeflationSpace(ksp,&U);
 47:       KSPReset(ksp);
 48:       KSPSetOperators(ksp,A,A);
 49:       KSPSetFromOptions(ksp);
 50:       KSPSetUp(ksp);
 51:       if (U) {
 52:         KSPHPDDMSetDeflationSpace(ksp,U);
 53:         MatDestroy(&U);
 54:       }
 55:     }
 56: #endif
 57:   }
 58:   VecDestroy(&x);
 59:   VecDestroy(&b);
 60:   MatDestroy(&A);
 61:   KSPDestroy(&ksp);
 62:   PetscFinalize();
 63:   return ierr;
 64: }

 66: /*TEST

 68:    test:
 69:       suffix: 1
 70:       nsize: 1
 71:       requires: hpddm datafilespath double !complex !define(PETSC_USE_64BIT_INDICES)
 72:       args: -nmat 1 -pc_type none -ksp_converged_reason -ksp_type {{gmres hpddm}shared output} -ksp_max_it 1000 -ksp_gmres_restart 1000 -ksp_rtol 1e-10 -ksp_hpddm_type {{gmres bgmres}shared output} -options_left no -load_dir ${DATAFILESPATH}/matrices/hpddm/GCRODR

 74:    test:
 75:       requires: hpddm datafilespath double !complex !define(PETSC_USE_64BIT_INDICES)
 76:       suffix: 1_icc
 77:       nsize: 1
 78:       args: -nmat 1 -pc_type icc -ksp_converged_reason -ksp_type {{gmres hpddm}shared output} -ksp_max_it 1000 -ksp_gmres_restart 1000 -ksp_rtol 1e-10 -load_dir ${DATAFILESPATH}/matrices/hpddm/GCRODR

 80:    testset:
 81:       requires: hpddm datafilespath double !complex !define(PETSC_USE_64BIT_INDICES)
 82:       args: -nmat 3 -pc_type none -ksp_converged_reason -ksp_type hpddm -ksp_max_it 1000 -ksp_gmres_restart 40 -ksp_rtol 1e-10 -ksp_hpddm_type {{gcrodr bgcrodr}shared output} -ksp_hpddm_recycle 20 -load_dir ${DATAFILESPATH}/matrices/hpddm/GCRODR
 83:       test:
 84:         nsize: 1
 85:         suffix: 2_seq
 86:         output_file: output/ex75_2.out
 87:       test:
 88:         nsize: 2
 89:         suffix: 2_par
 90:         output_file: output/ex75_2.out

 92:    test:
 93:       requires: hpddm datafilespath double !complex !define(PETSC_USE_64BIT_INDICES)
 94:       suffix: 2_icc
 95:       nsize: 1
 96:       args: -nmat 3 -pc_type icc -ksp_converged_reason -ksp_type hpddm -ksp_max_it 1000 -ksp_gmres_restart 40 -ksp_rtol 1e-10 -ksp_hpddm_type gcrodr -ksp_hpddm_recycle 20 -reset {{false true}shared output} -load_dir ${DATAFILESPATH}/matrices/hpddm/GCRODR

 98: TEST*/