1 #include <../src/ksp/pc/impls/bddc/bddc.h> 2 #include <../src/ksp/pc/impls/bddc/bddcprivate.h> 3 4 /* prototypes for deluxe functions */ 5 static PetscErrorCode PCBDDCScalingCreate_Deluxe(PC); 6 static PetscErrorCode PCBDDCScalingDestroy_Deluxe(PC); 7 static PetscErrorCode PCBDDCScalingSetUp_Deluxe(PC); 8 static PetscErrorCode PCBDDCScalingSetUp_Deluxe_Par(PC,PetscInt,PetscInt,PetscInt[],PetscInt[]); 9 static PetscErrorCode PCBDDCScalingSetUp_Deluxe_Seq(PC); 10 static PetscErrorCode PCBDDCScalingReset_Deluxe_Solvers(PCBDDCDeluxeScaling); 11 12 #undef __FUNCT__ 13 #define __FUNCT__ "PCBDDCScalingExtension_Basic" 14 static PetscErrorCode PCBDDCScalingExtension_Basic(PC pc, Vec local_interface_vector, Vec global_vector) 15 { 16 PC_IS* pcis = (PC_IS*)pc->data; 17 PC_BDDC* pcbddc = (PC_BDDC*)pc->data; 18 PetscErrorCode ierr; 19 20 PetscFunctionBegin; 21 /* Apply partition of unity */ 22 ierr = VecPointwiseMult(pcbddc->work_scaling,pcis->D,local_interface_vector);CHKERRQ(ierr); 23 ierr = VecSet(global_vector,0.0);CHKERRQ(ierr); 24 ierr = VecScatterBegin(pcis->global_to_B,pcbddc->work_scaling,global_vector,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 25 ierr = VecScatterEnd(pcis->global_to_B,pcbddc->work_scaling,global_vector,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 26 PetscFunctionReturn(0); 27 } 28 29 #undef __FUNCT__ 30 #define __FUNCT__ "PCBDDCScalingExtension_Deluxe" 31 static PetscErrorCode PCBDDCScalingExtension_Deluxe(PC pc, Vec x, Vec y) 32 { 33 PC_IS* pcis=(PC_IS*)pc->data; 34 PC_BDDC* pcbddc=(PC_BDDC*)pc->data; 35 PCBDDCDeluxeScaling deluxe_ctx = pcbddc->deluxe_ctx; 36 PCBDDCSubSchurs sub_schurs = pcbddc->sub_schurs; 37 PetscInt i; 38 PetscErrorCode ierr; 39 40 /* TODO CHECK STUFF RELATED WITH FAKE WORK */ 41 PetscFunctionBegin; 42 ierr = VecSet(pcbddc->work_scaling,0.0);CHKERRQ(ierr); /* needed by the fake work below */ 43 if (deluxe_ctx->n_simple) { 44 /* scale deluxe vertices using diagonal scaling */ 45 PetscScalar *array_x,*array_D,*array; 46 ierr = VecGetArray(x,&array_x);CHKERRQ(ierr); 47 ierr = VecGetArray(pcis->D,&array_D);CHKERRQ(ierr); 48 ierr = VecGetArray(pcbddc->work_scaling,&array);CHKERRQ(ierr); 49 for (i=0;i<deluxe_ctx->n_simple;i++) { 50 array[deluxe_ctx->idx_simple_B[i]] = array_x[deluxe_ctx->idx_simple_B[i]]*array_D[deluxe_ctx->idx_simple_B[i]]; 51 } 52 ierr = VecRestoreArray(pcbddc->work_scaling,&array);CHKERRQ(ierr); 53 ierr = VecRestoreArray(pcis->D,&array_D);CHKERRQ(ierr); 54 ierr = VecRestoreArray(x,&array_x);CHKERRQ(ierr); 55 } 56 /* sequential part : all problems and Schur applications collapsed into a single matrix vector multiplication and ksp solution */ 57 if (deluxe_ctx->seq_ksp) { 58 ierr = VecScatterBegin(deluxe_ctx->seq_scctx,x,deluxe_ctx->seq_work1,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 59 ierr = VecScatterEnd(deluxe_ctx->seq_scctx,x,deluxe_ctx->seq_work1,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 60 ierr = MatMultTranspose(sub_schurs->S_Ej_all,deluxe_ctx->seq_work1,deluxe_ctx->seq_work2);CHKERRQ(ierr); 61 ierr = KSPSolveTranspose(deluxe_ctx->seq_ksp,deluxe_ctx->seq_work2,deluxe_ctx->seq_work1);CHKERRQ(ierr); 62 /* fake work due to final ADD VALUES and vertices scaling needed? TODO: check it */ 63 ierr = VecScatterBegin(deluxe_ctx->seq_scctx,deluxe_ctx->seq_work1,pcbddc->work_scaling,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 64 ierr = VecScatterEnd(deluxe_ctx->seq_scctx,deluxe_ctx->seq_work1,pcbddc->work_scaling,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 65 } 66 /* parallel part */ 67 for (i=0;i<deluxe_ctx->par_colors;i++) { 68 if (deluxe_ctx->par_ksp[i]) { 69 PetscMPIInt color_rank; 70 PetscInt subidx = deluxe_ctx->par_col2sub[i]; 71 /* restrict on subset */ 72 ierr = VecScatterBegin(deluxe_ctx->par_scctx_s[i],x,deluxe_ctx->par_work1[i],INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 73 ierr = VecScatterEnd(deluxe_ctx->par_scctx_s[i],x,deluxe_ctx->par_work1[i],INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 74 /* S_Ej */ 75 ierr = MatMult(sub_schurs->S_Ej[subidx],deluxe_ctx->par_work1[i],deluxe_ctx->par_work2[i]);CHKERRQ(ierr); 76 /* (\sum_j S_Ej)^-1 */ 77 ierr = VecSet(deluxe_ctx->par_vec[i],0.0);CHKERRQ(ierr); 78 ierr = VecScatterBegin(deluxe_ctx->par_scctx_p[i],deluxe_ctx->par_work2[i],deluxe_ctx->par_vec[i],ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 79 ierr = VecScatterEnd(deluxe_ctx->par_scctx_p[i],deluxe_ctx->par_work2[i],deluxe_ctx->par_vec[i],ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 80 ierr = KSPSolve(deluxe_ctx->par_ksp[i],deluxe_ctx->par_vec[i],deluxe_ctx->par_vec[i]);CHKERRQ(ierr); 81 ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)deluxe_ctx->par_ksp[i]),&color_rank);CHKERRQ(ierr); 82 /* get back solution on subset */ 83 ierr = VecScatterBegin(deluxe_ctx->par_scctx_p[i],deluxe_ctx->par_vec[i],deluxe_ctx->par_work1[i],INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 84 ierr = VecScatterEnd(deluxe_ctx->par_scctx_p[i],deluxe_ctx->par_vec[i],deluxe_ctx->par_work1[i],INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 85 if (!color_rank) { /* only the master process in coloured comm copies the computed values */ 86 ierr = VecScatterBegin(deluxe_ctx->par_scctx_s[i],deluxe_ctx->par_work1[i],pcbddc->work_scaling,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 87 ierr = VecScatterEnd(deluxe_ctx->par_scctx_s[i],deluxe_ctx->par_work1[i],pcbddc->work_scaling,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 88 } 89 } 90 } 91 /* put local boundary part in global vector */ 92 ierr = VecSet(y,0.0);CHKERRQ(ierr); 93 ierr = VecScatterBegin(pcis->global_to_B,pcbddc->work_scaling,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 94 ierr = VecScatterEnd(pcis->global_to_B,pcbddc->work_scaling,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 95 PetscFunctionReturn(0); 96 } 97 98 #undef __FUNCT__ 99 #define __FUNCT__ "PCBDDCScalingExtension" 100 PetscErrorCode PCBDDCScalingExtension(PC pc, Vec local_interface_vector, Vec global_vector) 101 { 102 PC_BDDC *pcbddc=(PC_BDDC*)pc->data; 103 PetscErrorCode ierr; 104 105 PetscFunctionBegin; 106 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 107 PetscValidHeaderSpecific(local_interface_vector,VEC_CLASSID,2); 108 PetscValidHeaderSpecific(global_vector,VEC_CLASSID,3); 109 if (local_interface_vector == pcbddc->work_scaling) { 110 SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"Local vector cannot be pcbddc->work_scaling!\n"); 111 } 112 ierr = PetscTryMethod(pc,"PCBDDCScalingExtension_C",(PC,Vec,Vec),(pc,local_interface_vector,global_vector));CHKERRQ(ierr); 113 PetscFunctionReturn(0); 114 } 115 116 #undef __FUNCT__ 117 #define __FUNCT__ "PCBDDCScalingRestriction_Basic" 118 static PetscErrorCode PCBDDCScalingRestriction_Basic(PC pc, Vec global_vector, Vec local_interface_vector) 119 { 120 PetscErrorCode ierr; 121 PC_IS* pcis = (PC_IS*)pc->data; 122 123 PetscFunctionBegin; 124 ierr = VecScatterBegin(pcis->global_to_B,global_vector,local_interface_vector,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 125 ierr = VecScatterEnd(pcis->global_to_B,global_vector,local_interface_vector,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 126 /* Apply partition of unity */ 127 ierr = VecPointwiseMult(local_interface_vector,pcis->D,local_interface_vector);CHKERRQ(ierr); 128 PetscFunctionReturn(0); 129 } 130 131 #undef __FUNCT__ 132 #define __FUNCT__ "PCBDDCScalingRestriction_Deluxe" 133 static PetscErrorCode PCBDDCScalingRestriction_Deluxe(PC pc, Vec x, Vec y) 134 { 135 PC_IS* pcis=(PC_IS*)pc->data; 136 PC_BDDC* pcbddc=(PC_BDDC*)pc->data; 137 PCBDDCDeluxeScaling deluxe_ctx = pcbddc->deluxe_ctx; 138 PCBDDCSubSchurs sub_schurs = pcbddc->sub_schurs; 139 PetscInt i; 140 PetscErrorCode ierr; 141 142 PetscFunctionBegin; 143 /* get local boundary part of global vector */ 144 ierr = VecScatterBegin(pcis->global_to_B,x,y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 145 ierr = VecScatterEnd(pcis->global_to_B,x,y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 146 if (deluxe_ctx->n_simple) { 147 /* scale deluxe vertices using diagonal scaling */ 148 PetscScalar *array_y,*array_D; 149 ierr = VecGetArray(y,&array_y);CHKERRQ(ierr); 150 ierr = VecGetArray(pcis->D,&array_D);CHKERRQ(ierr); 151 for (i=0;i<deluxe_ctx->n_simple;i++) { 152 array_y[deluxe_ctx->idx_simple_B[i]] *= array_D[deluxe_ctx->idx_simple_B[i]]; 153 } 154 ierr = VecRestoreArray(pcis->D,&array_D);CHKERRQ(ierr); 155 ierr = VecRestoreArray(y,&array_y);CHKERRQ(ierr); 156 } 157 /* sequential part : all problems and Schur applications collapsed into a single matrix vector multiplication and ksp solution */ 158 if (deluxe_ctx->seq_ksp) { 159 ierr = VecScatterBegin(deluxe_ctx->seq_scctx,y,deluxe_ctx->seq_work1,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 160 ierr = VecScatterEnd(deluxe_ctx->seq_scctx,y,deluxe_ctx->seq_work1,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 161 ierr = KSPSolve(deluxe_ctx->seq_ksp,deluxe_ctx->seq_work1,deluxe_ctx->seq_work2);CHKERRQ(ierr); 162 ierr = MatMult(sub_schurs->S_Ej_all,deluxe_ctx->seq_work2,deluxe_ctx->seq_work1);CHKERRQ(ierr); 163 ierr = VecScatterBegin(deluxe_ctx->seq_scctx,deluxe_ctx->seq_work1,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 164 ierr = VecScatterEnd(deluxe_ctx->seq_scctx,deluxe_ctx->seq_work1,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 165 } 166 /* parallel part */ 167 for (i=0;i<deluxe_ctx->par_colors;i++) { 168 if (deluxe_ctx->par_ksp[i]) { 169 PetscInt subidx = deluxe_ctx->par_col2sub[i]; 170 /* restrict on subset */ 171 ierr = VecScatterBegin(deluxe_ctx->par_scctx_s[i],y,deluxe_ctx->par_work1[i],INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 172 ierr = VecScatterEnd(deluxe_ctx->par_scctx_s[i],y,deluxe_ctx->par_work1[i],INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 173 /* (\sum_j S_Ej)^-T */ 174 ierr = VecSet(deluxe_ctx->par_vec[i],0.0);CHKERRQ(ierr); 175 ierr = VecScatterBegin(deluxe_ctx->par_scctx_p[i],deluxe_ctx->par_work1[i],deluxe_ctx->par_vec[i],INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 176 ierr = VecScatterEnd(deluxe_ctx->par_scctx_p[i],deluxe_ctx->par_work1[i],deluxe_ctx->par_vec[i],INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 177 ierr = KSPSolveTranspose(deluxe_ctx->par_ksp[i],deluxe_ctx->par_vec[i],deluxe_ctx->par_vec[i]);CHKERRQ(ierr); 178 ierr = VecScatterBegin(deluxe_ctx->par_scctx_p[i],deluxe_ctx->par_vec[i],deluxe_ctx->par_work1[i],INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 179 ierr = VecScatterEnd(deluxe_ctx->par_scctx_p[i],deluxe_ctx->par_vec[i],deluxe_ctx->par_work1[i],INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 180 /* S_Ej^T */ 181 ierr = MatMultTranspose(sub_schurs->S_Ej[subidx],deluxe_ctx->par_work1[i],deluxe_ctx->par_work2[i]);CHKERRQ(ierr); 182 /* extend to boundary */ 183 ierr = VecScatterBegin(deluxe_ctx->par_scctx_s[i],deluxe_ctx->par_work2[i],y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 184 ierr = VecScatterEnd(deluxe_ctx->par_scctx_s[i],deluxe_ctx->par_work2[i],y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 185 } 186 } 187 PetscFunctionReturn(0); 188 } 189 190 #undef __FUNCT__ 191 #define __FUNCT__ "PCBDDCScalingRestriction" 192 PetscErrorCode PCBDDCScalingRestriction(PC pc, Vec global_vector, Vec local_interface_vector) 193 { 194 PC_BDDC *pcbddc=(PC_BDDC*)pc->data; 195 PetscErrorCode ierr; 196 197 PetscFunctionBegin; 198 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 199 PetscValidHeaderSpecific(global_vector,VEC_CLASSID,2); 200 PetscValidHeaderSpecific(local_interface_vector,VEC_CLASSID,3); 201 if (local_interface_vector == pcbddc->work_scaling) { 202 SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"Local vector should cannot be pcbddc->work_scaling!\n"); 203 } 204 ierr = PetscTryMethod(pc,"PCBDDCScalingRestriction_C",(PC,Vec,Vec),(pc,global_vector,local_interface_vector));CHKERRQ(ierr); 205 PetscFunctionReturn(0); 206 } 207 208 #undef __FUNCT__ 209 #define __FUNCT__ "PCBDDCScalingSetUp" 210 PetscErrorCode PCBDDCScalingSetUp(PC pc) 211 { 212 PC_IS* pcis=(PC_IS*)pc->data; 213 PC_BDDC* pcbddc=(PC_BDDC*)pc->data; 214 PetscErrorCode ierr; 215 216 PetscFunctionBegin; 217 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 218 /* create work vector for the operator */ 219 ierr = VecDestroy(&pcbddc->work_scaling);CHKERRQ(ierr); 220 ierr = VecDuplicate(pcis->vec1_B,&pcbddc->work_scaling);CHKERRQ(ierr); 221 /* always rebuild pcis->D */ 222 if (pcis->use_stiffness_scaling) { 223 ierr = MatGetDiagonal(pcbddc->local_mat,pcis->vec1_N);CHKERRQ(ierr); 224 ierr = VecScatterBegin(pcis->N_to_B,pcis->vec1_N,pcis->D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 225 ierr = VecScatterEnd(pcis->N_to_B,pcis->vec1_N,pcis->D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 226 } 227 ierr = VecCopy(pcis->D,pcis->vec1_B);CHKERRQ(ierr); 228 ierr = VecSet(pcis->vec1_global,0.0);CHKERRQ(ierr); 229 ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_B,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 230 ierr = VecScatterEnd(pcis->global_to_B,pcis->vec1_B,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 231 ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_global,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 232 ierr = VecScatterEnd(pcis->global_to_B,pcis->vec1_global,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 233 ierr = VecPointwiseDivide(pcis->D,pcis->D,pcis->vec1_B);CHKERRQ(ierr); 234 /* now setup */ 235 if (pcbddc->use_deluxe_scaling) { 236 if (!pcbddc->deluxe_ctx) { 237 ierr = PCBDDCScalingCreate_Deluxe(pc);CHKERRQ(ierr); 238 } 239 ierr = PCBDDCScalingSetUp_Deluxe(pc);CHKERRQ(ierr); 240 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCScalingRestriction_C",PCBDDCScalingRestriction_Deluxe);CHKERRQ(ierr); 241 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCScalingExtension_C",PCBDDCScalingExtension_Deluxe);CHKERRQ(ierr); 242 } else { 243 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCScalingRestriction_C",PCBDDCScalingRestriction_Basic);CHKERRQ(ierr); 244 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCScalingExtension_C",PCBDDCScalingExtension_Basic);CHKERRQ(ierr); 245 } 246 247 /* test */ 248 if (pcbddc->dbg_flag) { 249 Vec vec2_global; 250 PetscViewer viewer=pcbddc->dbg_viewer; 251 PetscReal error; 252 253 /* extension -> from local to parallel */ 254 ierr = VecSet(pcis->vec1_global,0.0);CHKERRQ(ierr); 255 ierr = VecSetRandom(pcis->vec1_B,NULL);CHKERRQ(ierr); 256 ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_B,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 257 ierr = VecScatterEnd(pcis->global_to_B,pcis->vec1_B,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 258 ierr = VecDuplicate(pcis->vec1_global,&vec2_global);CHKERRQ(ierr); 259 ierr = VecCopy(pcis->vec1_global,vec2_global);CHKERRQ(ierr); 260 261 ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_global,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 262 ierr = VecScatterEnd(pcis->global_to_B,pcis->vec1_global,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 263 ierr = PCBDDCScalingExtension(pc,pcis->vec1_B,pcis->vec1_global);CHKERRQ(ierr); 264 ierr = VecAXPY(pcis->vec1_global,-1.0,vec2_global);CHKERRQ(ierr); 265 ierr = VecNorm(pcis->vec1_global,NORM_INFINITY,&error);CHKERRQ(ierr); 266 ierr = PetscViewerASCIIPrintf(viewer,"Error scaling extension %1.14e\n",error);CHKERRQ(ierr); 267 if (error>1.e-8 && pcbddc->dbg_flag>1) { 268 ierr = VecView(pcis->vec1_global,viewer);CHKERRQ(ierr); 269 } 270 ierr = VecDestroy(&vec2_global);CHKERRQ(ierr); 271 272 /* restriction -> from parallel to local */ 273 ierr = VecSet(pcis->vec1_global,0.0);CHKERRQ(ierr); 274 ierr = VecSetRandom(pcis->vec1_B,NULL);CHKERRQ(ierr); 275 ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_B,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 276 ierr = VecScatterEnd(pcis->global_to_B,pcis->vec1_B,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 277 278 ierr = PCBDDCScalingRestriction(pc,pcis->vec1_global,pcis->vec1_B);CHKERRQ(ierr); 279 ierr = VecScale(pcis->vec1_B,-1.0);CHKERRQ(ierr); 280 ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_B,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 281 ierr = VecScatterEnd(pcis->global_to_B,pcis->vec1_B,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 282 ierr = VecNorm(pcis->vec1_global,NORM_INFINITY,&error);CHKERRQ(ierr); 283 ierr = PetscViewerASCIIPrintf(viewer,"Error scaling restriction %1.14e\n",error);CHKERRQ(ierr); 284 if (error>1.e-8 && pcbddc->dbg_flag>1) { 285 ierr = VecView(pcis->vec1_global,viewer);CHKERRQ(ierr); 286 } 287 } 288 PetscFunctionReturn(0); 289 } 290 291 #undef __FUNCT__ 292 #define __FUNCT__ "PCBDDCScalingDestroy" 293 PetscErrorCode PCBDDCScalingDestroy(PC pc) 294 { 295 PC_BDDC* pcbddc=(PC_BDDC*)pc->data; 296 PetscErrorCode ierr; 297 298 PetscFunctionBegin; 299 if (pcbddc->deluxe_ctx) { 300 ierr = PCBDDCScalingDestroy_Deluxe(pc);CHKERRQ(ierr); 301 } 302 ierr = VecDestroy(&pcbddc->work_scaling);CHKERRQ(ierr); 303 /* remove functions */ 304 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCScalingRestriction_C",NULL);CHKERRQ(ierr); 305 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCScalingExtension_C",NULL);CHKERRQ(ierr); 306 PetscFunctionReturn(0); 307 } 308 309 #undef __FUNCT__ 310 #define __FUNCT__ "PCBDDCScalingCreate_Deluxe" 311 static PetscErrorCode PCBDDCScalingCreate_Deluxe(PC pc) 312 { 313 PC_BDDC* pcbddc=(PC_BDDC*)pc->data; 314 PCBDDCDeluxeScaling deluxe_ctx; 315 PetscErrorCode ierr; 316 317 PetscFunctionBegin; 318 ierr = PetscNew(&deluxe_ctx);CHKERRQ(ierr); 319 pcbddc->deluxe_ctx = deluxe_ctx; 320 PetscFunctionReturn(0); 321 } 322 323 #undef __FUNCT__ 324 #define __FUNCT__ "PCBDDCScalingDestroy_Deluxe" 325 static PetscErrorCode PCBDDCScalingDestroy_Deluxe(PC pc) 326 { 327 PC_BDDC* pcbddc=(PC_BDDC*)pc->data; 328 PetscErrorCode ierr; 329 330 PetscFunctionBegin; 331 ierr = PCBDDCScalingReset_Deluxe_Solvers(pcbddc->deluxe_ctx);CHKERRQ(ierr); 332 ierr = PetscFree(pcbddc->deluxe_ctx);CHKERRQ(ierr); 333 PetscFunctionReturn(0); 334 } 335 336 #undef __FUNCT__ 337 #define __FUNCT__ "PCBDDCScalingReset_Deluxe_Solvers" 338 static PetscErrorCode PCBDDCScalingReset_Deluxe_Solvers(PCBDDCDeluxeScaling deluxe_ctx) 339 { 340 PetscErrorCode ierr; 341 342 PetscFunctionBegin; 343 ierr = PetscFree(deluxe_ctx->idx_simple_B);CHKERRQ(ierr); 344 deluxe_ctx->n_simple = 0; 345 if (deluxe_ctx->seq_ksp) { 346 ierr = VecScatterDestroy(&deluxe_ctx->seq_scctx);CHKERRQ(ierr); 347 ierr = VecDestroy(&deluxe_ctx->seq_work1);CHKERRQ(ierr); 348 ierr = VecDestroy(&deluxe_ctx->seq_work2);CHKERRQ(ierr); 349 ierr = KSPDestroy(&deluxe_ctx->seq_ksp);CHKERRQ(ierr); 350 } 351 if (deluxe_ctx->par_colors) { 352 PetscInt i; 353 for (i=0;i<deluxe_ctx->par_colors;i++) { 354 ierr = VecScatterDestroy(&deluxe_ctx->par_scctx_s[i]);CHKERRQ(ierr); 355 ierr = VecScatterDestroy(&deluxe_ctx->par_scctx_p[i]);CHKERRQ(ierr); 356 ierr = VecDestroy(&deluxe_ctx->par_vec[i]);CHKERRQ(ierr); 357 ierr = VecDestroy(&deluxe_ctx->par_work1[i]);CHKERRQ(ierr); 358 ierr = VecDestroy(&deluxe_ctx->par_work2[i]);CHKERRQ(ierr); 359 ierr = KSPDestroy(&deluxe_ctx->par_ksp[i]);CHKERRQ(ierr); 360 } 361 ierr = PetscFree7(deluxe_ctx->par_ksp, 362 deluxe_ctx->par_scctx_s, 363 deluxe_ctx->par_scctx_p, 364 deluxe_ctx->par_vec, 365 deluxe_ctx->par_work1, 366 deluxe_ctx->par_work2, 367 deluxe_ctx->par_col2sub);CHKERRQ(ierr); 368 } 369 deluxe_ctx->par_colors = 0; 370 PetscFunctionReturn(0); 371 } 372 373 #undef __FUNCT__ 374 #define __FUNCT__ "PCBDDCScalingSetUp_Deluxe" 375 static PetscErrorCode PCBDDCScalingSetUp_Deluxe(PC pc) 376 { 377 PC_IS *pcis=(PC_IS*)pc->data; 378 PC_BDDC *pcbddc=(PC_BDDC*)pc->data; 379 PCBDDCDeluxeScaling deluxe_ctx=pcbddc->deluxe_ctx; 380 PCBDDCSubSchurs sub_schurs=pcbddc->sub_schurs; 381 PetscBT bitmask; 382 PetscInt i; 383 const PetscInt *idxs; 384 PetscErrorCode ierr; 385 386 PetscFunctionBegin; 387 /* (TODO: reuse) throw away the solvers */ 388 ierr = PCBDDCScalingReset_Deluxe_Solvers(deluxe_ctx);CHKERRQ(ierr); 389 390 /* Compute data structures to solve parallel problems */ 391 ierr = PCBDDCScalingSetUp_Deluxe_Par(pc,sub_schurs->n_subs_par,sub_schurs->n_subs_par_g, 392 sub_schurs->auxglobal_parallel, 393 sub_schurs->index_parallel);CHKERRQ(ierr); 394 395 /* Compute data structures to solve sequential problems */ 396 ierr = PCBDDCScalingSetUp_Deluxe_Seq(pc);CHKERRQ(ierr); 397 398 /* diagonal scaling on interface dofs not contained in cc */ 399 ierr = PetscBTCreate(pcis->n,&bitmask);CHKERRQ(ierr); 400 ierr = ISGetIndices(pcis->is_I_local,&idxs);CHKERRQ(ierr); 401 for (i=0;i<pcis->n-pcis->n_B;i++) { 402 ierr = PetscBTSet(bitmask,idxs[i]);CHKERRQ(ierr); 403 } 404 ierr = ISRestoreIndices(pcis->is_I_local,&idxs);CHKERRQ(ierr); 405 406 for (i=0;i<sub_schurs->n_subs;i++) { 407 PetscInt subset_size,j; 408 409 ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr); 410 ierr = ISGetIndices(sub_schurs->is_subs[i],&idxs);CHKERRQ(ierr); 411 for (j=0;j<subset_size;j++) { 412 ierr = PetscBTSet(bitmask,idxs[j]);CHKERRQ(ierr); 413 } 414 ierr = ISRestoreIndices(sub_schurs->is_subs[i],&idxs);CHKERRQ(ierr); 415 } 416 417 deluxe_ctx->n_simple = 0; 418 for (i=0;i<pcis->n;i++) { 419 if (!PetscBTLookup(bitmask,i)) { 420 deluxe_ctx->n_simple++; 421 } 422 } 423 ierr = PetscMalloc1(deluxe_ctx->n_simple,&deluxe_ctx->idx_simple_B);CHKERRQ(ierr); 424 deluxe_ctx->n_simple = 0; 425 for (i=0;i<pcis->n;i++) { 426 if (!PetscBTLookup(bitmask,i)) { 427 deluxe_ctx->idx_simple_B[deluxe_ctx->n_simple++] = i; 428 } 429 } 430 ierr = ISGlobalToLocalMappingApply(pcis->BtoNmap,IS_GTOLM_DROP,deluxe_ctx->n_simple,deluxe_ctx->idx_simple_B,&i,deluxe_ctx->idx_simple_B);CHKERRQ(ierr); 431 if (i != deluxe_ctx->n_simple) { 432 SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Error when mapping simple scaling dofs! %d != %d",i,deluxe_ctx->n_simple); 433 } 434 ierr = PetscBTDestroy(&bitmask);CHKERRQ(ierr); 435 PetscFunctionReturn(0); 436 } 437 438 #undef __FUNCT__ 439 #define __FUNCT__ "PCBDDCScalingSetUp_Deluxe_Par" 440 static PetscErrorCode PCBDDCScalingSetUp_Deluxe_Par(PC pc, PetscInt n_local_parallel_problems,PetscInt n_parallel_problems,PetscInt global_parallel[],PetscInt index_parallel[]) 441 { 442 PC_BDDC *pcbddc=(PC_BDDC*)pc->data; 443 PCBDDCDeluxeScaling deluxe_ctx=pcbddc->deluxe_ctx; 444 /* coloring */ 445 Mat parallel_problems; 446 MatColoring coloring_obj; 447 ISColoring coloring_parallel_problems; 448 IS *par_is_colors,*is_colors; 449 /* working stuff */ 450 PetscInt i,j; 451 PetscErrorCode ierr; 452 453 PetscFunctionBegin; 454 if (!n_parallel_problems) { 455 PetscFunctionReturn(0); 456 } 457 /* Color parallel subproblems */ 458 ierr = MatCreate(PetscObjectComm((PetscObject)pc),¶llel_problems);CHKERRQ(ierr); 459 ierr = MatSetSizes(parallel_problems,PETSC_DECIDE,PETSC_DECIDE,n_parallel_problems,n_parallel_problems);CHKERRQ(ierr); 460 ierr = MatSetType(parallel_problems,MATAIJ);CHKERRQ(ierr); 461 ierr = MatSetUp(parallel_problems);CHKERRQ(ierr); 462 ierr = MatSetOption(parallel_problems,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_FALSE);CHKERRQ(ierr); 463 ierr = MatSetOption(parallel_problems,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_FALSE);CHKERRQ(ierr); 464 for (i=0;i<n_local_parallel_problems;i++) { 465 PetscInt row = global_parallel[i]; 466 for (j=0;j<n_local_parallel_problems;j++) { 467 PetscInt col = global_parallel[j]; 468 if (row != col) { 469 ierr = MatSetValue(parallel_problems,row,col,1.0,INSERT_VALUES);CHKERRQ(ierr); 470 } 471 } 472 } 473 ierr = MatAssemblyBegin(parallel_problems,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 474 ierr = MatAssemblyEnd(parallel_problems,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 475 if (pcbddc->dbg_flag > 1) { 476 ierr = PetscViewerFlush(pcbddc->dbg_viewer);CHKERRQ(ierr); 477 ierr = PetscViewerASCIIPrintf(pcbddc->dbg_viewer,"Adj matrix for deluxe parallel problems\n");CHKERRQ(ierr); 478 ierr = MatView(parallel_problems,pcbddc->dbg_viewer);CHKERRQ(ierr); 479 } 480 ierr = MatColoringCreate(parallel_problems,&coloring_obj);CHKERRQ(ierr); 481 ierr = MatColoringSetDistance(coloring_obj,1);CHKERRQ(ierr); 482 ierr = MatColoringSetType(coloring_obj,MATCOLORINGJP);CHKERRQ(ierr); 483 ierr = MatColoringApply(coloring_obj,&coloring_parallel_problems);CHKERRQ(ierr); 484 ierr = ISColoringGetIS(coloring_parallel_problems,&deluxe_ctx->par_colors,&par_is_colors);CHKERRQ(ierr); 485 if (pcbddc->dbg_flag) { 486 ierr = PetscViewerFlush(pcbddc->dbg_viewer);CHKERRQ(ierr); 487 ierr = PetscViewerASCIIPrintf(pcbddc->dbg_viewer,"Number of colors %d for parallel part of deluxe\n",deluxe_ctx->par_colors);CHKERRQ(ierr); 488 } 489 490 /* all procs should know the color distribution */ 491 ierr = PetscMalloc1(deluxe_ctx->par_colors,&is_colors);CHKERRQ(ierr); 492 for (i=0;i<deluxe_ctx->par_colors;i++) { 493 if (pcbddc->dbg_flag) { 494 ierr = PetscViewerASCIIPrintf(pcbddc->dbg_viewer,"Global problem indexes for color %d\n",i);CHKERRQ(ierr); 495 ierr = ISView(par_is_colors[i],pcbddc->dbg_viewer);CHKERRQ(ierr); 496 ierr = PetscViewerFlush(pcbddc->dbg_viewer);CHKERRQ(ierr); 497 } 498 ierr = ISAllGather(par_is_colors[i],&is_colors[i]);CHKERRQ(ierr); 499 } 500 501 /* free unneeded objects */ 502 ierr = ISColoringRestoreIS(coloring_parallel_problems,&par_is_colors);CHKERRQ(ierr); 503 ierr = ISColoringDestroy(&coloring_parallel_problems);CHKERRQ(ierr); 504 ierr = MatColoringDestroy(&coloring_obj);CHKERRQ(ierr); 505 ierr = MatDestroy(¶llel_problems);CHKERRQ(ierr); 506 507 /* allocate deluxe arrays for parallel problems */ 508 ierr = PetscCalloc7(deluxe_ctx->par_colors,&deluxe_ctx->par_ksp, 509 deluxe_ctx->par_colors,&deluxe_ctx->par_scctx_s, 510 deluxe_ctx->par_colors,&deluxe_ctx->par_scctx_p, 511 deluxe_ctx->par_colors,&deluxe_ctx->par_vec, 512 deluxe_ctx->par_colors,&deluxe_ctx->par_work1, 513 deluxe_ctx->par_colors,&deluxe_ctx->par_work2, 514 deluxe_ctx->par_colors,&deluxe_ctx->par_col2sub);CHKERRQ(ierr); 515 516 /* cycle on colors */ 517 for (i=0;i<deluxe_ctx->par_colors;i++) { 518 PetscSubcomm par_subcomm; 519 const PetscInt* idxs_subproblems; 520 PetscInt color_size; 521 PetscMPIInt rank,active_color; 522 523 /* get local index of i-th parallel colored problem */ 524 ierr = ISGetLocalSize(is_colors[i],&color_size);CHKERRQ(ierr); 525 ierr = ISGetIndices(is_colors[i],&idxs_subproblems);CHKERRQ(ierr); 526 /* split comm for computing parallel problems for this color */ 527 /* Processes not partecipating at this stage will have color = color_size */ 528 /* because PetscCommDuplicate does not handle MPI_COMM_NULL */ 529 active_color = color_size; 530 deluxe_ctx->par_col2sub[i] = -1; 531 for (j=0;j<n_local_parallel_problems;j++) { 532 PetscInt local_idx; 533 ierr = PetscFindInt(global_parallel[j],color_size,idxs_subproblems,&local_idx);CHKERRQ(ierr); 534 if (local_idx > -1) { 535 ierr = PetscMPIIntCast(local_idx,&active_color);CHKERRQ(ierr); 536 deluxe_ctx->par_col2sub[i] = index_parallel[j]; 537 break; 538 } 539 } 540 ierr = ISRestoreIndices(is_colors[i],&idxs_subproblems);CHKERRQ(ierr); 541 ierr = PetscSubcommCreate(PetscObjectComm((PetscObject)pc),&par_subcomm);CHKERRQ(ierr); 542 ierr = PetscSubcommSetNumber(par_subcomm,color_size+1);CHKERRQ(ierr); 543 ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)pc),&rank);CHKERRQ(ierr); 544 ierr = PetscSubcommSetTypeGeneral(par_subcomm,active_color,rank);CHKERRQ(ierr); 545 /* print debug info */ 546 if (pcbddc->dbg_flag) { 547 PetscMPIInt crank,csize; 548 ierr = MPI_Comm_rank(par_subcomm->comm,&crank);CHKERRQ(ierr); 549 ierr = MPI_Comm_size(par_subcomm->comm,&csize);CHKERRQ(ierr); 550 ierr = PetscViewerASCIIPrintf(pcbddc->dbg_viewer,"Color %d: size %d, details follows.\n",i,color_size);CHKERRQ(ierr); 551 ierr = PetscViewerFlush(pcbddc->dbg_viewer);CHKERRQ(ierr); 552 ierr = PetscViewerASCIISynchronizedAllow(pcbddc->dbg_viewer,PETSC_TRUE);CHKERRQ(ierr); 553 ierr = PetscViewerASCIISynchronizedPrintf(pcbddc->dbg_viewer," Subdomain %d: color in subcomm %d (rank %d out of %d) (lidx %d)\n",PetscGlobalRank,par_subcomm->color,crank,csize,deluxe_ctx->par_col2sub[i]);CHKERRQ(ierr); 554 ierr = PetscViewerFlush(pcbddc->dbg_viewer);CHKERRQ(ierr); 555 } 556 557 if (deluxe_ctx->par_col2sub[i] >= 0) { 558 PC pctemp; 559 PC_IS *pcis=(PC_IS*)pc->data; 560 Mat color_mat,color_mat_is,temp_mat; 561 ISLocalToGlobalMapping WtoNmap,l2gmap_subset; 562 IS is_local_numbering,isB_local,isW_local,isW; 563 PCBDDCSubSchurs sub_schurs = pcbddc->sub_schurs; 564 PetscInt subidx,n_local_dofs,n_global_dofs; 565 PetscInt *global_numbering,*local_numbering; 566 char ksp_prefix[256]; 567 size_t len; 568 569 /* Local index for schur complement on subset */ 570 subidx = deluxe_ctx->par_col2sub[i]; 571 572 /* Parallel numbering for dofs in colored subset */ 573 ierr = ISSum(sub_schurs->is_I_layer,sub_schurs->is_subs[subidx],&is_local_numbering);CHKERRQ(ierr); 574 ierr = ISGetLocalSize(is_local_numbering,&n_local_dofs);CHKERRQ(ierr); 575 ierr = ISGetIndices(is_local_numbering,(const PetscInt **)&local_numbering);CHKERRQ(ierr); 576 ierr = PCBDDCSubsetNumbering(par_subcomm->comm,pcbddc->mat_graph->l2gmap,n_local_dofs,local_numbering,PETSC_NULL,&n_global_dofs,&global_numbering);CHKERRQ(ierr); 577 ierr = ISRestoreIndices(is_local_numbering,(const PetscInt **)&local_numbering);CHKERRQ(ierr); 578 579 /* L2Gmap from relevant dofs to local dofs */ 580 ierr = ISLocalToGlobalMappingCreateIS(is_local_numbering,&WtoNmap);CHKERRQ(ierr); 581 582 /* L2Gmap from local to global dofs */ 583 ierr = ISLocalToGlobalMappingCreate(par_subcomm->comm,1,n_local_dofs,global_numbering,PETSC_COPY_VALUES,&l2gmap_subset);CHKERRQ(ierr); 584 585 /* compute parallel matrix (extended dirichlet problem on subset) */ 586 ierr = MatCreateIS(par_subcomm->comm,1,PETSC_DECIDE,PETSC_DECIDE,n_global_dofs,n_global_dofs,l2gmap_subset,&color_mat_is);CHKERRQ(ierr); 587 ierr = MatGetSubMatrix(pcbddc->local_mat,is_local_numbering,is_local_numbering,MAT_INITIAL_MATRIX,&temp_mat);CHKERRQ(ierr); 588 ierr = MatISSetLocalMat(color_mat_is,temp_mat);CHKERRQ(ierr); 589 ierr = MatDestroy(&temp_mat);CHKERRQ(ierr); 590 ierr = MatISGetMPIXAIJ(color_mat_is,MAT_INITIAL_MATRIX,&color_mat);CHKERRQ(ierr); 591 ierr = MatDestroy(&color_mat_is);CHKERRQ(ierr); 592 593 /* work vector for (parallel) extended dirichlet problem */ 594 ierr = MatCreateVecs(color_mat,&deluxe_ctx->par_vec[i],NULL);CHKERRQ(ierr); 595 596 /* compute scatters */ 597 /* deluxe_ctx->par_scctx_p[i] extension from local subset to extended dirichlet problem 598 deluxe_ctx->par_scctx_s[i] restriction from local boundary to subset -> simple copy of selected values */ 599 ierr = ISGlobalToLocalMappingApplyIS(pcis->BtoNmap,IS_GTOLM_DROP,sub_schurs->is_subs[subidx],&isB_local);CHKERRQ(ierr); 600 ierr = MatCreateVecs(sub_schurs->S_Ej[subidx],&deluxe_ctx->par_work1[i],&deluxe_ctx->par_work2[i]);CHKERRQ(ierr); 601 ierr = VecScatterCreate(pcbddc->work_scaling,isB_local,deluxe_ctx->par_work1[i],NULL,&deluxe_ctx->par_scctx_s[i]);CHKERRQ(ierr); 602 ierr = ISGlobalToLocalMappingApplyIS(WtoNmap,IS_GTOLM_DROP,sub_schurs->is_subs[subidx],&isW_local);CHKERRQ(ierr); 603 ierr = ISLocalToGlobalMappingApplyIS(l2gmap_subset,isW_local,&isW);CHKERRQ(ierr); 604 ierr = VecScatterCreate(deluxe_ctx->par_work1[i],NULL,deluxe_ctx->par_vec[i],isW,&deluxe_ctx->par_scctx_p[i]);CHKERRQ(ierr); 605 606 /* free objects no longer neeeded */ 607 ierr = ISDestroy(&isW);CHKERRQ(ierr); 608 ierr = ISDestroy(&isW_local);CHKERRQ(ierr); 609 ierr = ISDestroy(&isB_local);CHKERRQ(ierr); 610 ierr = ISLocalToGlobalMappingDestroy(&WtoNmap);CHKERRQ(ierr); 611 ierr = ISLocalToGlobalMappingDestroy(&l2gmap_subset);CHKERRQ(ierr); 612 ierr = ISDestroy(&is_local_numbering);CHKERRQ(ierr); 613 ierr = PetscFree(global_numbering);CHKERRQ(ierr); 614 615 /* KSP for extended dirichlet problem */ 616 ierr = KSPCreate(par_subcomm->comm,&deluxe_ctx->par_ksp[i]);CHKERRQ(ierr); 617 ierr = KSPSetOperators(deluxe_ctx->par_ksp[i],color_mat,color_mat);CHKERRQ(ierr); 618 ierr = KSPSetTolerances(deluxe_ctx->par_ksp[i],1.e-12,1.e-12,1.e10,10000);CHKERRQ(ierr); 619 ierr = KSPSetType(deluxe_ctx->par_ksp[i],KSPPREONLY);CHKERRQ(ierr); 620 ierr = KSPGetPC(deluxe_ctx->par_ksp[i],&pctemp);CHKERRQ(ierr); 621 ierr = PCSetType(pctemp,PCREDUNDANT);CHKERRQ(ierr); 622 ierr = PetscStrlen(((PetscObject)(pcbddc->ksp_D))->prefix,&len);CHKERRQ(ierr); 623 len -= 10; /* remove "dirichlet_" */ 624 ierr = PetscStrncpy(ksp_prefix,((PetscObject)(pcbddc->ksp_D))->prefix,len+1);CHKERRQ(ierr); /* PetscStrncpy puts a terminating char at the end */ 625 ierr = PetscStrcat(ksp_prefix,"deluxe_par_");CHKERRQ(ierr); 626 ierr = KSPSetOptionsPrefix(deluxe_ctx->par_ksp[i],ksp_prefix);CHKERRQ(ierr); 627 ierr = KSPSetFromOptions(deluxe_ctx->par_ksp[i]);CHKERRQ(ierr); 628 ierr = KSPSetUp(deluxe_ctx->par_ksp[i]);CHKERRQ(ierr); 629 ierr = MatDestroy(&color_mat);CHKERRQ(ierr); 630 } 631 ierr = PetscSubcommDestroy(&par_subcomm);CHKERRQ(ierr); 632 } 633 for (i=0;i<deluxe_ctx->par_colors;i++) { 634 ierr = ISDestroy(&is_colors[i]);CHKERRQ(ierr); 635 } 636 ierr = PetscFree(is_colors);CHKERRQ(ierr); 637 638 if (pcbddc->dbg_flag) { 639 Vec test_vec; 640 PetscReal error; 641 PCBDDCSubSchurs sub_schurs = pcbddc->sub_schurs; 642 /* test partition of unity of coloured schur complements */ 643 for (i=0;i<deluxe_ctx->par_colors;i++) { 644 PetscInt subidx = deluxe_ctx->par_col2sub[i]; 645 PetscBool error_found = PETSC_FALSE; 646 ierr = PetscViewerASCIISynchronizedAllow(pcbddc->dbg_viewer,PETSC_TRUE);CHKERRQ(ierr); 647 648 if (deluxe_ctx->par_ksp[i]) { 649 /* create random test vec being zero on internal nodes of the extende dirichlet problem */ 650 ierr = VecDuplicate(deluxe_ctx->par_vec[i],&test_vec);CHKERRQ(ierr); 651 ierr = VecSetRandom(deluxe_ctx->par_work1[i],PETSC_NULL);CHKERRQ(ierr); 652 ierr = VecSet(test_vec,0.0);CHKERRQ(ierr); 653 ierr = VecScatterBegin(deluxe_ctx->par_scctx_p[i],deluxe_ctx->par_work1[i],test_vec,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 654 ierr = VecScatterEnd(deluxe_ctx->par_scctx_p[i],deluxe_ctx->par_work1[i],test_vec,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 655 /* w_j */ 656 ierr = VecScatterBegin(deluxe_ctx->par_scctx_p[i],test_vec,deluxe_ctx->par_work1[i],INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 657 ierr = VecScatterEnd(deluxe_ctx->par_scctx_p[i],test_vec,deluxe_ctx->par_work1[i],INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 658 /* S_j*w_j */ 659 ierr = MatMult(sub_schurs->S_Ej[subidx],deluxe_ctx->par_work1[i],deluxe_ctx->par_work2[i]);CHKERRQ(ierr); 660 /* \sum_j S_j*w_j */ 661 ierr = VecSet(deluxe_ctx->par_vec[i],0.0);CHKERRQ(ierr); 662 ierr = VecScatterBegin(deluxe_ctx->par_scctx_p[i],deluxe_ctx->par_work2[i],deluxe_ctx->par_vec[i],ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 663 ierr = VecScatterEnd(deluxe_ctx->par_scctx_p[i],deluxe_ctx->par_work2[i],deluxe_ctx->par_vec[i],ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 664 /* (\sum_j S_j)^(-1)(\sum_j S_j*w_j) */ 665 ierr = KSPSolve(deluxe_ctx->par_ksp[i],deluxe_ctx->par_vec[i],deluxe_ctx->par_vec[i]);CHKERRQ(ierr); 666 ierr = VecScatterBegin(deluxe_ctx->par_scctx_p[i],deluxe_ctx->par_vec[i],deluxe_ctx->par_work1[i],INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 667 ierr = VecScatterEnd(deluxe_ctx->par_scctx_p[i],deluxe_ctx->par_vec[i],deluxe_ctx->par_work1[i],INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 668 ierr = VecSet(deluxe_ctx->par_vec[i],0.0);CHKERRQ(ierr); 669 ierr = VecScatterBegin(deluxe_ctx->par_scctx_p[i],deluxe_ctx->par_work1[i],deluxe_ctx->par_vec[i],INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 670 ierr = VecScatterEnd(deluxe_ctx->par_scctx_p[i],deluxe_ctx->par_work1[i],deluxe_ctx->par_vec[i],INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 671 /* test partition of unity */ 672 ierr = VecAXPY(test_vec,-1.0,deluxe_ctx->par_vec[i]);CHKERRQ(ierr); 673 ierr = VecNorm(test_vec,NORM_INFINITY,&error);CHKERRQ(ierr); 674 if (PetscAbsReal(error) > 1.e-2) { 675 /* ierr = VecView(test_vec,0);CHKERRQ(ierr); */ 676 error_found = PETSC_TRUE; 677 } 678 ierr = VecDestroy(&test_vec);CHKERRQ(ierr); 679 } 680 if (error_found) { 681 ierr = PetscViewerASCIISynchronizedPrintf(pcbddc->dbg_viewer,"Error testing local schur for color %d and subdomain %d\n",i,PetscGlobalRank);CHKERRQ(ierr); 682 } 683 ierr = PetscViewerFlush(pcbddc->dbg_viewer);CHKERRQ(ierr); 684 } 685 } 686 PetscFunctionReturn(0); 687 } 688 689 #undef __FUNCT__ 690 #define __FUNCT__ "PCBDDCScalingSetUp_Deluxe_Seq" 691 static PetscErrorCode PCBDDCScalingSetUp_Deluxe_Seq(PC pc) 692 { 693 PC_BDDC *pcbddc=(PC_BDDC*)pc->data; 694 PCBDDCDeluxeScaling deluxe_ctx=pcbddc->deluxe_ctx; 695 PCBDDCSubSchurs sub_schurs = pcbddc->sub_schurs; 696 PC pc_temp; 697 MatSolverPackage solver=NULL; 698 char ksp_prefix[256]; 699 size_t len; 700 PetscInt local_size; 701 PetscErrorCode ierr; 702 703 PetscFunctionBegin; 704 if (!sub_schurs->n_subs_seq_g) { 705 PetscFunctionReturn(0); 706 } 707 708 /* Create work vectors for sequential part of deluxe */ 709 ierr = MatCreateVecs(sub_schurs->sum_S_Ej_all,&deluxe_ctx->seq_work1,&deluxe_ctx->seq_work2);CHKERRQ(ierr); 710 711 /* Compute deluxe sequential scatter */ 712 ierr = VecScatterCreate(pcbddc->work_scaling,sub_schurs->is_Ej_all,deluxe_ctx->seq_work1,NULL,&deluxe_ctx->seq_scctx);CHKERRQ(ierr); 713 714 /* Create KSP object for sequential part of deluxe scaling */ 715 ierr = KSPCreate(PETSC_COMM_SELF,&deluxe_ctx->seq_ksp);CHKERRQ(ierr); 716 ierr = KSPSetOperators(deluxe_ctx->seq_ksp,sub_schurs->sum_S_Ej_all,sub_schurs->sum_S_Ej_all);CHKERRQ(ierr); 717 ierr = KSPSetType(deluxe_ctx->seq_ksp,KSPPREONLY);CHKERRQ(ierr); 718 ierr = KSPGetPC(deluxe_ctx->seq_ksp,&pc_temp);CHKERRQ(ierr); 719 ierr = PCSetType(pc_temp,PCLU);CHKERRQ(ierr); 720 ierr = KSPGetPC(pcbddc->ksp_D,&pc_temp);CHKERRQ(ierr); 721 ierr = PCFactorGetMatSolverPackage(pc_temp,(const MatSolverPackage*)&solver);CHKERRQ(ierr); 722 ierr = MatGetSize(sub_schurs->sum_S_Ej_all,&local_size,NULL);CHKERRQ(ierr); 723 if (solver && local_size) { /* if local_size is null, some external packages will report errors */ 724 PC new_pc; 725 PCType type; 726 ierr = PCGetType(pc_temp,&type);CHKERRQ(ierr); 727 ierr = KSPGetPC(deluxe_ctx->seq_ksp,&new_pc);CHKERRQ(ierr); 728 ierr = PCSetType(new_pc,type);CHKERRQ(ierr); 729 ierr = PCFactorSetMatSolverPackage(new_pc,solver);CHKERRQ(ierr); 730 } 731 ierr = PetscStrlen(((PetscObject)(pcbddc->ksp_D))->prefix,&len);CHKERRQ(ierr); 732 len -= 10; /* remove "dirichlet_" */ 733 ierr = PetscStrncpy(ksp_prefix,((PetscObject)(pcbddc->ksp_D))->prefix,len+1);CHKERRQ(ierr); 734 ierr = PetscStrcat(ksp_prefix,"deluxe_seq_");CHKERRQ(ierr); 735 ierr = KSPSetOptionsPrefix(deluxe_ctx->seq_ksp,ksp_prefix);CHKERRQ(ierr); 736 if (local_size) { 737 ierr = KSPSetFromOptions(deluxe_ctx->seq_ksp);CHKERRQ(ierr); 738 } 739 ierr = KSPSetUp(deluxe_ctx->seq_ksp);CHKERRQ(ierr); 740 PetscFunctionReturn(0); 741 } 742