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