xref: /petsc/src/ksp/pc/impls/telescope/telescope.h (revision ccb4e88a40f0b86eaeca07ff64c64e4de2fae686)
1 
2 #ifndef PETSCPC_TELESCOPE_H
3 #define PETSCPC_TELESCOPE_H
4 
5 /* Telescope */
6 typedef enum { TELESCOPE_DEFAULT = 0, TELESCOPE_DMDA, TELESCOPE_DMPLEX, TELESCOPE_COARSEDM } PCTelescopeType;
7 
8 typedef struct _PC_Telescope *PC_Telescope;
9 struct _PC_Telescope {
10   PetscSubcomm      psubcomm;
11   PetscSubcommType  subcommtype;
12   MPI_Comm          subcomm;
13   PetscInt          redfactor; /* factor to reduce comm size by */
14   KSP               ksp;
15   IS                isin;
16   VecScatter        scatter;
17   Vec               xred,yred,xtmp;
18   Mat               Bred;
19   PetscBool         ignore_dm,ignore_kspcomputeoperators,use_coarse_dm;
20   PCTelescopeType   sr_type;
21   void              *dm_ctx;
22   PetscErrorCode    (*pctelescope_setup_type)(PC,PC_Telescope);
23   PetscErrorCode    (*pctelescope_matcreate_type)(PC,PC_Telescope,MatReuse,Mat*);
24   PetscErrorCode    (*pctelescope_matnullspacecreate_type)(PC,PC_Telescope,Mat);
25   PetscErrorCode    (*pctelescope_reset_type)(PC);
26 };
27 
28 /* DMDA */
29 typedef struct {
30   DM              dmrepart;
31   Mat             permutation;
32   Vec             xp;
33   PetscInt        Mp_re,Np_re,Pp_re;
34   PetscInt        *range_i_re,*range_j_re,*range_k_re;
35   PetscInt        *start_i_re,*start_j_re,*start_k_re;
36 } PC_Telescope_DMDACtx;
37 
38 PETSC_STATIC_INLINE PetscBool PetscSubcomm_isActiveRank(PetscSubcomm scomm)
39 {
40   if (scomm->color == 0) return(PETSC_TRUE);
41   else return(PETSC_FALSE);
42 }
43 
44 PETSC_STATIC_INLINE PetscBool PCTelescope_isActiveRank(PC_Telescope sred)
45 {
46   if (sred->psubcomm) return(PetscSubcomm_isActiveRank(sred->psubcomm));
47   else {
48     if (sred->subcomm != MPI_COMM_NULL) return(PETSC_TRUE);
49     else return (PETSC_FALSE);
50   }
51 }
52 
53 PetscErrorCode PCTelescopeSetUp_dmda(PC,PC_Telescope);
54 PetscErrorCode PCTelescopeMatCreate_dmda(PC,PC_Telescope,MatReuse,Mat*);
55 PetscErrorCode PCTelescopeMatNullSpaceCreate_dmda(PC,PC_Telescope,Mat);
56 PetscErrorCode PCApply_Telescope_dmda(PC,Vec,Vec);
57 PetscErrorCode PCApplyRichardson_Telescope_dmda(PC,Vec,Vec,Vec,PetscReal,PetscReal,PetscReal,PetscInt,PetscBool,PetscInt*,PCRichardsonConvergedReason*);
58 PetscErrorCode PCReset_Telescope_dmda(PC);
59 PetscErrorCode PCTelescopeSetUp_CoarseDM(PC,PC_Telescope);
60 PetscErrorCode PCApply_Telescope_CoarseDM(PC,Vec,Vec);
61 PetscErrorCode PCTelescopeMatNullSpaceCreate_CoarseDM(PC,PC_Telescope,Mat);
62 PetscErrorCode PCReset_Telescope_CoarseDM(PC);
63 PetscErrorCode PCApplyRichardson_Telescope_CoarseDM(PC,Vec,Vec,Vec,PetscReal,PetscReal,PetscReal,PetscInt,PetscBool,PetscInt*,PCRichardsonConvergedReason*);
64 PetscErrorCode DMView_DA_Short(DM,PetscViewer);
65 
66 #endif
67