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