xref: /petsc/src/ksp/pc/impls/telescope/telescope.h (revision 8d9f7141f5118b8e7f737f854c827afb4e7eecc0)
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  PetscBool isActiveRank(PC_Telescope);
29  DM private_PCTelescopeGetSubDM(PC_Telescope);
30 
31 /* DMDA */
32 typedef struct {
33   DM              dmrepart;
34   Mat             permutation;
35   Vec             xp;
36   PetscInt        Mp_re,Np_re,Pp_re;
37   PetscInt        *range_i_re,*range_j_re,*range_k_re;
38   PetscInt        *start_i_re,*start_j_re,*start_k_re;
39 } PC_Telescope_DMDACtx;
40 
41  PetscErrorCode _DMDADetermineRankFromGlobalIJK(PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,
42                                                PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*,
43                                                PetscMPIInt*,PetscMPIInt*,PetscMPIInt*,PetscMPIInt*);
44 
45  PetscErrorCode _DMDADetermineGlobalS0(PetscInt,PetscMPIInt,PetscInt,PetscInt,PetscInt,PetscInt*,PetscInt*,PetscInt*,PetscInt*);
46 
47  PetscErrorCode PCTelescopeSetUp_dmda(PC,PC_Telescope);
48  PetscErrorCode PCTelescopeMatCreate_dmda(PC,PC_Telescope,MatReuse,Mat*);
49  PetscErrorCode PCTelescopeMatNullSpaceCreate_dmda(PC,PC_Telescope,Mat);
50  PetscErrorCode PCApply_Telescope_dmda(PC,Vec,Vec);
51 PetscErrorCode PCApplyRichardson_Telescope_dmda(PC,Vec,Vec,Vec,PetscReal,PetscReal,PetscReal,PetscInt,PetscBool,PetscInt*,PCRichardsonConvergedReason*);
52 PetscErrorCode PCReset_Telescope_dmda(PC);
53 PetscErrorCode PCTelescopeSetUp_CoarseDM(PC,PC_Telescope);
54 PetscErrorCode PCApply_Telescope_CoarseDM(PC,Vec,Vec);
55 PetscErrorCode PCTelescopeMatNullSpaceCreate_CoarseDM(PC,PC_Telescope,Mat);
56 PetscErrorCode PCReset_Telescope_CoarseDM(PC);
57 PetscErrorCode PCApplyRichardson_Telescope_CoarseDM(PC,Vec,Vec,Vec,PetscReal,PetscReal,PetscReal,PetscInt,PetscBool,PetscInt*,PCRichardsonConvergedReason*);
58 PetscErrorCode DMView_DA_Short(DM,PetscViewer);
59 
60 #endif
61