xref: /petsc/src/ksp/pc/impls/telescope/telescope.h (revision 2b8d69ca7ea5fe9190df62c1dce3bbd66fce84dd)
1 
2 #ifndef __PETSCPC_TELESCOPE_H
3 #define __PETSCPC_TELESCOPE_H
4 
5 /* Telescope */
6 typedef enum { TELESCOPE_DEFAULT = 0, TELESCOPE_DMDA, TELESCOPE_DMPLEX } PCTelescopeType;
7 
8 typedef struct _PC_Telescope *PC_Telescope;
9 struct _PC_Telescope {
10   PetscSubcomm      psubcomm;
11   PetscInt          redfactor; /* factor to reduce comm size by */
12   KSP               ksp;
13   IS                isin;
14   VecScatter        scatter;
15   Vec               xred,yred,xtmp;
16   Mat               Bred;
17   PetscBool         ignore_dm;
18   PCTelescopeType   sr_type;
19   void              *dm_ctx;
20   PetscErrorCode    (*pctelescope_setup_type)(PC,PC_Telescope);
21   PetscErrorCode    (*pctelescope_matcreate_type)(PC,PC_Telescope,MatReuse,Mat*);
22   PetscErrorCode    (*pctelescope_matnullspacecreate_type)(PC,PC_Telescope,Mat);
23   PetscErrorCode    (*pctelescope_reset_type)(PC);
24 };
25 
26 PetscBool isActiveRank(PetscSubcomm);
27 DM private_PCTelescopeGetSubDM(PC_Telescope);
28 
29 /* DMDA */
30 typedef struct {
31   DM              dmrepart;
32   Mat             permutation;
33   Vec             xp;
34   PetscInt        Mp_re,Np_re,Pp_re;
35   PetscInt        *range_i_re,*range_j_re,*range_k_re;
36   PetscInt        *start_i_re,*start_j_re,*start_k_re;
37 } PC_Telescope_DMDACtx;
38 
39 PetscErrorCode _DMDADetermineRankFromGlobalIJK(PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,
40                                                PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*,
41                                                PetscMPIInt*,PetscMPIInt*,PetscMPIInt*,PetscMPIInt*);
42 
43 PetscErrorCode _DMDADetermineGlobalS0(PetscInt,PetscMPIInt,PetscInt,PetscInt,PetscInt,PetscInt*,PetscInt*,PetscInt*,PetscInt*);
44 
45 PetscErrorCode PCTelescopeSetUp_dmda(PC,PC_Telescope);
46 PetscErrorCode PCTelescopeMatCreate_dmda(PC,PC_Telescope,MatReuse,Mat*);
47 PetscErrorCode PCTelescopeMatNullSpaceCreate_dmda(PC,PC_Telescope,Mat);
48 PetscErrorCode PCApply_Telescope_dmda(PC,Vec,Vec);
49 PetscErrorCode PCReset_Telescope_dmda(PC);
50 PetscErrorCode DMView_DMDAShort(DM,PetscViewer);
51 
52 #endif
53