4: #include <../src/ksp/pc/impls/is/pcis.h> 5: #include <../src/ksp/pc/impls/bddc/bddcstructs.h> 7: /* Private context (data structure) for the BDDC preconditioner. */ 8: typedef struct { 9: /* First MUST come the folowing line, for the stuff that is common to FETI and Neumann-Neumann. */ 10: PC_IS pcis; 11: /* Coarse stuffs needed by BDDC application in KSP */ 12: Vec coarse_vec; 13: KSP coarse_ksp; 14: Mat coarse_phi_B; 15: Mat coarse_phi_D; 16: Mat coarse_psi_B; 17: Mat coarse_psi_D; 18: PetscInt local_primal_size; 19: PetscInt coarse_size; 20: PetscInt* global_primal_indices; 21: VecScatter coarse_loc_to_glob; 22: /* Local stuffs needed by BDDC application in KSP */ 23: Vec vec1_P; 24: Vec vec1_C; 25: Mat local_auxmat1; 26: Mat local_auxmat2; 27: Vec vec1_R; 28: Vec vec2_R; 29: IS is_R_local; 30: VecScatter R_to_B; 31: VecScatter R_to_D; 32: KSP ksp_R; 33: KSP ksp_D; 34: /* Quantities defining constraining details (local) of the preconditioner */ 35: /* These quantities define the preconditioner itself */ 36: PetscInt n_vertices; 37: Mat ConstraintMatrix; 38: PetscBool new_primal_space; 39: PetscBool new_primal_space_local; 40: PetscInt *primal_indices_local_idxs; 41: PetscInt local_primal_size_cc; 42: PetscInt *local_primal_ref_node; 43: PetscInt *local_primal_ref_mult; 44: PetscBool use_change_of_basis; 45: PetscBool use_change_on_faces; 46: PetscBool fake_change; 47: Mat ChangeOfBasisMatrix; 48: Mat user_ChangeOfBasisMatrix; 49: PetscBool change_interior; 50: Mat switch_static_change; 51: Vec work_change; 52: Vec original_rhs; 53: Vec temp_solution; 54: Mat local_mat; 55: PetscBool use_exact_dirichlet_trick; 56: PetscBool exact_dirichlet_trick_app; 57: PetscBool ksp_guess_nonzero; 58: PetscBool rhs_change; 59: PetscBool temp_solution_used; 60: /* benign subspace trick */ 61: PetscBool benign_saddle_point; 62: PetscBool benign_have_null; 63: PetscBool benign_skip_correction; 64: PetscBool benign_compute_correction; 65: Mat benign_change; 66: Mat benign_original_mat; 67: IS *benign_zerodiag_subs; 68: Vec benign_vec; 69: Mat benign_B0; 70: PetscSF benign_sf; 71: PetscScalar *benign_p0; 72: PetscInt benign_n; 73: PetscInt *benign_p0_lidx; 74: PetscInt *benign_p0_gidx; 75: PetscBool benign_null; 76: PetscBool benign_change_explicit; 77: PetscBool benign_apply_coarse_only; 79: /* Some defaults on selecting vertices and constraints*/ 80: PetscBool use_local_adj; 81: PetscBool use_vertices; 82: PetscBool use_faces; 83: PetscBool use_edges; 84: /* Some customization is possible */ 85: PetscBool recompute_topography; 86: PetscBool graphanalyzed; 87: PCBDDCGraph mat_graph; 88: PetscInt graphmaxcount; 89: MatNullSpace onearnullspace; 90: PetscObjectState *onearnullvecs_state; 91: PetscBool NullSpace_corr[4]; 92: IS user_primal_vertices; 93: IS user_primal_vertices_local; 94: PetscBool use_nnsp_true; 95: PetscBool use_qr_single; 96: PetscBool user_provided_isfordofs; 97: PetscInt n_ISForDofs; 98: PetscInt n_ISForDofsLocal; 99: IS *ISForDofs; 100: IS *ISForDofsLocal; 101: IS NeumannBoundaries; 102: IS NeumannBoundariesLocal; 103: IS DirichletBoundaries; 104: IS DirichletBoundariesLocal; 105: PetscBool eliminate_dirdofs; 106: PetscBool switch_static; 107: PetscInt coarsening_ratio; 108: PetscInt coarse_adj_red; 109: PetscInt current_level; 110: PetscInt max_levels; 111: PetscInt coarse_eqs_per_proc; 112: IS coarse_subassembling; 113: PetscBool use_coarse_estimates; 114: PetscBool symmetric_primal; 115: PetscInt vertex_size; 117: /* no-net-flux */ 118: PetscBool compute_nonetflux; 119: Mat divudotp; 120: PetscBool divudotp_trans; 121: IS divudotp_vl2l; 123: /* nedelec */ 124: Mat discretegradient; 125: PetscInt nedorder; 126: PetscBool conforming; 127: PetscInt nedfield; 128: PetscBool nedglobal; 129: Mat nedcG; 130: IS nedclocal; 132: /* local disconnected subdomains */ 133: PetscBool detect_disconnected; 134: PetscInt n_local_subs; 135: IS *local_subs; 137: /* scaling */ 138: Vec work_scaling; 139: PetscBool use_deluxe_scaling; 140: PCBDDCDeluxeScaling deluxe_ctx; 141: PetscBool deluxe_zerorows; 142: PetscBool deluxe_singlemat; 144: /* schur complements on interface's subsets */ 145: PCBDDCSubSchurs sub_schurs; 146: PetscBool sub_schurs_rebuild; 147: PetscBool sub_schurs_exact_schur; 148: PetscInt sub_schurs_layers; 149: PetscBool sub_schurs_use_useradj; 150: PetscBool computed_rowadj; 152: /* adaptive selection of constraints */ 153: PetscBool adaptive_selection; 154: PetscBool adaptive_userdefined; 155: PetscReal adaptive_threshold; 156: PetscInt adaptive_nmin; 157: PetscInt adaptive_nmax; 158: PetscInt* adaptive_constraints_n; 159: PetscInt* adaptive_constraints_idxs; 160: PetscInt* adaptive_constraints_idxs_ptr; 161: PetscScalar* adaptive_constraints_data; 162: PetscInt* adaptive_constraints_data_ptr; 164: /* For verbose output of some bddc data structures */ 165: PetscInt dbg_flag; 166: PetscViewer dbg_viewer; 167: } PC_BDDC; 170: #endif /* __pcbddc_h */