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