1 #include <../src/ksp/pc/impls/bddc/bddc.h> 2 #include <../src/ksp/pc/impls/bddc/bddcprivate.h> 3 #include <petscblaslapack.h> 4 5 static PetscErrorCode MatMult_BDdelta_deluxe_nonred(Mat A, Vec x, Vec y) 6 { 7 BDdelta_DN ctx; 8 9 PetscFunctionBegin; 10 PetscCall(MatShellGetContext(A,&ctx)); 11 PetscCall(MatMultTranspose(ctx->BD,x,ctx->work)); 12 PetscCall(KSPSolveTranspose(ctx->kBD,ctx->work,y)); 13 /* No PC so cannot propagate up failure in KSPSolveTranspose() */ 14 PetscFunctionReturn(0); 15 } 16 17 static PetscErrorCode MatMultTranspose_BDdelta_deluxe_nonred(Mat A, Vec x, Vec y) 18 { 19 BDdelta_DN ctx; 20 21 PetscFunctionBegin; 22 PetscCall(MatShellGetContext(A,&ctx)); 23 PetscCall(KSPSolve(ctx->kBD,x,ctx->work)); 24 /* No PC so cannot propagate up failure in KSPSolve() */ 25 PetscCall(MatMult(ctx->BD,ctx->work,y)); 26 PetscFunctionReturn(0); 27 } 28 29 static PetscErrorCode MatDestroy_BDdelta_deluxe_nonred(Mat A) 30 { 31 BDdelta_DN ctx; 32 33 PetscFunctionBegin; 34 PetscCall(MatShellGetContext(A,&ctx)); 35 PetscCall(MatDestroy(&ctx->BD)); 36 PetscCall(KSPDestroy(&ctx->kBD)); 37 PetscCall(VecDestroy(&ctx->work)); 38 PetscCall(PetscFree(ctx)); 39 PetscFunctionReturn(0); 40 } 41 42 PetscErrorCode PCBDDCCreateFETIDPMatContext(PC pc, FETIDPMat_ctx *fetidpmat_ctx) 43 { 44 FETIDPMat_ctx newctx; 45 46 PetscFunctionBegin; 47 PetscCall(PetscNew(&newctx)); 48 /* increase the reference count for BDDC preconditioner */ 49 PetscCall(PetscObjectReference((PetscObject)pc)); 50 newctx->pc = pc; 51 *fetidpmat_ctx = newctx; 52 PetscFunctionReturn(0); 53 } 54 55 PetscErrorCode PCBDDCCreateFETIDPPCContext(PC pc, FETIDPPC_ctx *fetidppc_ctx) 56 { 57 FETIDPPC_ctx newctx; 58 59 PetscFunctionBegin; 60 PetscCall(PetscNew(&newctx)); 61 /* increase the reference count for BDDC preconditioner */ 62 PetscCall(PetscObjectReference((PetscObject)pc)); 63 newctx->pc = pc; 64 *fetidppc_ctx = newctx; 65 PetscFunctionReturn(0); 66 } 67 68 PetscErrorCode PCBDDCDestroyFETIDPMat(Mat A) 69 { 70 FETIDPMat_ctx mat_ctx; 71 72 PetscFunctionBegin; 73 PetscCall(MatShellGetContext(A,&mat_ctx)); 74 PetscCall(VecDestroy(&mat_ctx->lambda_local)); 75 PetscCall(VecDestroy(&mat_ctx->temp_solution_D)); 76 PetscCall(VecDestroy(&mat_ctx->temp_solution_B)); 77 PetscCall(MatDestroy(&mat_ctx->B_delta)); 78 PetscCall(MatDestroy(&mat_ctx->B_Ddelta)); 79 PetscCall(MatDestroy(&mat_ctx->B_BB)); 80 PetscCall(MatDestroy(&mat_ctx->B_BI)); 81 PetscCall(MatDestroy(&mat_ctx->Bt_BB)); 82 PetscCall(MatDestroy(&mat_ctx->Bt_BI)); 83 PetscCall(MatDestroy(&mat_ctx->C)); 84 PetscCall(VecDestroy(&mat_ctx->rhs_flip)); 85 PetscCall(VecDestroy(&mat_ctx->vP)); 86 PetscCall(VecDestroy(&mat_ctx->xPg)); 87 PetscCall(VecDestroy(&mat_ctx->yPg)); 88 PetscCall(VecScatterDestroy(&mat_ctx->l2g_lambda)); 89 PetscCall(VecScatterDestroy(&mat_ctx->l2g_lambda_only)); 90 PetscCall(VecScatterDestroy(&mat_ctx->l2g_p)); 91 PetscCall(VecScatterDestroy(&mat_ctx->g2g_p)); 92 PetscCall(PCDestroy(&mat_ctx->pc)); /* decrease PCBDDC reference count */ 93 PetscCall(ISDestroy(&mat_ctx->pressure)); 94 PetscCall(ISDestroy(&mat_ctx->lagrange)); 95 PetscCall(PetscFree(mat_ctx)); 96 PetscFunctionReturn(0); 97 } 98 99 PetscErrorCode PCBDDCDestroyFETIDPPC(PC pc) 100 { 101 FETIDPPC_ctx pc_ctx; 102 103 PetscFunctionBegin; 104 PetscCall(PCShellGetContext(pc,&pc_ctx)); 105 PetscCall(VecDestroy(&pc_ctx->lambda_local)); 106 PetscCall(MatDestroy(&pc_ctx->B_Ddelta)); 107 PetscCall(VecScatterDestroy(&pc_ctx->l2g_lambda)); 108 PetscCall(MatDestroy(&pc_ctx->S_j)); 109 PetscCall(PCDestroy(&pc_ctx->pc)); /* decrease PCBDDC reference count */ 110 PetscCall(VecDestroy(&pc_ctx->xPg)); 111 PetscCall(VecDestroy(&pc_ctx->yPg)); 112 PetscCall(PetscFree(pc_ctx)); 113 PetscFunctionReturn(0); 114 } 115 116 PetscErrorCode PCBDDCSetupFETIDPMatContext(FETIDPMat_ctx fetidpmat_ctx) 117 { 118 PC_IS *pcis=(PC_IS*)fetidpmat_ctx->pc->data; 119 PC_BDDC *pcbddc=(PC_BDDC*)fetidpmat_ctx->pc->data; 120 PCBDDCGraph mat_graph=pcbddc->mat_graph; 121 Mat_IS *matis = (Mat_IS*)fetidpmat_ctx->pc->pmat->data; 122 MPI_Comm comm; 123 Mat ScalingMat,BD1,BD2; 124 Vec fetidp_global; 125 IS IS_l2g_lambda; 126 IS subset,subset_mult,subset_n,isvert; 127 PetscBool skip_node,fully_redundant; 128 PetscInt i,j,k,s,n_boundary_dofs,n_global_lambda,n_vertices,partial_sum; 129 PetscInt cum,n_local_lambda,n_lambda_for_dof,dual_size,n_neg_values,n_pos_values; 130 PetscMPIInt rank,size,buf_size,neigh; 131 PetscScalar scalar_value; 132 const PetscInt *vertex_indices; 133 PetscInt *dual_dofs_boundary_indices,*aux_local_numbering_1; 134 const PetscInt *aux_global_numbering; 135 PetscInt *aux_sums,*cols_B_delta,*l2g_indices; 136 PetscScalar *array,*scaling_factors,*vals_B_delta; 137 PetscScalar **all_factors; 138 PetscInt *aux_local_numbering_2; 139 PetscLayout llay; 140 141 /* saddlepoint */ 142 ISLocalToGlobalMapping l2gmap_p; 143 PetscLayout play; 144 IS gP,pP; 145 PetscInt nPl,nPg,nPgl; 146 147 PetscFunctionBegin; 148 PetscCall(PetscObjectGetComm((PetscObject)(fetidpmat_ctx->pc),&comm)); 149 PetscCallMPI(MPI_Comm_rank(comm,&rank)); 150 PetscCallMPI(MPI_Comm_size(comm,&size)); 151 152 /* saddlepoint */ 153 nPl = 0; 154 nPg = 0; 155 nPgl = 0; 156 gP = NULL; 157 pP = NULL; 158 l2gmap_p = NULL; 159 play = NULL; 160 PetscCall(PetscObjectQuery((PetscObject)fetidpmat_ctx->pc,"__KSPFETIDP_pP",(PetscObject*)&pP)); 161 if (pP) { /* saddle point */ 162 /* subdomain pressures in global numbering */ 163 PetscCall(PetscObjectQuery((PetscObject)fetidpmat_ctx->pc,"__KSPFETIDP_gP",(PetscObject*)&gP)); 164 PetscCheck(gP,PETSC_COMM_SELF,PETSC_ERR_PLIB,"gP not present"); 165 PetscCall(ISGetLocalSize(gP,&nPl)); 166 PetscCall(VecCreate(PETSC_COMM_SELF,&fetidpmat_ctx->vP)); 167 PetscCall(VecSetSizes(fetidpmat_ctx->vP,nPl,nPl)); 168 PetscCall(VecSetType(fetidpmat_ctx->vP,VECSTANDARD)); 169 PetscCall(VecSetUp(fetidpmat_ctx->vP)); 170 171 /* pressure matrix */ 172 PetscCall(PetscObjectQuery((PetscObject)fetidpmat_ctx->pc,"__KSPFETIDP_C",(PetscObject*)&fetidpmat_ctx->C)); 173 if (!fetidpmat_ctx->C) { /* null pressure block, compute layout and global numbering for pressures */ 174 IS Pg; 175 176 PetscCall(ISRenumber(gP,NULL,&nPg,&Pg)); 177 PetscCall(ISLocalToGlobalMappingCreateIS(Pg,&l2gmap_p)); 178 PetscCall(ISDestroy(&Pg)); 179 PetscCall(PetscLayoutCreate(comm,&play)); 180 PetscCall(PetscLayoutSetBlockSize(play,1)); 181 PetscCall(PetscLayoutSetSize(play,nPg)); 182 PetscCall(ISGetLocalSize(pP,&nPgl)); 183 PetscCall(PetscLayoutSetLocalSize(play,nPgl)); 184 PetscCall(PetscLayoutSetUp(play)); 185 } else { 186 PetscCall(PetscObjectReference((PetscObject)fetidpmat_ctx->C)); 187 PetscCall(MatISGetLocalToGlobalMapping(fetidpmat_ctx->C,&l2gmap_p,NULL)); 188 PetscCall(PetscObjectReference((PetscObject)l2gmap_p)); 189 PetscCall(MatGetSize(fetidpmat_ctx->C,&nPg,NULL)); 190 PetscCall(MatGetLocalSize(fetidpmat_ctx->C,NULL,&nPgl)); 191 PetscCall(MatGetLayouts(fetidpmat_ctx->C,NULL,&llay)); 192 PetscCall(PetscLayoutReference(llay,&play)); 193 } 194 PetscCall(VecCreateMPIWithArray(comm,1,nPgl,nPg,NULL,&fetidpmat_ctx->xPg)); 195 PetscCall(VecCreateMPIWithArray(comm,1,nPgl,nPg,NULL,&fetidpmat_ctx->yPg)); 196 197 /* import matrices for pressures coupling */ 198 PetscCall(PetscObjectQuery((PetscObject)fetidpmat_ctx->pc,"__KSPFETIDP_B_BI",(PetscObject*)&fetidpmat_ctx->B_BI)); 199 PetscCheck(fetidpmat_ctx->B_BI,PETSC_COMM_SELF,PETSC_ERR_PLIB,"B_BI not present"); 200 PetscCall(PetscObjectReference((PetscObject)fetidpmat_ctx->B_BI)); 201 202 PetscCall(PetscObjectQuery((PetscObject)fetidpmat_ctx->pc,"__KSPFETIDP_B_BB",(PetscObject*)&fetidpmat_ctx->B_BB)); 203 PetscCheck(fetidpmat_ctx->B_BB,PETSC_COMM_SELF,PETSC_ERR_PLIB,"B_BB not present"); 204 PetscCall(PetscObjectReference((PetscObject)fetidpmat_ctx->B_BB)); 205 206 PetscCall(PetscObjectQuery((PetscObject)fetidpmat_ctx->pc,"__KSPFETIDP_Bt_BI",(PetscObject*)&fetidpmat_ctx->Bt_BI)); 207 PetscCheck(fetidpmat_ctx->Bt_BI,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Bt_BI not present"); 208 PetscCall(PetscObjectReference((PetscObject)fetidpmat_ctx->Bt_BI)); 209 210 PetscCall(PetscObjectQuery((PetscObject)fetidpmat_ctx->pc,"__KSPFETIDP_Bt_BB",(PetscObject*)&fetidpmat_ctx->Bt_BB)); 211 PetscCheck(fetidpmat_ctx->Bt_BB,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Bt_BB not present"); 212 PetscCall(PetscObjectReference((PetscObject)fetidpmat_ctx->Bt_BB)); 213 214 PetscCall(PetscObjectQuery((PetscObject)fetidpmat_ctx->pc,"__KSPFETIDP_flip" ,(PetscObject*)&fetidpmat_ctx->rhs_flip)); 215 if (fetidpmat_ctx->rhs_flip) { 216 PetscCall(PetscObjectReference((PetscObject)fetidpmat_ctx->rhs_flip)); 217 } 218 } 219 220 /* Default type of lagrange multipliers is non-redundant */ 221 fully_redundant = fetidpmat_ctx->fully_redundant; 222 223 /* Evaluate local and global number of lagrange multipliers */ 224 PetscCall(VecSet(pcis->vec1_N,0.0)); 225 n_local_lambda = 0; 226 partial_sum = 0; 227 n_boundary_dofs = 0; 228 s = 0; 229 230 /* Get Vertices used to define the BDDC */ 231 PetscCall(PCBDDCGraphGetCandidatesIS(pcbddc->mat_graph,NULL,NULL,NULL,NULL,&isvert)); 232 PetscCall(ISGetLocalSize(isvert,&n_vertices)); 233 PetscCall(ISGetIndices(isvert,&vertex_indices)); 234 235 dual_size = pcis->n_B-n_vertices; 236 PetscCall(PetscMalloc1(dual_size,&dual_dofs_boundary_indices)); 237 PetscCall(PetscMalloc1(dual_size,&aux_local_numbering_1)); 238 PetscCall(PetscMalloc1(dual_size,&aux_local_numbering_2)); 239 240 PetscCall(VecGetArray(pcis->vec1_N,&array)); 241 for (i=0;i<pcis->n;i++) { 242 j = mat_graph->count[i]; /* RECALL: mat_graph->count[i] does not count myself */ 243 if (j > 0) n_boundary_dofs++; 244 skip_node = PETSC_FALSE; 245 if (s < n_vertices && vertex_indices[s] == i) { /* it works for a sorted set of vertices */ 246 skip_node = PETSC_TRUE; 247 s++; 248 } 249 if (j < 1) skip_node = PETSC_TRUE; 250 if (mat_graph->special_dof[i] == PCBDDCGRAPH_DIRICHLET_MARK) skip_node = PETSC_TRUE; 251 if (!skip_node) { 252 if (fully_redundant) { 253 /* fully redundant set of lagrange multipliers */ 254 n_lambda_for_dof = (j*(j+1))/2; 255 } else { 256 n_lambda_for_dof = j; 257 } 258 n_local_lambda += j; 259 /* needed to evaluate global number of lagrange multipliers */ 260 array[i]=(1.0*n_lambda_for_dof)/(j+1.0); /* already scaled for the next global sum */ 261 /* store some data needed */ 262 dual_dofs_boundary_indices[partial_sum] = n_boundary_dofs-1; 263 aux_local_numbering_1[partial_sum] = i; 264 aux_local_numbering_2[partial_sum] = n_lambda_for_dof; 265 partial_sum++; 266 } 267 } 268 PetscCall(VecRestoreArray(pcis->vec1_N,&array)); 269 PetscCall(ISRestoreIndices(isvert,&vertex_indices)); 270 PetscCall(PCBDDCGraphRestoreCandidatesIS(pcbddc->mat_graph,NULL,NULL,NULL,NULL,&isvert)); 271 dual_size = partial_sum; 272 273 /* compute global ordering of lagrange multipliers and associate l2g map */ 274 PetscCall(ISCreateGeneral(comm,partial_sum,aux_local_numbering_1,PETSC_COPY_VALUES,&subset_n)); 275 PetscCall(ISLocalToGlobalMappingApplyIS(pcis->mapping,subset_n,&subset)); 276 PetscCall(ISDestroy(&subset_n)); 277 PetscCall(ISCreateGeneral(comm,partial_sum,aux_local_numbering_2,PETSC_OWN_POINTER,&subset_mult)); 278 PetscCall(ISRenumber(subset,subset_mult,&fetidpmat_ctx->n_lambda,&subset_n)); 279 PetscCall(ISDestroy(&subset)); 280 281 if (PetscDefined(USE_DEBUG)) { 282 PetscCall(VecSet(pcis->vec1_global,0.0)); 283 PetscCall(VecScatterBegin(matis->rctx,pcis->vec1_N,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE)); 284 PetscCall(VecScatterEnd(matis->rctx,pcis->vec1_N,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE)); 285 PetscCall(VecSum(pcis->vec1_global,&scalar_value)); 286 i = (PetscInt)PetscRealPart(scalar_value); 287 PetscCheck(i == fetidpmat_ctx->n_lambda,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"Global number of multipliers mismatch! (%D != %D)",fetidpmat_ctx->n_lambda,i); 288 } 289 290 /* init data for scaling factors exchange */ 291 if (!pcbddc->use_deluxe_scaling) { 292 PetscInt *ptrs_buffer,neigh_position; 293 PetscScalar *send_buffer,*recv_buffer; 294 MPI_Request *send_reqs,*recv_reqs; 295 296 partial_sum = 0; 297 PetscCall(PetscMalloc1(pcis->n_neigh,&ptrs_buffer)); 298 PetscCall(PetscMalloc1(PetscMax(pcis->n_neigh-1,0),&send_reqs)); 299 PetscCall(PetscMalloc1(PetscMax(pcis->n_neigh-1,0),&recv_reqs)); 300 PetscCall(PetscMalloc1(pcis->n+1,&all_factors)); 301 if (pcis->n_neigh > 0) ptrs_buffer[0]=0; 302 for (i=1;i<pcis->n_neigh;i++) { 303 partial_sum += pcis->n_shared[i]; 304 ptrs_buffer[i] = ptrs_buffer[i-1]+pcis->n_shared[i]; 305 } 306 PetscCall(PetscMalloc1(partial_sum,&send_buffer)); 307 PetscCall(PetscMalloc1(partial_sum,&recv_buffer)); 308 PetscCall(PetscMalloc1(partial_sum,&all_factors[0])); 309 for (i=0;i<pcis->n-1;i++) { 310 j = mat_graph->count[i]; 311 all_factors[i+1]=all_factors[i]+j; 312 } 313 314 /* scatter B scaling to N vec */ 315 PetscCall(VecScatterBegin(pcis->N_to_B,pcis->D,pcis->vec1_N,INSERT_VALUES,SCATTER_REVERSE)); 316 PetscCall(VecScatterEnd(pcis->N_to_B,pcis->D,pcis->vec1_N,INSERT_VALUES,SCATTER_REVERSE)); 317 /* communications */ 318 PetscCall(VecGetArrayRead(pcis->vec1_N,(const PetscScalar**)&array)); 319 for (i=1;i<pcis->n_neigh;i++) { 320 for (j=0;j<pcis->n_shared[i];j++) { 321 send_buffer[ptrs_buffer[i-1]+j]=array[pcis->shared[i][j]]; 322 } 323 PetscCall(PetscMPIIntCast(ptrs_buffer[i]-ptrs_buffer[i-1],&buf_size)); 324 PetscCall(PetscMPIIntCast(pcis->neigh[i],&neigh)); 325 PetscCallMPI(MPI_Isend(&send_buffer[ptrs_buffer[i-1]],buf_size,MPIU_SCALAR,neigh,0,comm,&send_reqs[i-1])); 326 PetscCallMPI(MPI_Irecv(&recv_buffer[ptrs_buffer[i-1]],buf_size,MPIU_SCALAR,neigh,0,comm,&recv_reqs[i-1])); 327 } 328 PetscCall(VecRestoreArrayRead(pcis->vec1_N,(const PetscScalar**)&array)); 329 if (pcis->n_neigh > 0) { 330 PetscCallMPI(MPI_Waitall(pcis->n_neigh-1,recv_reqs,MPI_STATUSES_IGNORE)); 331 } 332 /* put values in correct places */ 333 for (i=1;i<pcis->n_neigh;i++) { 334 for (j=0;j<pcis->n_shared[i];j++) { 335 k = pcis->shared[i][j]; 336 neigh_position = 0; 337 while (mat_graph->neighbours_set[k][neigh_position] != pcis->neigh[i]) {neigh_position++;} 338 all_factors[k][neigh_position]=recv_buffer[ptrs_buffer[i-1]+j]; 339 } 340 } 341 if (pcis->n_neigh > 0) { 342 PetscCallMPI(MPI_Waitall(pcis->n_neigh-1,send_reqs,MPI_STATUSES_IGNORE)); 343 } 344 PetscCall(PetscFree(send_reqs)); 345 PetscCall(PetscFree(recv_reqs)); 346 PetscCall(PetscFree(send_buffer)); 347 PetscCall(PetscFree(recv_buffer)); 348 PetscCall(PetscFree(ptrs_buffer)); 349 } 350 351 /* Compute B and B_delta (local actions) */ 352 PetscCall(PetscMalloc1(pcis->n_neigh,&aux_sums)); 353 PetscCall(PetscMalloc1(n_local_lambda,&l2g_indices)); 354 PetscCall(PetscMalloc1(n_local_lambda,&vals_B_delta)); 355 PetscCall(PetscMalloc1(n_local_lambda,&cols_B_delta)); 356 if (!pcbddc->use_deluxe_scaling) { 357 PetscCall(PetscMalloc1(n_local_lambda,&scaling_factors)); 358 } else { 359 scaling_factors = NULL; 360 all_factors = NULL; 361 } 362 PetscCall(ISGetIndices(subset_n,&aux_global_numbering)); 363 partial_sum=0; 364 cum = 0; 365 for (i=0;i<dual_size;i++) { 366 n_global_lambda = aux_global_numbering[cum]; 367 j = mat_graph->count[aux_local_numbering_1[i]]; 368 aux_sums[0]=0; 369 for (s=1;s<j;s++) { 370 aux_sums[s]=aux_sums[s-1]+j-s+1; 371 } 372 if (all_factors) array = all_factors[aux_local_numbering_1[i]]; 373 n_neg_values = 0; 374 while (n_neg_values < j && mat_graph->neighbours_set[aux_local_numbering_1[i]][n_neg_values] < rank) {n_neg_values++;} 375 n_pos_values = j - n_neg_values; 376 if (fully_redundant) { 377 for (s=0;s<n_neg_values;s++) { 378 l2g_indices [partial_sum+s]=aux_sums[s]+n_neg_values-s-1+n_global_lambda; 379 cols_B_delta [partial_sum+s]=dual_dofs_boundary_indices[i]; 380 vals_B_delta [partial_sum+s]=-1.0; 381 if (!pcbddc->use_deluxe_scaling) scaling_factors[partial_sum+s]=array[s]; 382 } 383 for (s=0;s<n_pos_values;s++) { 384 l2g_indices [partial_sum+s+n_neg_values]=aux_sums[n_neg_values]+s+n_global_lambda; 385 cols_B_delta [partial_sum+s+n_neg_values]=dual_dofs_boundary_indices[i]; 386 vals_B_delta [partial_sum+s+n_neg_values]=1.0; 387 if (!pcbddc->use_deluxe_scaling) scaling_factors[partial_sum+s+n_neg_values]=array[s+n_neg_values]; 388 } 389 partial_sum += j; 390 } else { 391 /* l2g_indices and default cols and vals of B_delta */ 392 for (s=0;s<j;s++) { 393 l2g_indices [partial_sum+s]=n_global_lambda+s; 394 cols_B_delta [partial_sum+s]=dual_dofs_boundary_indices[i]; 395 vals_B_delta [partial_sum+s]=0.0; 396 } 397 /* B_delta */ 398 if (n_neg_values > 0) { /* there's a rank next to me to the left */ 399 vals_B_delta [partial_sum+n_neg_values-1]=-1.0; 400 } 401 if (n_neg_values < j) { /* there's a rank next to me to the right */ 402 vals_B_delta [partial_sum+n_neg_values]=1.0; 403 } 404 /* scaling as in Klawonn-Widlund 1999 */ 405 if (!pcbddc->use_deluxe_scaling) { 406 for (s=0;s<n_neg_values;s++) { 407 scalar_value = 0.0; 408 for (k=0;k<s+1;k++) scalar_value += array[k]; 409 scaling_factors[partial_sum+s] = -scalar_value; 410 } 411 for (s=0;s<n_pos_values;s++) { 412 scalar_value = 0.0; 413 for (k=s+n_neg_values;k<j;k++) scalar_value += array[k]; 414 scaling_factors[partial_sum+s+n_neg_values] = scalar_value; 415 } 416 } 417 partial_sum += j; 418 } 419 cum += aux_local_numbering_2[i]; 420 } 421 PetscCall(ISRestoreIndices(subset_n,&aux_global_numbering)); 422 PetscCall(ISDestroy(&subset_mult)); 423 PetscCall(ISDestroy(&subset_n)); 424 PetscCall(PetscFree(aux_sums)); 425 PetscCall(PetscFree(aux_local_numbering_1)); 426 PetscCall(PetscFree(dual_dofs_boundary_indices)); 427 if (all_factors) { 428 PetscCall(PetscFree(all_factors[0])); 429 PetscCall(PetscFree(all_factors)); 430 } 431 432 /* Create local part of B_delta */ 433 PetscCall(MatCreate(PETSC_COMM_SELF,&fetidpmat_ctx->B_delta)); 434 PetscCall(MatSetSizes(fetidpmat_ctx->B_delta,n_local_lambda,pcis->n_B,n_local_lambda,pcis->n_B)); 435 PetscCall(MatSetType(fetidpmat_ctx->B_delta,MATSEQAIJ)); 436 PetscCall(MatSeqAIJSetPreallocation(fetidpmat_ctx->B_delta,1,NULL)); 437 PetscCall(MatSetOption(fetidpmat_ctx->B_delta,MAT_IGNORE_ZERO_ENTRIES,PETSC_TRUE)); 438 for (i=0;i<n_local_lambda;i++) { 439 PetscCall(MatSetValue(fetidpmat_ctx->B_delta,i,cols_B_delta[i],vals_B_delta[i],INSERT_VALUES)); 440 } 441 PetscCall(PetscFree(vals_B_delta)); 442 PetscCall(MatAssemblyBegin(fetidpmat_ctx->B_delta,MAT_FINAL_ASSEMBLY)); 443 PetscCall(MatAssemblyEnd(fetidpmat_ctx->B_delta,MAT_FINAL_ASSEMBLY)); 444 445 BD1 = NULL; 446 BD2 = NULL; 447 if (fully_redundant) { 448 PetscCheck(!pcbddc->use_deluxe_scaling,comm,PETSC_ERR_SUP,"Deluxe FETIDP with fully-redundant multipliers to be implemented"); 449 PetscCall(MatCreate(PETSC_COMM_SELF,&ScalingMat)); 450 PetscCall(MatSetSizes(ScalingMat,n_local_lambda,n_local_lambda,n_local_lambda,n_local_lambda)); 451 PetscCall(MatSetType(ScalingMat,MATSEQAIJ)); 452 PetscCall(MatSeqAIJSetPreallocation(ScalingMat,1,NULL)); 453 for (i=0;i<n_local_lambda;i++) { 454 PetscCall(MatSetValue(ScalingMat,i,i,scaling_factors[i],INSERT_VALUES)); 455 } 456 PetscCall(MatAssemblyBegin(ScalingMat,MAT_FINAL_ASSEMBLY)); 457 PetscCall(MatAssemblyEnd(ScalingMat,MAT_FINAL_ASSEMBLY)); 458 PetscCall(MatMatMult(ScalingMat,fetidpmat_ctx->B_delta,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&fetidpmat_ctx->B_Ddelta)); 459 PetscCall(MatDestroy(&ScalingMat)); 460 } else { 461 PetscCall(MatCreate(PETSC_COMM_SELF,&fetidpmat_ctx->B_Ddelta)); 462 PetscCall(MatSetSizes(fetidpmat_ctx->B_Ddelta,n_local_lambda,pcis->n_B,n_local_lambda,pcis->n_B)); 463 if (!pcbddc->use_deluxe_scaling || !pcbddc->sub_schurs) { 464 PetscCall(MatSetType(fetidpmat_ctx->B_Ddelta,MATSEQAIJ)); 465 PetscCall(MatSeqAIJSetPreallocation(fetidpmat_ctx->B_Ddelta,1,NULL)); 466 for (i=0;i<n_local_lambda;i++) { 467 PetscCall(MatSetValue(fetidpmat_ctx->B_Ddelta,i,cols_B_delta[i],scaling_factors[i],INSERT_VALUES)); 468 } 469 PetscCall(MatAssemblyBegin(fetidpmat_ctx->B_Ddelta,MAT_FINAL_ASSEMBLY)); 470 PetscCall(MatAssemblyEnd(fetidpmat_ctx->B_Ddelta,MAT_FINAL_ASSEMBLY)); 471 } else { 472 /* scaling as in Klawonn-Widlund 1999 */ 473 PCBDDCDeluxeScaling deluxe_ctx = pcbddc->deluxe_ctx; 474 PCBDDCSubSchurs sub_schurs = pcbddc->sub_schurs; 475 Mat T; 476 PetscScalar *W,lwork,*Bwork; 477 const PetscInt *idxs = NULL; 478 PetscInt cum,mss,*nnz; 479 PetscBLASInt *pivots,B_lwork,B_N,B_ierr; 480 481 PetscCheck(pcbddc->deluxe_singlemat,comm,PETSC_ERR_USER,"Cannot compute B_Ddelta! rerun with -pc_bddc_deluxe_singlemat"); 482 mss = 0; 483 PetscCall(PetscCalloc1(pcis->n_B,&nnz)); 484 if (sub_schurs->is_Ej_all) { 485 PetscCall(ISGetIndices(sub_schurs->is_Ej_all,&idxs)); 486 for (i=0,cum=0;i<sub_schurs->n_subs;i++) { 487 PetscInt subset_size; 488 489 PetscCall(ISGetLocalSize(sub_schurs->is_subs[i],&subset_size)); 490 for (j=0;j<subset_size;j++) nnz[idxs[j+cum]] = subset_size; 491 mss = PetscMax(mss,subset_size); 492 cum += subset_size; 493 } 494 } 495 PetscCall(MatCreate(PETSC_COMM_SELF,&T)); 496 PetscCall(MatSetSizes(T,pcis->n_B,pcis->n_B,pcis->n_B,pcis->n_B)); 497 PetscCall(MatSetType(T,MATSEQAIJ)); 498 PetscCall(MatSeqAIJSetPreallocation(T,0,nnz)); 499 PetscCall(PetscFree(nnz)); 500 501 /* workspace allocation */ 502 B_lwork = 0; 503 if (mss) { 504 PetscScalar dummy = 1; 505 506 B_lwork = -1; 507 PetscCall(PetscBLASIntCast(mss,&B_N)); 508 PetscCall(PetscFPTrapPush(PETSC_FP_TRAP_OFF)); 509 PetscStackCallBLAS("LAPACKgetri",LAPACKgetri_(&B_N,&dummy,&B_N,&B_N,&lwork,&B_lwork,&B_ierr)); 510 PetscCall(PetscFPTrapPop()); 511 PetscCheck(!B_ierr,PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in query to GETRI Lapack routine %d",(int)B_ierr); 512 PetscCall(PetscBLASIntCast((PetscInt)PetscRealPart(lwork),&B_lwork)); 513 } 514 PetscCall(PetscMalloc3(mss*mss,&W,mss,&pivots,B_lwork,&Bwork)); 515 516 for (i=0,cum=0;i<sub_schurs->n_subs;i++) { 517 const PetscScalar *M; 518 PetscInt subset_size; 519 520 PetscCall(ISGetLocalSize(sub_schurs->is_subs[i],&subset_size)); 521 PetscCall(PetscBLASIntCast(subset_size,&B_N)); 522 PetscCall(MatDenseGetArrayRead(deluxe_ctx->seq_mat[i],&M)); 523 PetscCall(PetscArraycpy(W,M,subset_size*subset_size)); 524 PetscCall(MatDenseRestoreArrayRead(deluxe_ctx->seq_mat[i],&M)); 525 PetscCall(PetscFPTrapPush(PETSC_FP_TRAP_OFF)); 526 PetscStackCallBLAS("LAPACKgetrf",LAPACKgetrf_(&B_N,&B_N,W,&B_N,pivots,&B_ierr)); 527 PetscCheck(!B_ierr,PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRF Lapack routine %d",(int)B_ierr); 528 PetscStackCallBLAS("LAPACKgetri",LAPACKgetri_(&B_N,W,&B_N,pivots,Bwork,&B_lwork,&B_ierr)); 529 PetscCheck(!B_ierr,PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRI Lapack routine %d",(int)B_ierr); 530 PetscCall(PetscFPTrapPop()); 531 /* silent static analyzer */ 532 PetscCheck(idxs,PETSC_COMM_SELF,PETSC_ERR_PLIB,"IDXS not present"); 533 PetscCall(MatSetValues(T,subset_size,idxs+cum,subset_size,idxs+cum,W,INSERT_VALUES)); 534 cum += subset_size; 535 } 536 PetscCall(MatAssemblyBegin(T,MAT_FINAL_ASSEMBLY)); 537 PetscCall(MatAssemblyEnd(T,MAT_FINAL_ASSEMBLY)); 538 PetscCall(MatMatTransposeMult(T,fetidpmat_ctx->B_delta,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&BD1)); 539 PetscCall(MatMatMult(fetidpmat_ctx->B_delta,BD1,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&BD2)); 540 PetscCall(MatDestroy(&T)); 541 PetscCall(PetscFree3(W,pivots,Bwork)); 542 if (sub_schurs->is_Ej_all) { 543 PetscCall(ISRestoreIndices(sub_schurs->is_Ej_all,&idxs)); 544 } 545 } 546 } 547 PetscCall(PetscFree(scaling_factors)); 548 PetscCall(PetscFree(cols_B_delta)); 549 550 /* Layout of multipliers */ 551 PetscCall(PetscLayoutCreate(comm,&llay)); 552 PetscCall(PetscLayoutSetBlockSize(llay,1)); 553 PetscCall(PetscLayoutSetSize(llay,fetidpmat_ctx->n_lambda)); 554 PetscCall(PetscLayoutSetUp(llay)); 555 PetscCall(PetscLayoutGetLocalSize(llay,&fetidpmat_ctx->n)); 556 557 /* Local work vector of multipliers */ 558 PetscCall(VecCreate(PETSC_COMM_SELF,&fetidpmat_ctx->lambda_local)); 559 PetscCall(VecSetSizes(fetidpmat_ctx->lambda_local,n_local_lambda,n_local_lambda)); 560 PetscCall(VecSetType(fetidpmat_ctx->lambda_local,VECSEQ)); 561 562 if (BD2) { 563 ISLocalToGlobalMapping l2g; 564 Mat T,TA,*pT; 565 IS is; 566 PetscInt nl,N; 567 BDdelta_DN ctx; 568 569 PetscCall(PetscLayoutGetLocalSize(llay,&nl)); 570 PetscCall(PetscLayoutGetSize(llay,&N)); 571 PetscCall(MatCreate(comm,&T)); 572 PetscCall(MatSetSizes(T,nl,nl,N,N)); 573 PetscCall(MatSetType(T,MATIS)); 574 PetscCall(ISLocalToGlobalMappingCreate(comm,1,n_local_lambda,l2g_indices,PETSC_COPY_VALUES,&l2g)); 575 PetscCall(MatSetLocalToGlobalMapping(T,l2g,l2g)); 576 PetscCall(ISLocalToGlobalMappingDestroy(&l2g)); 577 PetscCall(MatISSetLocalMat(T,BD2)); 578 PetscCall(MatAssemblyBegin(T,MAT_FINAL_ASSEMBLY)); 579 PetscCall(MatAssemblyEnd(T,MAT_FINAL_ASSEMBLY)); 580 PetscCall(MatDestroy(&BD2)); 581 PetscCall(MatConvert(T,MATAIJ,MAT_INITIAL_MATRIX,&TA)); 582 PetscCall(MatDestroy(&T)); 583 PetscCall(ISCreateGeneral(comm,n_local_lambda,l2g_indices,PETSC_USE_POINTER,&is)); 584 PetscCall(MatCreateSubMatrices(TA,1,&is,&is,MAT_INITIAL_MATRIX,&pT)); 585 PetscCall(MatDestroy(&TA)); 586 PetscCall(ISDestroy(&is)); 587 BD2 = pT[0]; 588 PetscCall(PetscFree(pT)); 589 590 /* B_Ddelta for non-redundant multipliers with deluxe scaling */ 591 PetscCall(PetscNew(&ctx)); 592 PetscCall(MatSetType(fetidpmat_ctx->B_Ddelta,MATSHELL)); 593 PetscCall(MatShellSetContext(fetidpmat_ctx->B_Ddelta,ctx)); 594 PetscCall(MatShellSetOperation(fetidpmat_ctx->B_Ddelta,MATOP_MULT,(void (*)(void))MatMult_BDdelta_deluxe_nonred)); 595 PetscCall(MatShellSetOperation(fetidpmat_ctx->B_Ddelta,MATOP_MULT_TRANSPOSE,(void (*)(void))MatMultTranspose_BDdelta_deluxe_nonred)); 596 PetscCall(MatShellSetOperation(fetidpmat_ctx->B_Ddelta,MATOP_DESTROY,(void (*)(void))MatDestroy_BDdelta_deluxe_nonred)); 597 PetscCall(MatSetUp(fetidpmat_ctx->B_Ddelta)); 598 599 PetscCall(PetscObjectReference((PetscObject)BD1)); 600 ctx->BD = BD1; 601 PetscCall(KSPCreate(PETSC_COMM_SELF,&ctx->kBD)); 602 PetscCall(KSPSetOperators(ctx->kBD,BD2,BD2)); 603 PetscCall(VecDuplicate(fetidpmat_ctx->lambda_local,&ctx->work)); 604 fetidpmat_ctx->deluxe_nonred = PETSC_TRUE; 605 } 606 PetscCall(MatDestroy(&BD1)); 607 PetscCall(MatDestroy(&BD2)); 608 609 /* fetidpmat sizes */ 610 fetidpmat_ctx->n += nPgl; 611 fetidpmat_ctx->N = fetidpmat_ctx->n_lambda+nPg; 612 613 /* Global vector for FETI-DP linear system */ 614 PetscCall(VecCreate(comm,&fetidp_global)); 615 PetscCall(VecSetSizes(fetidp_global,fetidpmat_ctx->n,fetidpmat_ctx->N)); 616 PetscCall(VecSetType(fetidp_global,VECMPI)); 617 PetscCall(VecSetUp(fetidp_global)); 618 619 /* Decide layout for fetidp dofs: if it is a saddle point problem 620 pressure is ordered first in the local part of the global vector 621 of the FETI-DP linear system */ 622 if (nPg) { 623 Vec v; 624 IS IS_l2g_p,ais; 625 PetscLayout alay; 626 const PetscInt *idxs,*pranges,*aranges,*lranges; 627 PetscInt *l2g_indices_p,rst; 628 PetscMPIInt rank; 629 630 PetscCall(PetscMalloc1(nPl,&l2g_indices_p)); 631 PetscCall(VecGetLayout(fetidp_global,&alay)); 632 PetscCall(PetscLayoutGetRanges(alay,&aranges)); 633 PetscCall(PetscLayoutGetRanges(play,&pranges)); 634 PetscCall(PetscLayoutGetRanges(llay,&lranges)); 635 636 PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)fetidp_global),&rank)); 637 PetscCall(ISCreateStride(PetscObjectComm((PetscObject)fetidp_global),pranges[rank+1]-pranges[rank],aranges[rank],1,&fetidpmat_ctx->pressure)); 638 PetscCall(PetscObjectSetName((PetscObject)fetidpmat_ctx->pressure,"F_P")); 639 PetscCall(ISCreateStride(PetscObjectComm((PetscObject)fetidp_global),lranges[rank+1]-lranges[rank],aranges[rank]+pranges[rank+1]-pranges[rank],1,&fetidpmat_ctx->lagrange)); 640 PetscCall(PetscObjectSetName((PetscObject)fetidpmat_ctx->lagrange,"F_L")); 641 PetscCall(ISLocalToGlobalMappingGetIndices(l2gmap_p,&idxs)); 642 /* shift local to global indices for pressure */ 643 for (i=0;i<nPl;i++) { 644 PetscMPIInt owner; 645 646 PetscCall(PetscLayoutFindOwner(play,idxs[i],&owner)); 647 l2g_indices_p[i] = idxs[i]-pranges[owner]+aranges[owner]; 648 } 649 PetscCall(ISLocalToGlobalMappingRestoreIndices(l2gmap_p,&idxs)); 650 PetscCall(ISCreateGeneral(comm,nPl,l2g_indices_p,PETSC_OWN_POINTER,&IS_l2g_p)); 651 652 /* local to global scatter for pressure */ 653 PetscCall(VecScatterCreate(fetidpmat_ctx->vP,NULL,fetidp_global,IS_l2g_p,&fetidpmat_ctx->l2g_p)); 654 PetscCall(ISDestroy(&IS_l2g_p)); 655 656 /* scatter for lagrange multipliers only */ 657 PetscCall(VecCreate(comm,&v)); 658 PetscCall(VecSetType(v,VECSTANDARD)); 659 PetscCall(VecSetLayout(v,llay)); 660 PetscCall(VecSetUp(v)); 661 PetscCall(ISCreateGeneral(comm,n_local_lambda,l2g_indices,PETSC_COPY_VALUES,&ais)); 662 PetscCall(VecScatterCreate(fetidpmat_ctx->lambda_local,NULL,v,ais,&fetidpmat_ctx->l2g_lambda_only)); 663 PetscCall(ISDestroy(&ais)); 664 PetscCall(VecDestroy(&v)); 665 666 /* shift local to global indices for multipliers */ 667 for (i=0;i<n_local_lambda;i++) { 668 PetscInt ps; 669 PetscMPIInt owner; 670 671 PetscCall(PetscLayoutFindOwner(llay,l2g_indices[i],&owner)); 672 ps = pranges[owner+1]-pranges[owner]; 673 l2g_indices[i] = l2g_indices[i]-lranges[owner]+aranges[owner]+ps; 674 } 675 676 /* scatter from alldofs to pressures global fetidp vector */ 677 PetscCall(PetscLayoutGetRange(alay,&rst,NULL)); 678 PetscCall(ISCreateStride(comm,nPgl,rst,1,&ais)); 679 PetscCall(VecScatterCreate(pcis->vec1_global,pP,fetidp_global,ais,&fetidpmat_ctx->g2g_p)); 680 PetscCall(ISDestroy(&ais)); 681 } 682 PetscCall(PetscLayoutDestroy(&llay)); 683 PetscCall(PetscLayoutDestroy(&play)); 684 PetscCall(ISCreateGeneral(comm,n_local_lambda,l2g_indices,PETSC_OWN_POINTER,&IS_l2g_lambda)); 685 686 /* scatter from local to global multipliers */ 687 PetscCall(VecScatterCreate(fetidpmat_ctx->lambda_local,NULL,fetidp_global,IS_l2g_lambda,&fetidpmat_ctx->l2g_lambda)); 688 PetscCall(ISDestroy(&IS_l2g_lambda)); 689 PetscCall(ISLocalToGlobalMappingDestroy(&l2gmap_p)); 690 PetscCall(VecDestroy(&fetidp_global)); 691 692 /* Create some work vectors needed by fetidp */ 693 PetscCall(VecDuplicate(pcis->vec1_B,&fetidpmat_ctx->temp_solution_B)); 694 PetscCall(VecDuplicate(pcis->vec1_D,&fetidpmat_ctx->temp_solution_D)); 695 PetscFunctionReturn(0); 696 } 697 698 PetscErrorCode PCBDDCSetupFETIDPPCContext(Mat fetimat, FETIDPPC_ctx fetidppc_ctx) 699 { 700 FETIDPMat_ctx mat_ctx; 701 PC_BDDC *pcbddc = (PC_BDDC*)fetidppc_ctx->pc->data; 702 PC_IS *pcis = (PC_IS*)fetidppc_ctx->pc->data; 703 PetscBool lumped = PETSC_FALSE; 704 705 PetscFunctionBegin; 706 PetscCall(MatShellGetContext(fetimat,&mat_ctx)); 707 /* get references from objects created when setting up feti mat context */ 708 PetscCall(PetscObjectReference((PetscObject)mat_ctx->lambda_local)); 709 fetidppc_ctx->lambda_local = mat_ctx->lambda_local; 710 PetscCall(PetscObjectReference((PetscObject)mat_ctx->B_Ddelta)); 711 fetidppc_ctx->B_Ddelta = mat_ctx->B_Ddelta; 712 if (mat_ctx->deluxe_nonred) { 713 PC pc,mpc; 714 BDdelta_DN ctx; 715 MatSolverType solver; 716 const char *prefix; 717 718 PetscCall(MatShellGetContext(mat_ctx->B_Ddelta,&ctx)); 719 PetscCall(KSPSetType(ctx->kBD,KSPPREONLY)); 720 PetscCall(KSPGetPC(ctx->kBD,&mpc)); 721 PetscCall(KSPGetPC(pcbddc->ksp_D,&pc)); 722 PetscCall(PCSetType(mpc,PCLU)); 723 PetscCall(PCFactorGetMatSolverType(pc,(MatSolverType*)&solver)); 724 if (solver) { 725 PetscCall(PCFactorSetMatSolverType(mpc,solver)); 726 } 727 PetscCall(MatGetOptionsPrefix(fetimat,&prefix)); 728 PetscCall(KSPSetOptionsPrefix(ctx->kBD,prefix)); 729 PetscCall(KSPAppendOptionsPrefix(ctx->kBD,"bddelta_")); 730 PetscCall(KSPSetFromOptions(ctx->kBD)); 731 } 732 733 if (mat_ctx->l2g_lambda_only) { 734 PetscCall(PetscObjectReference((PetscObject)mat_ctx->l2g_lambda_only)); 735 fetidppc_ctx->l2g_lambda = mat_ctx->l2g_lambda_only; 736 } else { 737 PetscCall(PetscObjectReference((PetscObject)mat_ctx->l2g_lambda)); 738 fetidppc_ctx->l2g_lambda = mat_ctx->l2g_lambda; 739 } 740 /* Dirichlet preconditioner */ 741 PetscCall(PetscOptionsGetBool(NULL,((PetscObject)fetimat)->prefix,"-pc_lumped",&lumped,NULL)); 742 if (!lumped) { 743 IS iV; 744 PetscBool discrete_harmonic = PETSC_FALSE; 745 746 PetscCall(PetscObjectQuery((PetscObject)fetidppc_ctx->pc,"__KSPFETIDP_iV",(PetscObject*)&iV)); 747 if (iV) { 748 PetscCall(PetscOptionsGetBool(NULL,((PetscObject)fetimat)->prefix,"-pc_discrete_harmonic",&discrete_harmonic,NULL)); 749 } 750 if (discrete_harmonic) { 751 KSP sksp; 752 PC pc; 753 PCBDDCSubSchurs sub_schurs = pcbddc->sub_schurs; 754 Mat A_II,A_IB,A_BI; 755 IS iP = NULL; 756 PetscBool isshell,reuse = PETSC_FALSE; 757 KSPType ksptype; 758 const char *prefix; 759 760 /* 761 We constructs a Schur complement for 762 763 | A_II A_ID | 764 | A_DI A_DD | 765 766 instead of 767 768 | A_II B^t_II A_ID | 769 | B_II -C_II B_ID | 770 | A_DI B^t_ID A_DD | 771 772 */ 773 if (sub_schurs && sub_schurs->reuse_solver) { 774 PetscCall(PetscObjectQuery((PetscObject)sub_schurs->A,"__KSPFETIDP_iP",(PetscObject*)&iP)); 775 if (iP) reuse = PETSC_TRUE; 776 } 777 if (!reuse) { 778 IS aB; 779 PetscInt nb; 780 PetscCall(ISGetLocalSize(pcis->is_B_local,&nb)); 781 PetscCall(ISCreateStride(PetscObjectComm((PetscObject)pcis->A_II),nb,0,1,&aB)); 782 PetscCall(MatCreateSubMatrix(pcis->A_II,iV,iV,MAT_INITIAL_MATRIX,&A_II)); 783 PetscCall(MatCreateSubMatrix(pcis->A_IB,iV,aB,MAT_INITIAL_MATRIX,&A_IB)); 784 PetscCall(MatCreateSubMatrix(pcis->A_BI,aB,iV,MAT_INITIAL_MATRIX,&A_BI)); 785 PetscCall(ISDestroy(&aB)); 786 } else { 787 PetscCall(MatCreateSubMatrix(sub_schurs->A,pcis->is_I_local,pcis->is_B_local,MAT_INITIAL_MATRIX,&A_IB)); 788 PetscCall(MatCreateSubMatrix(sub_schurs->A,pcis->is_B_local,pcis->is_I_local,MAT_INITIAL_MATRIX,&A_BI)); 789 PetscCall(PetscObjectReference((PetscObject)pcis->A_II)); 790 A_II = pcis->A_II; 791 } 792 PetscCall(MatCreateSchurComplement(A_II,A_II,A_IB,A_BI,pcis->A_BB,&fetidppc_ctx->S_j)); 793 794 /* propagate settings of solver */ 795 PetscCall(MatSchurComplementGetKSP(fetidppc_ctx->S_j,&sksp)); 796 PetscCall(KSPGetType(pcis->ksp_D,&ksptype)); 797 PetscCall(KSPSetType(sksp,ksptype)); 798 PetscCall(KSPGetPC(pcis->ksp_D,&pc)); 799 PetscCall(PetscObjectTypeCompare((PetscObject)pc,PCSHELL,&isshell)); 800 if (!isshell) { 801 MatSolverType solver; 802 PCType pctype; 803 804 PetscCall(PCGetType(pc,&pctype)); 805 PetscCall(PCFactorGetMatSolverType(pc,(MatSolverType*)&solver)); 806 PetscCall(KSPGetPC(sksp,&pc)); 807 PetscCall(PCSetType(pc,pctype)); 808 if (solver) { 809 PetscCall(PCFactorSetMatSolverType(pc,solver)); 810 } 811 } else { 812 PetscCall(KSPGetPC(sksp,&pc)); 813 PetscCall(PCSetType(pc,PCLU)); 814 } 815 PetscCall(MatDestroy(&A_II)); 816 PetscCall(MatDestroy(&A_IB)); 817 PetscCall(MatDestroy(&A_BI)); 818 PetscCall(MatGetOptionsPrefix(fetimat,&prefix)); 819 PetscCall(KSPSetOptionsPrefix(sksp,prefix)); 820 PetscCall(KSPAppendOptionsPrefix(sksp,"harmonic_")); 821 PetscCall(KSPSetFromOptions(sksp)); 822 if (reuse) { 823 PetscCall(KSPSetPC(sksp,sub_schurs->reuse_solver->interior_solver)); 824 PetscCall(PetscObjectIncrementTabLevel((PetscObject)sub_schurs->reuse_solver->interior_solver,(PetscObject)sksp,0)); 825 } 826 } else { /* default Dirichlet preconditioner is pde-harmonic */ 827 PetscCall(MatCreateSchurComplement(pcis->A_II,pcis->A_II,pcis->A_IB,pcis->A_BI,pcis->A_BB,&fetidppc_ctx->S_j)); 828 PetscCall(MatSchurComplementSetKSP(fetidppc_ctx->S_j,pcis->ksp_D)); 829 } 830 } else { 831 PetscCall(PetscObjectReference((PetscObject)pcis->A_BB)); 832 fetidppc_ctx->S_j = pcis->A_BB; 833 } 834 /* saddle-point */ 835 if (mat_ctx->xPg) { 836 PetscCall(PetscObjectReference((PetscObject)mat_ctx->xPg)); 837 fetidppc_ctx->xPg = mat_ctx->xPg; 838 PetscCall(PetscObjectReference((PetscObject)mat_ctx->yPg)); 839 fetidppc_ctx->yPg = mat_ctx->yPg; 840 } 841 PetscFunctionReturn(0); 842 } 843 844 PetscErrorCode FETIDPMatMult_Kernel(Mat fetimat, Vec x, Vec y, PetscBool trans) 845 { 846 FETIDPMat_ctx mat_ctx; 847 PC_BDDC *pcbddc; 848 PC_IS *pcis; 849 850 PetscFunctionBegin; 851 PetscCall(MatShellGetContext(fetimat,&mat_ctx)); 852 pcis = (PC_IS*)mat_ctx->pc->data; 853 pcbddc = (PC_BDDC*)mat_ctx->pc->data; 854 /* Application of B_delta^T */ 855 PetscCall(VecSet(pcis->vec1_B,0.)); 856 PetscCall(VecScatterBegin(mat_ctx->l2g_lambda,x,mat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE)); 857 PetscCall(VecScatterEnd(mat_ctx->l2g_lambda,x,mat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE)); 858 PetscCall(MatMultTranspose(mat_ctx->B_delta,mat_ctx->lambda_local,pcis->vec1_B)); 859 860 /* Add contribution from saddle point */ 861 if (mat_ctx->l2g_p) { 862 PetscCall(VecScatterBegin(mat_ctx->l2g_p,x,mat_ctx->vP,INSERT_VALUES,SCATTER_REVERSE)); 863 PetscCall(VecScatterEnd(mat_ctx->l2g_p,x,mat_ctx->vP,INSERT_VALUES,SCATTER_REVERSE)); 864 if (pcbddc->switch_static) { 865 if (trans) { 866 PetscCall(MatMultTranspose(mat_ctx->B_BI,mat_ctx->vP,pcis->vec1_D)); 867 } else { 868 PetscCall(MatMult(mat_ctx->Bt_BI,mat_ctx->vP,pcis->vec1_D)); 869 } 870 } 871 if (trans) { 872 PetscCall(MatMultTransposeAdd(mat_ctx->B_BB,mat_ctx->vP,pcis->vec1_B,pcis->vec1_B)); 873 } else { 874 PetscCall(MatMultAdd(mat_ctx->Bt_BB,mat_ctx->vP,pcis->vec1_B,pcis->vec1_B)); 875 } 876 } else { 877 if (pcbddc->switch_static) { 878 PetscCall(VecSet(pcis->vec1_D,0.0)); 879 } 880 } 881 /* Application of \widetilde{S}^-1 */ 882 PetscCall(PetscArrayzero(pcbddc->benign_p0,pcbddc->benign_n)); 883 PetscCall(PCBDDCApplyInterfacePreconditioner(mat_ctx->pc,trans)); 884 PetscCall(PetscArrayzero(pcbddc->benign_p0,pcbddc->benign_n)); 885 PetscCall(VecSet(y,0.0)); 886 /* Application of B_delta */ 887 PetscCall(MatMult(mat_ctx->B_delta,pcis->vec1_B,mat_ctx->lambda_local)); 888 /* Contribution from boundary pressures */ 889 if (mat_ctx->C) { 890 const PetscScalar *lx; 891 PetscScalar *ly; 892 893 /* pressure ordered first in the local part of x and y */ 894 PetscCall(VecGetArrayRead(x,&lx)); 895 PetscCall(VecGetArray(y,&ly)); 896 PetscCall(VecPlaceArray(mat_ctx->xPg,lx)); 897 PetscCall(VecPlaceArray(mat_ctx->yPg,ly)); 898 if (trans) { 899 PetscCall(MatMultTranspose(mat_ctx->C,mat_ctx->xPg,mat_ctx->yPg)); 900 } else { 901 PetscCall(MatMult(mat_ctx->C,mat_ctx->xPg,mat_ctx->yPg)); 902 } 903 PetscCall(VecResetArray(mat_ctx->xPg)); 904 PetscCall(VecResetArray(mat_ctx->yPg)); 905 PetscCall(VecRestoreArrayRead(x,&lx)); 906 PetscCall(VecRestoreArray(y,&ly)); 907 } 908 /* Add contribution from saddle point */ 909 if (mat_ctx->l2g_p) { 910 if (trans) { 911 PetscCall(MatMultTranspose(mat_ctx->Bt_BB,pcis->vec1_B,mat_ctx->vP)); 912 } else { 913 PetscCall(MatMult(mat_ctx->B_BB,pcis->vec1_B,mat_ctx->vP)); 914 } 915 if (pcbddc->switch_static) { 916 if (trans) { 917 PetscCall(MatMultTransposeAdd(mat_ctx->Bt_BI,pcis->vec1_D,mat_ctx->vP,mat_ctx->vP)); 918 } else { 919 PetscCall(MatMultAdd(mat_ctx->B_BI,pcis->vec1_D,mat_ctx->vP,mat_ctx->vP)); 920 } 921 } 922 PetscCall(VecScatterBegin(mat_ctx->l2g_p,mat_ctx->vP,y,ADD_VALUES,SCATTER_FORWARD)); 923 PetscCall(VecScatterEnd(mat_ctx->l2g_p,mat_ctx->vP,y,ADD_VALUES,SCATTER_FORWARD)); 924 } 925 PetscCall(VecScatterBegin(mat_ctx->l2g_lambda,mat_ctx->lambda_local,y,ADD_VALUES,SCATTER_FORWARD)); 926 PetscCall(VecScatterEnd(mat_ctx->l2g_lambda,mat_ctx->lambda_local,y,ADD_VALUES,SCATTER_FORWARD)); 927 PetscFunctionReturn(0); 928 } 929 930 PetscErrorCode FETIDPMatMult(Mat fetimat, Vec x, Vec y) 931 { 932 PetscFunctionBegin; 933 PetscCall(FETIDPMatMult_Kernel(fetimat,x,y,PETSC_FALSE)); 934 PetscFunctionReturn(0); 935 } 936 937 PetscErrorCode FETIDPMatMultTranspose(Mat fetimat, Vec x, Vec y) 938 { 939 PetscFunctionBegin; 940 PetscCall(FETIDPMatMult_Kernel(fetimat,x,y,PETSC_TRUE)); 941 PetscFunctionReturn(0); 942 } 943 944 PetscErrorCode FETIDPPCApply_Kernel(PC fetipc, Vec x, Vec y, PetscBool trans) 945 { 946 FETIDPPC_ctx pc_ctx; 947 PC_IS *pcis; 948 949 PetscFunctionBegin; 950 PetscCall(PCShellGetContext(fetipc,&pc_ctx)); 951 pcis = (PC_IS*)pc_ctx->pc->data; 952 /* Application of B_Ddelta^T */ 953 PetscCall(VecScatterBegin(pc_ctx->l2g_lambda,x,pc_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE)); 954 PetscCall(VecScatterEnd(pc_ctx->l2g_lambda,x,pc_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE)); 955 PetscCall(VecSet(pcis->vec2_B,0.0)); 956 PetscCall(MatMultTranspose(pc_ctx->B_Ddelta,pc_ctx->lambda_local,pcis->vec2_B)); 957 /* Application of local Schur complement */ 958 if (trans) { 959 PetscCall(MatMultTranspose(pc_ctx->S_j,pcis->vec2_B,pcis->vec1_B)); 960 } else { 961 PetscCall(MatMult(pc_ctx->S_j,pcis->vec2_B,pcis->vec1_B)); 962 } 963 /* Application of B_Ddelta */ 964 PetscCall(MatMult(pc_ctx->B_Ddelta,pcis->vec1_B,pc_ctx->lambda_local)); 965 PetscCall(VecSet(y,0.0)); 966 PetscCall(VecScatterBegin(pc_ctx->l2g_lambda,pc_ctx->lambda_local,y,ADD_VALUES,SCATTER_FORWARD)); 967 PetscCall(VecScatterEnd(pc_ctx->l2g_lambda,pc_ctx->lambda_local,y,ADD_VALUES,SCATTER_FORWARD)); 968 PetscFunctionReturn(0); 969 } 970 971 PetscErrorCode FETIDPPCApply(PC pc, Vec x, Vec y) 972 { 973 PetscFunctionBegin; 974 PetscCall(FETIDPPCApply_Kernel(pc,x,y,PETSC_FALSE)); 975 PetscFunctionReturn(0); 976 } 977 978 PetscErrorCode FETIDPPCApplyTranspose(PC pc, Vec x, Vec y) 979 { 980 PetscFunctionBegin; 981 PetscCall(FETIDPPCApply_Kernel(pc,x,y,PETSC_TRUE)); 982 PetscFunctionReturn(0); 983 } 984 985 PetscErrorCode FETIDPPCView(PC pc, PetscViewer viewer) 986 { 987 FETIDPPC_ctx pc_ctx; 988 PetscBool iascii; 989 PetscViewer sviewer; 990 991 PetscFunctionBegin; 992 PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii)); 993 if (iascii) { 994 PetscMPIInt rank; 995 PetscBool isschur,isshell; 996 997 PetscCall(PCShellGetContext(pc,&pc_ctx)); 998 PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)pc),&rank)); 999 PetscCall(PetscObjectTypeCompare((PetscObject)pc_ctx->S_j,MATSCHURCOMPLEMENT,&isschur)); 1000 if (isschur) { 1001 PetscCall(PetscViewerASCIIPrintf(viewer," Dirichlet preconditioner (just from rank 0)\n")); 1002 } else { 1003 PetscCall(PetscViewerASCIIPrintf(viewer," Lumped preconditioner (just from rank 0)\n")); 1004 } 1005 PetscCall(PetscViewerGetSubViewer(viewer,PetscObjectComm((PetscObject)pc_ctx->S_j),&sviewer)); 1006 if (rank == 0) { 1007 PetscCall(PetscViewerPushFormat(sviewer,PETSC_VIEWER_ASCII_INFO)); 1008 PetscCall(PetscViewerASCIIPushTab(sviewer)); 1009 PetscCall(MatView(pc_ctx->S_j,sviewer)); 1010 PetscCall(PetscViewerASCIIPopTab(sviewer)); 1011 PetscCall(PetscViewerPopFormat(sviewer)); 1012 } 1013 PetscCall(PetscViewerRestoreSubViewer(viewer,PetscObjectComm((PetscObject)pc_ctx->S_j),&sviewer)); 1014 PetscCall(PetscObjectTypeCompare((PetscObject)pc_ctx->B_Ddelta,MATSHELL,&isshell)); 1015 if (isshell) { 1016 BDdelta_DN ctx; 1017 PetscCall(PetscViewerASCIIPrintf(viewer," FETI-DP BDdelta: DB^t * (B D^-1 B^t)^-1 for deluxe scaling (just from rank 0)\n")); 1018 PetscCall(MatShellGetContext(pc_ctx->B_Ddelta,&ctx)); 1019 PetscCall(PetscViewerGetSubViewer(viewer,PetscObjectComm((PetscObject)pc_ctx->S_j),&sviewer)); 1020 if (rank == 0) { 1021 PetscInt tl; 1022 1023 PetscCall(PetscViewerASCIIGetTab(sviewer,&tl)); 1024 PetscCall(PetscObjectSetTabLevel((PetscObject)ctx->kBD,tl)); 1025 PetscCall(KSPView(ctx->kBD,sviewer)); 1026 PetscCall(PetscViewerPushFormat(sviewer,PETSC_VIEWER_ASCII_INFO)); 1027 PetscCall(MatView(ctx->BD,sviewer)); 1028 PetscCall(PetscViewerPopFormat(sviewer)); 1029 } 1030 PetscCall(PetscViewerRestoreSubViewer(viewer,PetscObjectComm((PetscObject)pc_ctx->S_j),&sviewer)); 1031 } 1032 PetscCall(PetscViewerFlush(viewer)); 1033 } 1034 PetscFunctionReturn(0); 1035 } 1036