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