Actual source code: ex78.c

  1: #include <petsc.h>

  3: static char help[] = "Exercises switching back and forth between different KSP and KSPHPDDM types.\n\n";

  5: int main(int argc,char **args)
  6: {
  7:   KSP            ksp;
  8: #if defined(PETSC_HAVE_HPDDM)
  9:   KSPHPDDMType   type;
 10:   PetscBool      flg;
 11: #endif
 12:   PetscInt       i;
 13:   const char     *common[] = {KSPGMRES,KSPCG,KSPPREONLY};

 16:   PetscInitialize(&argc,&args,NULL,help);if (ierr) return ierr;
 17:   KSPCreate(PETSC_COMM_WORLD,&ksp);
 18:   for (i=0; i<3; i++) {
 19:     KSPSetType(ksp,common[i]);
 20:     KSPSetType(ksp,KSPHPDDM);
 21: #if defined(PETSC_HAVE_HPDDM)
 22:     KSPHPDDMGetType(ksp,&type);
 23:     PetscStrcmp(KSPHPDDMTypes[type],common[i],&flg);
 24:     if (!flg) SETERRQ2(PetscObjectComm((PetscObject)ksp),PETSC_ERR_PLIB,"KSPType and KSPHPDDMType do not match: %s != %s", common[i], type);
 25:     KSPSetFromOptions(ksp);
 26:     KSPHPDDMGetType(ksp,&type);
 27:     if (type != KSP_HPDDM_TYPE_GCRODR) SETERRQ1(PetscObjectComm((PetscObject)ksp),PETSC_ERR_PLIB,"-ksp_hpddm_type gcrodr and KSPHPDDMType do not match: gcrodr != %s", KSPHPDDMTypes[type]);
 28:     KSPHPDDMSetType(ksp,KSP_HPDDM_TYPE_BGMRES);
 29: #endif
 30:   }
 31:   KSPDestroy(&ksp);
 32:   PetscFinalize();
 33:   return ierr;
 34: }

 36: /*TEST

 38:    test:
 39:       requires: hpddm
 40:       nsize: 1
 41:       suffix: 1
 42:       output_file: output/ex77_preonly.out
 43:       args: -ksp_type hpddm -ksp_hpddm_type gcrodr

 45: TEST*/