1 #include <../src/ksp/pc/impls/bddc/bddc.h> 2 #include <../src/ksp/pc/impls/bddc/bddcprivate.h> 3 #include <petscblaslapack.h> 4 #include <../src/mat/impls/dense/seq/dense.h> 5 6 /* prototypes for deluxe functions */ 7 static PetscErrorCode PCBDDCScalingCreate_Deluxe(PC); 8 static PetscErrorCode PCBDDCScalingDestroy_Deluxe(PC); 9 static PetscErrorCode PCBDDCScalingSetUp_Deluxe(PC); 10 static PetscErrorCode PCBDDCScalingSetUp_Deluxe_Private(PC); 11 static PetscErrorCode PCBDDCScalingReset_Deluxe_Solvers(PCBDDCDeluxeScaling); 12 13 static PetscErrorCode PCBDDCMatTransposeMatSolve_SeqDense(Mat A,Mat B,Mat X) 14 { 15 Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 16 PetscErrorCode ierr; 17 PetscScalar *b,*x; 18 PetscInt n; 19 PetscBLASInt nrhs,info,m; 20 PetscBool flg; 21 22 PetscFunctionBegin; 23 ierr = PetscBLASIntCast(A->rmap->n,&m);CHKERRQ(ierr); 24 ierr = PetscObjectTypeCompareAny((PetscObject)B,&flg,MATSEQDENSE,MATMPIDENSE,NULL);CHKERRQ(ierr); 25 if (!flg) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONG,"Matrix B must be MATDENSE matrix"); 26 ierr = PetscObjectTypeCompareAny((PetscObject)X,&flg,MATSEQDENSE,MATMPIDENSE,NULL);CHKERRQ(ierr); 27 if (!flg) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONG,"Matrix X must be MATDENSE matrix"); 28 29 ierr = MatGetSize(B,NULL,&n);CHKERRQ(ierr); 30 ierr = PetscBLASIntCast(n,&nrhs);CHKERRQ(ierr); 31 ierr = MatDenseGetArray(B,&b);CHKERRQ(ierr); 32 ierr = MatDenseGetArray(X,&x);CHKERRQ(ierr); 33 34 ierr = PetscMemcpy(x,b,m*nrhs*sizeof(PetscScalar));CHKERRQ(ierr); 35 36 if (A->factortype == MAT_FACTOR_LU) { 37 #if defined(PETSC_MISSING_LAPACK_GETRS) 38 SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"GETRS - Lapack routine is unavailable."); 39 #else 40 PetscStackCallBLAS("LAPACKgetrs",LAPACKgetrs_("T",&m,&nrhs,mat->v,&mat->lda,mat->pivots,x,&m,&info)); 41 if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"GETRS - Bad solve"); 42 #endif 43 } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Only LU factor supported"); 44 45 ierr = MatDenseRestoreArray(B,&b);CHKERRQ(ierr); 46 ierr = MatDenseRestoreArray(X,&x);CHKERRQ(ierr); 47 ierr = PetscLogFlops(nrhs*(2.0*m*m - m));CHKERRQ(ierr); 48 PetscFunctionReturn(0); 49 } 50 51 static PetscErrorCode PCBDDCScalingExtension_Basic(PC pc, Vec local_interface_vector, Vec global_vector) 52 { 53 PC_IS* pcis = (PC_IS*)pc->data; 54 PC_BDDC* pcbddc = (PC_BDDC*)pc->data; 55 PetscErrorCode ierr; 56 57 PetscFunctionBegin; 58 /* Apply partition of unity */ 59 ierr = VecPointwiseMult(pcbddc->work_scaling,pcis->D,local_interface_vector);CHKERRQ(ierr); 60 ierr = VecSet(global_vector,0.0);CHKERRQ(ierr); 61 ierr = VecScatterBegin(pcis->global_to_B,pcbddc->work_scaling,global_vector,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 62 ierr = VecScatterEnd(pcis->global_to_B,pcbddc->work_scaling,global_vector,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 63 PetscFunctionReturn(0); 64 } 65 66 static PetscErrorCode PCBDDCScalingExtension_Deluxe(PC pc, Vec x, Vec y) 67 { 68 PC_IS* pcis=(PC_IS*)pc->data; 69 PC_BDDC* pcbddc=(PC_BDDC*)pc->data; 70 PCBDDCDeluxeScaling deluxe_ctx = pcbddc->deluxe_ctx; 71 PetscErrorCode ierr; 72 73 PetscFunctionBegin; 74 ierr = VecSet(pcbddc->work_scaling,0.0);CHKERRQ(ierr); 75 ierr = VecSet(y,0.0);CHKERRQ(ierr); 76 if (deluxe_ctx->n_simple) { /* scale deluxe vertices using diagonal scaling */ 77 PetscInt i; 78 const PetscScalar *array_x,*array_D; 79 PetscScalar *array; 80 ierr = VecGetArrayRead(x,&array_x);CHKERRQ(ierr); 81 ierr = VecGetArrayRead(pcis->D,&array_D);CHKERRQ(ierr); 82 ierr = VecGetArray(pcbddc->work_scaling,&array);CHKERRQ(ierr); 83 for (i=0;i<deluxe_ctx->n_simple;i++) { 84 array[deluxe_ctx->idx_simple_B[i]] = array_x[deluxe_ctx->idx_simple_B[i]]*array_D[deluxe_ctx->idx_simple_B[i]]; 85 } 86 ierr = VecRestoreArray(pcbddc->work_scaling,&array);CHKERRQ(ierr); 87 ierr = VecRestoreArrayRead(pcis->D,&array_D);CHKERRQ(ierr); 88 ierr = VecRestoreArrayRead(x,&array_x);CHKERRQ(ierr); 89 } 90 /* sequential part : all problems and Schur applications collapsed into a single matrix vector multiplication or a matvec and a solve */ 91 if (deluxe_ctx->seq_mat) { 92 PetscInt i; 93 for (i=0;i<deluxe_ctx->seq_n;i++) { 94 if (deluxe_ctx->change) { 95 ierr = VecScatterBegin(deluxe_ctx->seq_scctx[i],x,deluxe_ctx->seq_work2[i],INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 96 ierr = VecScatterEnd(deluxe_ctx->seq_scctx[i],x,deluxe_ctx->seq_work2[i],INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 97 if (deluxe_ctx->change_with_qr) { 98 Mat change; 99 100 ierr = KSPGetOperators(deluxe_ctx->change[i],&change,NULL);CHKERRQ(ierr); 101 ierr = MatMultTranspose(change,deluxe_ctx->seq_work2[i],deluxe_ctx->seq_work1[i]);CHKERRQ(ierr); 102 } else { 103 ierr = KSPSolve(deluxe_ctx->change[i],deluxe_ctx->seq_work2[i],deluxe_ctx->seq_work1[i]);CHKERRQ(ierr); 104 } 105 } else { 106 ierr = VecScatterBegin(deluxe_ctx->seq_scctx[i],x,deluxe_ctx->seq_work1[i],INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 107 ierr = VecScatterEnd(deluxe_ctx->seq_scctx[i],x,deluxe_ctx->seq_work1[i],INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 108 } 109 ierr = MatMultTranspose(deluxe_ctx->seq_mat[i],deluxe_ctx->seq_work1[i],deluxe_ctx->seq_work2[i]);CHKERRQ(ierr); 110 if (deluxe_ctx->seq_mat_inv_sum[i]) { 111 PetscScalar *x; 112 113 ierr = VecGetArray(deluxe_ctx->seq_work2[i],&x);CHKERRQ(ierr); 114 ierr = VecPlaceArray(deluxe_ctx->seq_work1[i],x);CHKERRQ(ierr); 115 ierr = VecRestoreArray(deluxe_ctx->seq_work2[i],&x);CHKERRQ(ierr); 116 ierr = MatSolveTranspose(deluxe_ctx->seq_mat_inv_sum[i],deluxe_ctx->seq_work1[i],deluxe_ctx->seq_work2[i]);CHKERRQ(ierr); 117 ierr = VecResetArray(deluxe_ctx->seq_work1[i]);CHKERRQ(ierr); 118 } 119 if (deluxe_ctx->change) { 120 Mat change; 121 122 ierr = KSPGetOperators(deluxe_ctx->change[i],&change,NULL);CHKERRQ(ierr); 123 ierr = MatMult(change,deluxe_ctx->seq_work2[i],deluxe_ctx->seq_work1[i]);CHKERRQ(ierr); 124 ierr = VecScatterBegin(deluxe_ctx->seq_scctx[i],deluxe_ctx->seq_work1[i],pcbddc->work_scaling,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 125 ierr = VecScatterEnd(deluxe_ctx->seq_scctx[i],deluxe_ctx->seq_work1[i],pcbddc->work_scaling,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 126 } else { 127 ierr = VecScatterBegin(deluxe_ctx->seq_scctx[i],deluxe_ctx->seq_work2[i],pcbddc->work_scaling,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 128 ierr = VecScatterEnd(deluxe_ctx->seq_scctx[i],deluxe_ctx->seq_work2[i],pcbddc->work_scaling,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 129 } 130 } 131 } 132 /* put local boundary part in global vector */ 133 ierr = VecScatterBegin(pcis->global_to_B,pcbddc->work_scaling,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 134 ierr = VecScatterEnd(pcis->global_to_B,pcbddc->work_scaling,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 135 PetscFunctionReturn(0); 136 } 137 138 PetscErrorCode PCBDDCScalingExtension(PC pc, Vec local_interface_vector, Vec global_vector) 139 { 140 PC_BDDC *pcbddc=(PC_BDDC*)pc->data; 141 PetscErrorCode ierr; 142 143 PetscFunctionBegin; 144 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 145 PetscValidHeaderSpecific(local_interface_vector,VEC_CLASSID,2); 146 PetscValidHeaderSpecific(global_vector,VEC_CLASSID,3); 147 if (local_interface_vector == pcbddc->work_scaling) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Local vector cannot be pcbddc->work_scaling!\n"); 148 ierr = PetscUseMethod(pc,"PCBDDCScalingExtension_C",(PC,Vec,Vec),(pc,local_interface_vector,global_vector));CHKERRQ(ierr); 149 PetscFunctionReturn(0); 150 } 151 152 static PetscErrorCode PCBDDCScalingRestriction_Basic(PC pc, Vec global_vector, Vec local_interface_vector) 153 { 154 PetscErrorCode ierr; 155 PC_IS *pcis = (PC_IS*)pc->data; 156 157 PetscFunctionBegin; 158 ierr = VecScatterBegin(pcis->global_to_B,global_vector,local_interface_vector,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 159 ierr = VecScatterEnd(pcis->global_to_B,global_vector,local_interface_vector,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 160 /* Apply partition of unity */ 161 ierr = VecPointwiseMult(local_interface_vector,pcis->D,local_interface_vector);CHKERRQ(ierr); 162 PetscFunctionReturn(0); 163 } 164 165 static PetscErrorCode PCBDDCScalingRestriction_Deluxe(PC pc, Vec x, Vec y) 166 { 167 PC_IS* pcis=(PC_IS*)pc->data; 168 PC_BDDC* pcbddc=(PC_BDDC*)pc->data; 169 PCBDDCDeluxeScaling deluxe_ctx = pcbddc->deluxe_ctx; 170 PetscErrorCode ierr; 171 172 PetscFunctionBegin; 173 /* get local boundary part of global vector */ 174 ierr = VecScatterBegin(pcis->global_to_B,x,y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 175 ierr = VecScatterEnd(pcis->global_to_B,x,y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 176 if (deluxe_ctx->n_simple) { /* scale deluxe vertices using diagonal scaling */ 177 PetscInt i; 178 PetscScalar *array_y; 179 const PetscScalar *array_D; 180 ierr = VecGetArray(y,&array_y);CHKERRQ(ierr); 181 ierr = VecGetArrayRead(pcis->D,&array_D);CHKERRQ(ierr); 182 for (i=0;i<deluxe_ctx->n_simple;i++) { 183 array_y[deluxe_ctx->idx_simple_B[i]] *= array_D[deluxe_ctx->idx_simple_B[i]]; 184 } 185 ierr = VecRestoreArrayRead(pcis->D,&array_D);CHKERRQ(ierr); 186 ierr = VecRestoreArray(y,&array_y);CHKERRQ(ierr); 187 } 188 /* sequential part : all problems and Schur applications collapsed into a single matrix vector multiplication or a matvec and a solve */ 189 if (deluxe_ctx->seq_mat) { 190 PetscInt i; 191 for (i=0;i<deluxe_ctx->seq_n;i++) { 192 if (deluxe_ctx->change) { 193 Mat change; 194 195 ierr = VecScatterBegin(deluxe_ctx->seq_scctx[i],y,deluxe_ctx->seq_work2[i],INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 196 ierr = VecScatterEnd(deluxe_ctx->seq_scctx[i],y,deluxe_ctx->seq_work2[i],INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 197 ierr = KSPGetOperators(deluxe_ctx->change[i],&change,NULL);CHKERRQ(ierr); 198 ierr = MatMultTranspose(change,deluxe_ctx->seq_work2[i],deluxe_ctx->seq_work1[i]);CHKERRQ(ierr); 199 } else { 200 ierr = VecScatterBegin(deluxe_ctx->seq_scctx[i],y,deluxe_ctx->seq_work1[i],INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 201 ierr = VecScatterEnd(deluxe_ctx->seq_scctx[i],y,deluxe_ctx->seq_work1[i],INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 202 } 203 if (deluxe_ctx->seq_mat_inv_sum[i]) { 204 PetscScalar *x; 205 206 ierr = VecGetArray(deluxe_ctx->seq_work1[i],&x);CHKERRQ(ierr); 207 ierr = VecPlaceArray(deluxe_ctx->seq_work2[i],x);CHKERRQ(ierr); 208 ierr = VecRestoreArray(deluxe_ctx->seq_work1[i],&x);CHKERRQ(ierr); 209 ierr = MatSolve(deluxe_ctx->seq_mat_inv_sum[i],deluxe_ctx->seq_work1[i],deluxe_ctx->seq_work2[i]);CHKERRQ(ierr); 210 ierr = VecResetArray(deluxe_ctx->seq_work2[i]);CHKERRQ(ierr); 211 } 212 ierr = MatMult(deluxe_ctx->seq_mat[i],deluxe_ctx->seq_work1[i],deluxe_ctx->seq_work2[i]);CHKERRQ(ierr); 213 if (deluxe_ctx->change) { 214 if (deluxe_ctx->change_with_qr) { 215 Mat change; 216 217 ierr = KSPGetOperators(deluxe_ctx->change[i],&change,NULL);CHKERRQ(ierr); 218 ierr = MatMult(change,deluxe_ctx->seq_work2[i],deluxe_ctx->seq_work1[i]);CHKERRQ(ierr); 219 } else { 220 ierr = KSPSolveTranspose(deluxe_ctx->change[i],deluxe_ctx->seq_work2[i],deluxe_ctx->seq_work1[i]);CHKERRQ(ierr); 221 } 222 ierr = VecScatterBegin(deluxe_ctx->seq_scctx[i],deluxe_ctx->seq_work1[i],y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 223 ierr = VecScatterEnd(deluxe_ctx->seq_scctx[i],deluxe_ctx->seq_work1[i],y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 224 } else { 225 ierr = VecScatterBegin(deluxe_ctx->seq_scctx[i],deluxe_ctx->seq_work2[i],y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 226 ierr = VecScatterEnd(deluxe_ctx->seq_scctx[i],deluxe_ctx->seq_work2[i],y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 227 } 228 } 229 } 230 PetscFunctionReturn(0); 231 } 232 233 PetscErrorCode PCBDDCScalingRestriction(PC pc, Vec global_vector, Vec local_interface_vector) 234 { 235 PC_BDDC *pcbddc=(PC_BDDC*)pc->data; 236 PetscErrorCode ierr; 237 238 PetscFunctionBegin; 239 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 240 PetscValidHeaderSpecific(global_vector,VEC_CLASSID,2); 241 PetscValidHeaderSpecific(local_interface_vector,VEC_CLASSID,3); 242 if (local_interface_vector == pcbddc->work_scaling) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Local vector cannot be pcbddc->work_scaling!\n"); 243 ierr = PetscUseMethod(pc,"PCBDDCScalingRestriction_C",(PC,Vec,Vec),(pc,global_vector,local_interface_vector));CHKERRQ(ierr); 244 PetscFunctionReturn(0); 245 } 246 247 PetscErrorCode PCBDDCScalingSetUp(PC pc) 248 { 249 PC_IS* pcis=(PC_IS*)pc->data; 250 PC_BDDC* pcbddc=(PC_BDDC*)pc->data; 251 PetscErrorCode ierr; 252 253 PetscFunctionBegin; 254 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 255 /* create work vector for the operator */ 256 ierr = VecDestroy(&pcbddc->work_scaling);CHKERRQ(ierr); 257 ierr = VecDuplicate(pcis->vec1_B,&pcbddc->work_scaling);CHKERRQ(ierr); 258 /* always rebuild pcis->D */ 259 if (pcis->use_stiffness_scaling) { 260 ierr = MatGetDiagonal(pcbddc->local_mat,pcis->vec1_N);CHKERRQ(ierr); 261 ierr = VecScatterBegin(pcis->N_to_B,pcis->vec1_N,pcis->D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 262 ierr = VecScatterEnd(pcis->N_to_B,pcis->vec1_N,pcis->D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 263 } 264 ierr = VecCopy(pcis->D,pcis->vec1_B);CHKERRQ(ierr); 265 ierr = VecSet(pcis->vec1_global,0.0);CHKERRQ(ierr); 266 ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_B,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 267 ierr = VecScatterEnd(pcis->global_to_B,pcis->vec1_B,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 268 ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_global,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 269 ierr = VecScatterEnd(pcis->global_to_B,pcis->vec1_global,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 270 ierr = VecPointwiseDivide(pcis->D,pcis->D,pcis->vec1_B);CHKERRQ(ierr); 271 /* now setup */ 272 if (pcbddc->use_deluxe_scaling) { 273 if (!pcbddc->deluxe_ctx) { 274 ierr = PCBDDCScalingCreate_Deluxe(pc);CHKERRQ(ierr); 275 } 276 ierr = PCBDDCScalingSetUp_Deluxe(pc);CHKERRQ(ierr); 277 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCScalingRestriction_C",PCBDDCScalingRestriction_Deluxe);CHKERRQ(ierr); 278 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCScalingExtension_C",PCBDDCScalingExtension_Deluxe);CHKERRQ(ierr); 279 } else { 280 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCScalingRestriction_C",PCBDDCScalingRestriction_Basic);CHKERRQ(ierr); 281 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCScalingExtension_C",PCBDDCScalingExtension_Basic);CHKERRQ(ierr); 282 } 283 284 /* test */ 285 if (pcbddc->dbg_flag) { 286 Mat B0_B = NULL; 287 Vec B0_Bv = NULL, B0_Bv2 = NULL; 288 Vec vec2_global; 289 PetscViewer viewer = pcbddc->dbg_viewer; 290 PetscReal error; 291 292 /* extension -> from local to parallel */ 293 ierr = VecSet(pcis->vec1_global,0.0);CHKERRQ(ierr); 294 ierr = VecSetRandom(pcis->vec1_B,NULL);CHKERRQ(ierr); 295 ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_B,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 296 ierr = VecScatterEnd(pcis->global_to_B,pcis->vec1_B,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 297 ierr = VecDuplicate(pcis->vec1_global,&vec2_global);CHKERRQ(ierr); 298 ierr = VecCopy(pcis->vec1_global,vec2_global);CHKERRQ(ierr); 299 ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_global,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 300 ierr = VecScatterEnd(pcis->global_to_B,pcis->vec1_global,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 301 if (pcbddc->benign_n) { 302 IS is_dummy; 303 304 ierr = ISCreateStride(PETSC_COMM_SELF,pcbddc->benign_n,0,1,&is_dummy);CHKERRQ(ierr); 305 ierr = MatCreateSubMatrix(pcbddc->benign_B0,is_dummy,pcis->is_B_local,MAT_INITIAL_MATRIX,&B0_B);CHKERRQ(ierr); 306 ierr = ISDestroy(&is_dummy);CHKERRQ(ierr); 307 ierr = MatCreateVecs(B0_B,NULL,&B0_Bv);CHKERRQ(ierr); 308 ierr = VecDuplicate(B0_Bv,&B0_Bv2);CHKERRQ(ierr); 309 ierr = MatMult(B0_B,pcis->vec1_B,B0_Bv);CHKERRQ(ierr); 310 } 311 ierr = PCBDDCScalingExtension(pc,pcis->vec1_B,pcis->vec1_global);CHKERRQ(ierr); 312 if (pcbddc->benign_saddle_point) { 313 PetscReal errorl = 0.; 314 ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_global,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 315 ierr = VecScatterEnd(pcis->global_to_B,pcis->vec1_global,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 316 if (pcbddc->benign_n) { 317 ierr = MatMult(B0_B,pcis->vec1_B,B0_Bv2);CHKERRQ(ierr); 318 ierr = VecAXPY(B0_Bv,-1.0,B0_Bv2);CHKERRQ(ierr); 319 ierr = VecNorm(B0_Bv,NORM_INFINITY,&errorl);CHKERRQ(ierr); 320 } 321 ierr = MPI_Allreduce(&errorl,&error,1,MPIU_REAL,MPI_SUM,PetscObjectComm((PetscObject)pc));CHKERRQ(ierr); 322 ierr = PetscViewerASCIIPrintf(viewer,"Error benign extension %1.14e\n",error);CHKERRQ(ierr); 323 } 324 ierr = VecAXPY(pcis->vec1_global,-1.0,vec2_global);CHKERRQ(ierr); 325 ierr = VecNorm(pcis->vec1_global,NORM_INFINITY,&error);CHKERRQ(ierr); 326 ierr = PetscViewerASCIIPrintf(viewer,"Error scaling extension %1.14e\n",error);CHKERRQ(ierr); 327 ierr = VecDestroy(&vec2_global);CHKERRQ(ierr); 328 329 /* restriction -> from parallel to local */ 330 ierr = VecSet(pcis->vec1_global,0.0);CHKERRQ(ierr); 331 ierr = VecSetRandom(pcis->vec1_B,NULL);CHKERRQ(ierr); 332 ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_B,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 333 ierr = VecScatterEnd(pcis->global_to_B,pcis->vec1_B,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 334 ierr = PCBDDCScalingRestriction(pc,pcis->vec1_global,pcis->vec1_B);CHKERRQ(ierr); 335 ierr = VecScale(pcis->vec1_B,-1.0);CHKERRQ(ierr); 336 ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_B,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 337 ierr = VecScatterEnd(pcis->global_to_B,pcis->vec1_B,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 338 ierr = VecNorm(pcis->vec1_global,NORM_INFINITY,&error);CHKERRQ(ierr); 339 ierr = PetscViewerASCIIPrintf(viewer,"Error scaling restriction %1.14e\n",error);CHKERRQ(ierr); 340 ierr = MatDestroy(&B0_B);CHKERRQ(ierr); 341 ierr = VecDestroy(&B0_Bv);CHKERRQ(ierr); 342 ierr = VecDestroy(&B0_Bv2);CHKERRQ(ierr); 343 } 344 PetscFunctionReturn(0); 345 } 346 347 PetscErrorCode PCBDDCScalingDestroy(PC pc) 348 { 349 PC_BDDC* pcbddc=(PC_BDDC*)pc->data; 350 PetscErrorCode ierr; 351 352 PetscFunctionBegin; 353 if (pcbddc->deluxe_ctx) { 354 ierr = PCBDDCScalingDestroy_Deluxe(pc);CHKERRQ(ierr); 355 } 356 ierr = VecDestroy(&pcbddc->work_scaling);CHKERRQ(ierr); 357 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCScalingRestriction_C",NULL);CHKERRQ(ierr); 358 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCScalingExtension_C",NULL);CHKERRQ(ierr); 359 PetscFunctionReturn(0); 360 } 361 362 static PetscErrorCode PCBDDCScalingCreate_Deluxe(PC pc) 363 { 364 PC_BDDC* pcbddc=(PC_BDDC*)pc->data; 365 PCBDDCDeluxeScaling deluxe_ctx; 366 PetscErrorCode ierr; 367 368 PetscFunctionBegin; 369 ierr = PetscNew(&deluxe_ctx);CHKERRQ(ierr); 370 pcbddc->deluxe_ctx = deluxe_ctx; 371 PetscFunctionReturn(0); 372 } 373 374 static PetscErrorCode PCBDDCScalingDestroy_Deluxe(PC pc) 375 { 376 PC_BDDC* pcbddc=(PC_BDDC*)pc->data; 377 PetscErrorCode ierr; 378 379 PetscFunctionBegin; 380 ierr = PCBDDCScalingReset_Deluxe_Solvers(pcbddc->deluxe_ctx);CHKERRQ(ierr); 381 ierr = PetscFree(pcbddc->deluxe_ctx);CHKERRQ(ierr); 382 PetscFunctionReturn(0); 383 } 384 385 static PetscErrorCode PCBDDCScalingReset_Deluxe_Solvers(PCBDDCDeluxeScaling deluxe_ctx) 386 { 387 PetscInt i; 388 PetscErrorCode ierr; 389 390 PetscFunctionBegin; 391 ierr = PetscFree(deluxe_ctx->idx_simple_B);CHKERRQ(ierr); 392 deluxe_ctx->n_simple = 0; 393 for (i=0;i<deluxe_ctx->seq_n;i++) { 394 ierr = VecScatterDestroy(&deluxe_ctx->seq_scctx[i]);CHKERRQ(ierr); 395 ierr = VecDestroy(&deluxe_ctx->seq_work1[i]);CHKERRQ(ierr); 396 ierr = VecDestroy(&deluxe_ctx->seq_work2[i]);CHKERRQ(ierr); 397 ierr = MatDestroy(&deluxe_ctx->seq_mat[i]);CHKERRQ(ierr); 398 ierr = MatDestroy(&deluxe_ctx->seq_mat_inv_sum[i]);CHKERRQ(ierr); 399 } 400 ierr = PetscFree5(deluxe_ctx->seq_scctx,deluxe_ctx->seq_work1,deluxe_ctx->seq_work2,deluxe_ctx->seq_mat,deluxe_ctx->seq_mat_inv_sum);CHKERRQ(ierr); 401 ierr = PetscFree(deluxe_ctx->workspace);CHKERRQ(ierr); 402 deluxe_ctx->seq_n = 0; 403 PetscFunctionReturn(0); 404 } 405 406 static PetscErrorCode PCBDDCScalingSetUp_Deluxe(PC pc) 407 { 408 PC_IS *pcis=(PC_IS*)pc->data; 409 PC_BDDC *pcbddc=(PC_BDDC*)pc->data; 410 PCBDDCDeluxeScaling deluxe_ctx=pcbddc->deluxe_ctx; 411 PCBDDCSubSchurs sub_schurs=pcbddc->sub_schurs; 412 PetscErrorCode ierr; 413 414 PetscFunctionBegin; 415 /* reset data structures if the topology has changed */ 416 if (pcbddc->recompute_topography) { 417 ierr = PCBDDCScalingReset_Deluxe_Solvers(deluxe_ctx);CHKERRQ(ierr); 418 } 419 420 /* Compute data structures to solve sequential problems */ 421 ierr = PCBDDCScalingSetUp_Deluxe_Private(pc);CHKERRQ(ierr); 422 423 /* diagonal scaling on interface dofs not contained in cc */ 424 if (sub_schurs->is_vertices || sub_schurs->is_dir) { 425 PetscInt n_com,n_dir; 426 n_com = 0; 427 if (sub_schurs->is_vertices) { 428 ierr = ISGetLocalSize(sub_schurs->is_vertices,&n_com);CHKERRQ(ierr); 429 } 430 n_dir = 0; 431 if (sub_schurs->is_dir) { 432 ierr = ISGetLocalSize(sub_schurs->is_dir,&n_dir);CHKERRQ(ierr); 433 } 434 if (!deluxe_ctx->n_simple) { 435 deluxe_ctx->n_simple = n_dir + n_com; 436 ierr = PetscMalloc1(deluxe_ctx->n_simple,&deluxe_ctx->idx_simple_B);CHKERRQ(ierr); 437 if (sub_schurs->is_vertices) { 438 PetscInt nmap; 439 const PetscInt *idxs; 440 441 ierr = ISGetIndices(sub_schurs->is_vertices,&idxs);CHKERRQ(ierr); 442 ierr = ISGlobalToLocalMappingApply(pcis->BtoNmap,IS_GTOLM_DROP,n_com,idxs,&nmap,deluxe_ctx->idx_simple_B);CHKERRQ(ierr); 443 if (nmap != n_com) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Error when mapping simply scaled dofs (is_vertices)! %d != %d",nmap,n_com); 444 ierr = ISRestoreIndices(sub_schurs->is_vertices,&idxs);CHKERRQ(ierr); 445 } 446 if (sub_schurs->is_dir) { 447 PetscInt nmap; 448 const PetscInt *idxs; 449 450 ierr = ISGetIndices(sub_schurs->is_dir,&idxs);CHKERRQ(ierr); 451 ierr = ISGlobalToLocalMappingApply(pcis->BtoNmap,IS_GTOLM_DROP,n_dir,idxs,&nmap,deluxe_ctx->idx_simple_B+n_com);CHKERRQ(ierr); 452 if (nmap != n_dir) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Error when mapping simply scaled dofs (sub_schurs->is_dir)! %d != %d",nmap,n_dir); 453 ierr = ISRestoreIndices(sub_schurs->is_dir,&idxs);CHKERRQ(ierr); 454 } 455 ierr = PetscSortInt(deluxe_ctx->n_simple,deluxe_ctx->idx_simple_B);CHKERRQ(ierr); 456 } else { 457 if (deluxe_ctx->n_simple != n_dir + n_com) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Number of simply scaled dofs %d is different from the previous one computed %d",n_dir + n_com,deluxe_ctx->n_simple); 458 } 459 } else { 460 deluxe_ctx->n_simple = 0; 461 deluxe_ctx->idx_simple_B = 0; 462 } 463 PetscFunctionReturn(0); 464 } 465 466 static PetscErrorCode PCBDDCScalingSetUp_Deluxe_Private(PC pc) 467 { 468 PC_BDDC *pcbddc=(PC_BDDC*)pc->data; 469 PCBDDCDeluxeScaling deluxe_ctx=pcbddc->deluxe_ctx; 470 PCBDDCSubSchurs sub_schurs = pcbddc->sub_schurs; 471 PetscScalar *matdata,*matdata2; 472 PetscInt i,max_subset_size,cum,cum2; 473 const PetscInt *idxs; 474 PetscBool newsetup = PETSC_FALSE; 475 PetscErrorCode ierr; 476 477 PetscFunctionBegin; 478 if (!sub_schurs->n_subs) { 479 PetscFunctionReturn(0); 480 } 481 482 /* Allocate arrays for subproblems */ 483 if (!deluxe_ctx->seq_n) { 484 deluxe_ctx->seq_n = sub_schurs->n_subs; 485 ierr = PetscCalloc5(deluxe_ctx->seq_n,&deluxe_ctx->seq_scctx,deluxe_ctx->seq_n,&deluxe_ctx->seq_work1,deluxe_ctx->seq_n,&deluxe_ctx->seq_work2,deluxe_ctx->seq_n,&deluxe_ctx->seq_mat,deluxe_ctx->seq_n,&deluxe_ctx->seq_mat_inv_sum);CHKERRQ(ierr); 486 newsetup = PETSC_TRUE; 487 } else if (deluxe_ctx->seq_n != sub_schurs->n_subs) { 488 SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Number of deluxe subproblems %d is different from the sub_schurs %d",deluxe_ctx->seq_n,sub_schurs->n_subs); 489 } 490 /* the change of basis is just a reference to sub_schurs->change (if any) */ 491 deluxe_ctx->change = sub_schurs->change; 492 deluxe_ctx->change_with_qr = sub_schurs->change_with_qr; 493 494 /* Create objects for deluxe */ 495 max_subset_size = 0; 496 for (i=0;i<sub_schurs->n_subs;i++) { 497 PetscInt subset_size; 498 ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr); 499 max_subset_size = PetscMax(subset_size,max_subset_size); 500 } 501 if (newsetup) { 502 ierr = PetscMalloc1(2*max_subset_size,&deluxe_ctx->workspace);CHKERRQ(ierr); 503 } 504 cum = cum2 = 0; 505 ierr = ISGetIndices(sub_schurs->is_Ej_all,&idxs);CHKERRQ(ierr); 506 ierr = MatSeqAIJGetArray(sub_schurs->S_Ej_all,&matdata);CHKERRQ(ierr); 507 ierr = MatSeqAIJGetArray(sub_schurs->sum_S_Ej_all,&matdata2);CHKERRQ(ierr); 508 for (i=0;i<deluxe_ctx->seq_n;i++) { 509 PetscInt subset_size; 510 511 ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr); 512 if (newsetup) { 513 IS sub; 514 /* work vectors */ 515 ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,1,subset_size,deluxe_ctx->workspace,&deluxe_ctx->seq_work1[i]);CHKERRQ(ierr); 516 ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,1,subset_size,deluxe_ctx->workspace+subset_size,&deluxe_ctx->seq_work2[i]);CHKERRQ(ierr); 517 518 /* scatters */ 519 ierr = ISCreateGeneral(PETSC_COMM_SELF,subset_size,idxs+cum,PETSC_COPY_VALUES,&sub);CHKERRQ(ierr); 520 ierr = VecScatterCreate(pcbddc->work_scaling,sub,deluxe_ctx->seq_work1[i],NULL,&deluxe_ctx->seq_scctx[i]);CHKERRQ(ierr); 521 ierr = ISDestroy(&sub);CHKERRQ(ierr); 522 } 523 524 /* S_E_j */ 525 ierr = MatDestroy(&deluxe_ctx->seq_mat[i]);CHKERRQ(ierr); 526 ierr = MatCreateSeqDense(PETSC_COMM_SELF,subset_size,subset_size,matdata+cum2,&deluxe_ctx->seq_mat[i]);CHKERRQ(ierr); 527 528 /* \sum_k S^k_E_j */ 529 ierr = MatDestroy(&deluxe_ctx->seq_mat_inv_sum[i]);CHKERRQ(ierr); 530 ierr = MatCreateSeqDense(PETSC_COMM_SELF,subset_size,subset_size,matdata2+cum2,&deluxe_ctx->seq_mat_inv_sum[i]);CHKERRQ(ierr); 531 ierr = MatSetOption(deluxe_ctx->seq_mat_inv_sum[i],MAT_SPD,sub_schurs->is_posdef);CHKERRQ(ierr); 532 ierr = MatSetOption(deluxe_ctx->seq_mat_inv_sum[i],MAT_HERMITIAN,sub_schurs->is_hermitian);CHKERRQ(ierr); 533 if (sub_schurs->is_hermitian) { 534 ierr = MatCholeskyFactor(deluxe_ctx->seq_mat_inv_sum[i],NULL,NULL);CHKERRQ(ierr); 535 } else { 536 ierr = MatLUFactor(deluxe_ctx->seq_mat_inv_sum[i],NULL,NULL,NULL);CHKERRQ(ierr); 537 } 538 if (pcbddc->deluxe_singlemat) { 539 Mat X,Y; 540 if (!sub_schurs->is_hermitian) { 541 ierr = MatTranspose(deluxe_ctx->seq_mat[i],MAT_INITIAL_MATRIX,&X);CHKERRQ(ierr); 542 } else { 543 ierr = PetscObjectReference((PetscObject)deluxe_ctx->seq_mat[i]);CHKERRQ(ierr); 544 X = deluxe_ctx->seq_mat[i]; 545 } 546 ierr = MatDuplicate(X,MAT_DO_NOT_COPY_VALUES,&Y);CHKERRQ(ierr); 547 if (!sub_schurs->is_hermitian) { 548 ierr = PCBDDCMatTransposeMatSolve_SeqDense(deluxe_ctx->seq_mat_inv_sum[i],X,Y);CHKERRQ(ierr); 549 } else { 550 ierr = MatMatSolve(deluxe_ctx->seq_mat_inv_sum[i],X,Y);CHKERRQ(ierr); 551 } 552 553 ierr = MatDestroy(&deluxe_ctx->seq_mat_inv_sum[i]);CHKERRQ(ierr); 554 ierr = MatDestroy(&deluxe_ctx->seq_mat[i]);CHKERRQ(ierr); 555 ierr = MatDestroy(&X);CHKERRQ(ierr); 556 if (deluxe_ctx->change) { 557 Mat C,CY; 558 559 if (!deluxe_ctx->change_with_qr) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Only QR based change of basis"); 560 ierr = KSPGetOperators(deluxe_ctx->change[i],&C,NULL);CHKERRQ(ierr); 561 ierr = MatMatMult(C,Y,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&CY);CHKERRQ(ierr); 562 ierr = MatMatTransposeMult(CY,C,MAT_REUSE_MATRIX,PETSC_DEFAULT,&Y);CHKERRQ(ierr); 563 ierr = MatDestroy(&CY);CHKERRQ(ierr); 564 } 565 ierr = MatTranspose(Y,MAT_INPLACE_MATRIX,&Y);CHKERRQ(ierr); 566 deluxe_ctx->seq_mat[i] = Y; 567 } 568 cum += subset_size; 569 cum2 += subset_size*subset_size; 570 } 571 ierr = ISRestoreIndices(sub_schurs->is_Ej_all,&idxs);CHKERRQ(ierr); 572 ierr = MatSeqAIJRestoreArray(sub_schurs->S_Ej_all,&matdata);CHKERRQ(ierr); 573 ierr = MatSeqAIJRestoreArray(sub_schurs->sum_S_Ej_all,&matdata2);CHKERRQ(ierr); 574 if (pcbddc->deluxe_singlemat) { 575 deluxe_ctx->change = NULL; 576 deluxe_ctx->change_with_qr = PETSC_FALSE; 577 } 578 579 if (deluxe_ctx->change && !deluxe_ctx->change_with_qr) { 580 for (i=0;i<deluxe_ctx->seq_n;i++) { 581 if (newsetup) { 582 PC pc; 583 584 ierr = KSPGetPC(deluxe_ctx->change[i],&pc);CHKERRQ(ierr); 585 ierr = PCSetType(pc,PCLU);CHKERRQ(ierr); 586 ierr = KSPSetFromOptions(deluxe_ctx->change[i]);CHKERRQ(ierr); 587 } 588 ierr = KSPSetUp(deluxe_ctx->change[i]);CHKERRQ(ierr); 589 } 590 } 591 PetscFunctionReturn(0); 592 } 593