Actual source code: bddc.h

petsc-3.3-p7 2013-05-11

  5: #include <../src/ksp/pc/impls/is/pcis.h>

  7: /* BDDC requires metis 5.0.1 for multilevel */
  8: #include "metis.h"
  9: #define MetisInt    idx_t
 10: #define MetisScalar real_t


 13: /* Structure for graph partitioning (adapted from Metis) */
 14: struct _PCBDDCGraph {
 15:   PetscInt nvtxs;
 16:   PetscInt ncmps;
 17:   PetscInt *xadj;
 18:   PetscInt *adjncy;
 19:   PetscInt *where;
 20:   PetscInt *which_dof;
 21:   PetscInt *cptr;
 22:   PetscInt *queue;
 23:   PetscInt *count;
 24:   PetscMPIInt *where_ncmps;
 25:   PetscBool *touched;
 26: };

 28: typedef enum {SCATTERS_BDDC,GATHERS_BDDC} CoarseCommunicationsType;
 29: typedef struct _PCBDDCGraph *PCBDDCGraph;

 31: /* Private context (data structure) for the BDDC preconditioner.  */
 32: typedef struct {
 33:   /* First MUST come the folowing line, for the stuff that is common to FETI and Neumann-Neumann. */
 34:   PC_IS         pcis;
 35:   /* Coarse stuffs needed by BDDC application in KSP */
 36:   PetscInt      coarse_size;
 37:   Mat           coarse_mat;
 38:   Vec           coarse_vec;
 39:   Vec           coarse_rhs;
 40:   KSP           coarse_ksp;
 41:   Mat           coarse_phi_B;
 42:   Mat           coarse_phi_D;
 43:   PetscMPIInt   local_primal_size;
 44:   PetscMPIInt   *local_primal_indices;
 45:   PetscMPIInt   *local_primal_displacements;
 46:   PetscMPIInt   *local_primal_sizes;
 47:   PetscMPIInt   replicated_primal_size;
 48:   PetscMPIInt   *replicated_local_primal_indices;
 49:   PetscScalar   *replicated_local_primal_values;
 50:   VecScatter    coarse_loc_to_glob;
 51:   /* Local stuffs needed by BDDC application in KSP */
 52:   Vec           vec1_P;
 53:   Vec           vec1_C;
 54:   Mat           local_auxmat1;
 55:   Mat           local_auxmat2;
 56:   Vec           vec1_R;
 57:   Vec           vec2_R;
 58:   VecScatter    R_to_B;
 59:   VecScatter    R_to_D;
 60:   KSP           ksp_R;
 61:   KSP           ksp_D;
 62:   Vec           vec4_D;
 63:   /* Quantities defining constraining details (local) of the preconditioner */
 64:   /* These quantities define the preconditioner itself */
 65:   IS            *ISForFaces;
 66:   IS            *ISForEdges;
 67:   IS            ISForVertices;
 68:   PetscInt      n_ISForFaces;
 69:   PetscInt      n_ISForEdges;
 70:   PetscInt      n_constraints;
 71:   PetscInt      n_vertices;
 72:   Mat           ConstraintMatrix;
 73:   /* Some defaults on selecting vertices and constraints*/
 74:   PetscBool     vertices_flag;
 75:   PetscBool     constraints_flag;
 76:   PetscBool     faces_flag;
 77:   PetscBool     edges_flag;
 78:   /* Some customization is possible */
 79:   PetscBool                  use_nnsp_true;
 80:   PetscInt                   n_ISForDofs;
 81:   IS                         *ISForDofs;
 82:   IS                         NeumannBoundaries;
 83:   IS                         DirichletBoundaries;
 84:   PetscBool                  prec_type;
 85:   CoarseProblemType          coarse_problem_type;
 86:   CoarseCommunicationsType   coarse_communications_type;
 87:   PetscInt                   coarsening_ratio;
 88:   PetscInt                   active_procs;
 89:   /* For verbose output of some bddc data structures */
 90:   PetscBool                  dbg_flag;
 91:   PetscViewer                dbg_viewer;
 92: } PC_BDDC;

 94: /* In case of multilevel BDDC, this is the minimum number of procs for which it will be allowed */
 95: #define MIN_PROCS_FOR_BDDC 5 

 97: /* prototypes for functions contained in bddc.c */
 98: static PetscErrorCode PCBDDCCoarseSetUp(PC);
 99: static PetscErrorCode PCBDDCFindConnectedComponents(PCBDDCGraph,PetscInt);
100: static PetscErrorCode PCBDDCSetupCoarseEnvironment(PC,PetscScalar*);
101: static PetscErrorCode PCBDDCManageLocalBoundaries(PC);
102: static PetscErrorCode PCBDDCApplyInterfacePreconditioner(PC,Vec);
103: static PetscErrorCode PCBDDCSolveSaddlePoint(PC);
104: static PetscErrorCode PCBDDCScatterCoarseDataBegin(PC,Vec,Vec,InsertMode,ScatterMode);
105: static PetscErrorCode PCBDDCScatterCoarseDataEnd(PC,Vec,Vec,InsertMode,ScatterMode);
106: static PetscErrorCode PCBDDCCreateConstraintMatrix(PC);

108: #endif /* __pcbddc_h */