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