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; 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 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 /* For communication of scaling factors */ 115 PetscInt *ptrs_buffer,neigh_position; 116 PetscScalar **all_factors,*send_buffer,*recv_buffer; 117 MPI_Request *send_reqs,*recv_reqs; 118 /* tests */ 119 PetscBool test_fetidp; 120 PetscViewer viewer; 121 /* saddlepoint */ 122 ISLocalToGlobalMapping l2gmap_p; 123 PetscLayout play; 124 IS gP,pP; 125 PetscInt nPl,nPg,nPgl; 126 127 PetscFunctionBegin; 128 ierr = PetscObjectGetComm((PetscObject)(fetidpmat_ctx->pc),&comm);CHKERRQ(ierr); 129 ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 130 ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 131 132 /* saddlepoint */ 133 nPl = 0; 134 nPg = 0; 135 nPgl = 0; 136 gP = NULL; 137 pP = NULL; 138 l2gmap_p = NULL; 139 play = NULL; 140 ierr = PetscObjectQuery((PetscObject)fetidpmat_ctx->pc,"__KSPFETIDP_pP",(PetscObject*)&pP);CHKERRQ(ierr); 141 if (pP) { /* saddle point */ 142 /* subdomain pressures in global numbering */ 143 ierr = PetscObjectQuery((PetscObject)fetidpmat_ctx->pc,"__KSPFETIDP_gP",(PetscObject*)&gP);CHKERRQ(ierr); 144 if (!gP) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"gP not present"); 145 ierr = ISGetLocalSize(gP,&nPl);CHKERRQ(ierr); 146 ierr = VecCreate(PETSC_COMM_SELF,&fetidpmat_ctx->vP);CHKERRQ(ierr); 147 ierr = VecSetSizes(fetidpmat_ctx->vP,nPl,nPl);CHKERRQ(ierr); 148 ierr = VecSetType(fetidpmat_ctx->vP,VECSTANDARD);CHKERRQ(ierr); 149 ierr = VecSetUp(fetidpmat_ctx->vP);CHKERRQ(ierr); 150 151 /* interface pressure matrix */ 152 ierr = PetscObjectQuery((PetscObject)fetidpmat_ctx->pc,"__KSPFETIDP_C",(PetscObject*)&fetidpmat_ctx->C);CHKERRQ(ierr); 153 if (!fetidpmat_ctx->C) { /* null pressure block, compute layout and global numbering for interface pressures */ 154 IS Pg; 155 156 ierr = ISRenumber(gP,NULL,&nPg,&Pg);CHKERRQ(ierr); 157 ierr = ISLocalToGlobalMappingCreateIS(Pg,&l2gmap_p);CHKERRQ(ierr); 158 ierr = ISDestroy(&Pg);CHKERRQ(ierr); 159 ierr = PetscLayoutCreate(comm,&play);CHKERRQ(ierr); 160 ierr = PetscLayoutSetBlockSize(play,1);CHKERRQ(ierr); 161 ierr = PetscLayoutSetSize(play,nPg);CHKERRQ(ierr); 162 ierr = ISGetLocalSize(pP,&nPgl);CHKERRQ(ierr); 163 ierr = PetscLayoutSetLocalSize(play,nPgl);CHKERRQ(ierr); 164 ierr = PetscLayoutSetUp(play);CHKERRQ(ierr); 165 } else { 166 ierr = PetscObjectReference((PetscObject)fetidpmat_ctx->C);CHKERRQ(ierr); 167 ierr = MatGetLocalToGlobalMapping(fetidpmat_ctx->C,&l2gmap_p,NULL);CHKERRQ(ierr); 168 ierr = PetscObjectReference((PetscObject)l2gmap_p);CHKERRQ(ierr); 169 ierr = MatGetSize(fetidpmat_ctx->C,&nPg,NULL);CHKERRQ(ierr); 170 ierr = MatGetLocalSize(fetidpmat_ctx->C,NULL,&nPgl);CHKERRQ(ierr); 171 ierr = MatGetLayouts(fetidpmat_ctx->C,NULL,&llay);CHKERRQ(ierr); 172 ierr = PetscLayoutReference(llay,&play);CHKERRQ(ierr); 173 } 174 ierr = VecCreateMPIWithArray(comm,1,nPgl,nPg,NULL,&fetidpmat_ctx->xPg);CHKERRQ(ierr); 175 ierr = VecCreateMPIWithArray(comm,1,nPgl,nPg,NULL,&fetidpmat_ctx->yPg);CHKERRQ(ierr); 176 177 /* import matrices for interface pressures coupling */ 178 ierr = PetscObjectQuery((PetscObject)fetidpmat_ctx->pc,"__KSPFETIDP_B_BI",(PetscObject*)&fetidpmat_ctx->B_BI);CHKERRQ(ierr); 179 if (!fetidpmat_ctx->B_BI) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"B_BI not present"); 180 ierr = PetscObjectReference((PetscObject)fetidpmat_ctx->B_BI);CHKERRQ(ierr); 181 182 ierr = PetscObjectQuery((PetscObject)fetidpmat_ctx->pc,"__KSPFETIDP_B_BB",(PetscObject*)&fetidpmat_ctx->B_BB);CHKERRQ(ierr); 183 if (!fetidpmat_ctx->B_BB) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"B_BB not present"); 184 ierr = PetscObjectReference((PetscObject)fetidpmat_ctx->B_BB);CHKERRQ(ierr); 185 186 ierr = PetscObjectQuery((PetscObject)fetidpmat_ctx->pc,"__KSPFETIDP_Bt_BI",(PetscObject*)&fetidpmat_ctx->Bt_BI);CHKERRQ(ierr); 187 if (!fetidpmat_ctx->Bt_BI) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Bt_BI not present"); 188 ierr = PetscObjectReference((PetscObject)fetidpmat_ctx->Bt_BI);CHKERRQ(ierr); 189 190 ierr = PetscObjectQuery((PetscObject)fetidpmat_ctx->pc,"__KSPFETIDP_Bt_BB",(PetscObject*)&fetidpmat_ctx->Bt_BB);CHKERRQ(ierr); 191 if (!fetidpmat_ctx->Bt_BB) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Bt_BB not present"); 192 ierr = PetscObjectReference((PetscObject)fetidpmat_ctx->Bt_BB);CHKERRQ(ierr); 193 194 ierr = PetscObjectQuery((PetscObject)fetidpmat_ctx->pc,"__KSPFETIDP_flip" ,(PetscObject*)&fetidpmat_ctx->rhs_flip);CHKERRQ(ierr); 195 if (fetidpmat_ctx->rhs_flip) { 196 ierr = PetscObjectReference((PetscObject)fetidpmat_ctx->rhs_flip);CHKERRQ(ierr); 197 } 198 } 199 200 /* Default type of lagrange multipliers is non-redundant */ 201 fully_redundant = fetidpmat_ctx->fully_redundant; 202 203 /* Evaluate local and global number of lagrange multipliers */ 204 ierr = VecSet(pcis->vec1_N,0.0);CHKERRQ(ierr); 205 n_local_lambda = 0; 206 partial_sum = 0; 207 n_boundary_dofs = 0; 208 s = 0; 209 /* Get Vertices used to define the BDDC */ 210 n_vertices = pcbddc->n_vertices; 211 vertex_indices = pcbddc->local_primal_ref_node; 212 dual_size = pcis->n_B-n_vertices; 213 ierr = PetscMalloc1(dual_size,&dual_dofs_boundary_indices);CHKERRQ(ierr); 214 ierr = PetscMalloc1(dual_size,&aux_local_numbering_1);CHKERRQ(ierr); 215 ierr = PetscMalloc1(dual_size,&aux_local_numbering_2);CHKERRQ(ierr); 216 217 ierr = VecGetArray(pcis->vec1_N,&array);CHKERRQ(ierr); 218 for (i=0;i<pcis->n;i++){ 219 j = mat_graph->count[i]; /* RECALL: mat_graph->count[i] does not count myself */ 220 if (j > 0) n_boundary_dofs++; 221 skip_node = PETSC_FALSE; 222 if (s < n_vertices && vertex_indices[s] == i) { /* it works for a sorted set of vertices */ 223 skip_node = PETSC_TRUE; 224 s++; 225 } 226 if (j < 1) skip_node = PETSC_TRUE; 227 if (mat_graph->special_dof[i] == PCBDDCGRAPH_DIRICHLET_MARK) skip_node = PETSC_TRUE; 228 if (!skip_node) { 229 if (fully_redundant) { 230 /* fully redundant set of lagrange multipliers */ 231 n_lambda_for_dof = (j*(j+1))/2; 232 } else { 233 n_lambda_for_dof = j; 234 } 235 n_local_lambda += j; 236 /* needed to evaluate global number of lagrange multipliers */ 237 array[i]=(1.0*n_lambda_for_dof)/(j+1.0); /* already scaled for the next global sum */ 238 /* store some data needed */ 239 dual_dofs_boundary_indices[partial_sum] = n_boundary_dofs-1; 240 aux_local_numbering_1[partial_sum] = i; 241 aux_local_numbering_2[partial_sum] = n_lambda_for_dof; 242 partial_sum++; 243 } 244 } 245 ierr = VecRestoreArray(pcis->vec1_N,&array);CHKERRQ(ierr); 246 dual_size = partial_sum; 247 248 /* compute global ordering of lagrange multipliers and associate l2g map */ 249 ierr = ISCreateGeneral(comm,partial_sum,aux_local_numbering_1,PETSC_COPY_VALUES,&subset_n);CHKERRQ(ierr); 250 ierr = ISLocalToGlobalMappingApplyIS(pcis->mapping,subset_n,&subset);CHKERRQ(ierr); 251 ierr = ISDestroy(&subset_n);CHKERRQ(ierr); 252 ierr = ISCreateGeneral(comm,partial_sum,aux_local_numbering_2,PETSC_OWN_POINTER,&subset_mult);CHKERRQ(ierr); 253 ierr = ISRenumber(subset,subset_mult,&fetidpmat_ctx->n_lambda,&subset_n);CHKERRQ(ierr); 254 ierr = ISDestroy(&subset);CHKERRQ(ierr); 255 256 #if defined(PETSC_USE_DEBUG) 257 ierr = VecSet(pcis->vec1_global,0.0);CHKERRQ(ierr); 258 ierr = VecScatterBegin(matis->rctx,pcis->vec1_N,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 259 ierr = VecScatterEnd(matis->rctx,pcis->vec1_N,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 260 ierr = VecSum(pcis->vec1_global,&scalar_value);CHKERRQ(ierr); 261 i = (PetscInt)PetscRealPart(scalar_value); 262 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); 263 #endif 264 265 /* init data for scaling factors exchange */ 266 partial_sum = 0; 267 ierr = PetscMalloc1(pcis->n_neigh,&ptrs_buffer);CHKERRQ(ierr); 268 ierr = PetscMalloc1(PetscMax(pcis->n_neigh-1,0),&send_reqs);CHKERRQ(ierr); 269 ierr = PetscMalloc1(PetscMax(pcis->n_neigh-1,0),&recv_reqs);CHKERRQ(ierr); 270 ierr = PetscMalloc1(pcis->n,&all_factors);CHKERRQ(ierr); 271 if (pcis->n_neigh > 0) ptrs_buffer[0]=0; 272 for (i=1;i<pcis->n_neigh;i++) { 273 partial_sum += pcis->n_shared[i]; 274 ptrs_buffer[i] = ptrs_buffer[i-1]+pcis->n_shared[i]; 275 } 276 ierr = PetscMalloc1(partial_sum,&send_buffer);CHKERRQ(ierr); 277 ierr = PetscMalloc1(partial_sum,&recv_buffer);CHKERRQ(ierr); 278 ierr = PetscMalloc1(partial_sum,&all_factors[0]);CHKERRQ(ierr); 279 for (i=0;i<pcis->n-1;i++) { 280 j = mat_graph->count[i]; 281 all_factors[i+1]=all_factors[i]+j; 282 } 283 /* scatter B scaling to N vec */ 284 ierr = VecScatterBegin(pcis->N_to_B,pcis->D,pcis->vec1_N,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 285 ierr = VecScatterEnd(pcis->N_to_B,pcis->D,pcis->vec1_N,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 286 /* communications */ 287 ierr = VecGetArrayRead(pcis->vec1_N,(const PetscScalar**)&array);CHKERRQ(ierr); 288 for (i=1;i<pcis->n_neigh;i++) { 289 for (j=0;j<pcis->n_shared[i];j++) { 290 send_buffer[ptrs_buffer[i-1]+j]=array[pcis->shared[i][j]]; 291 } 292 ierr = PetscMPIIntCast(ptrs_buffer[i]-ptrs_buffer[i-1],&buf_size);CHKERRQ(ierr); 293 ierr = PetscMPIIntCast(pcis->neigh[i],&neigh);CHKERRQ(ierr); 294 ierr = MPI_Isend(&send_buffer[ptrs_buffer[i-1]],buf_size,MPIU_SCALAR,neigh,0,comm,&send_reqs[i-1]);CHKERRQ(ierr); 295 ierr = MPI_Irecv(&recv_buffer[ptrs_buffer[i-1]],buf_size,MPIU_SCALAR,neigh,0,comm,&recv_reqs[i-1]);CHKERRQ(ierr); 296 } 297 ierr = VecRestoreArrayRead(pcis->vec1_N,(const PetscScalar**)&array);CHKERRQ(ierr); 298 if (pcis->n_neigh > 0) { 299 ierr = MPI_Waitall(pcis->n_neigh-1,recv_reqs,MPI_STATUSES_IGNORE);CHKERRQ(ierr); 300 } 301 /* put values in correct places */ 302 for (i=1;i<pcis->n_neigh;i++) { 303 for (j=0;j<pcis->n_shared[i];j++) { 304 k = pcis->shared[i][j]; 305 neigh_position = 0; 306 while(mat_graph->neighbours_set[k][neigh_position] != pcis->neigh[i]) {neigh_position++;} 307 all_factors[k][neigh_position]=recv_buffer[ptrs_buffer[i-1]+j]; 308 } 309 } 310 if (pcis->n_neigh > 0) { 311 ierr = MPI_Waitall(pcis->n_neigh-1,send_reqs,MPI_STATUSES_IGNORE);CHKERRQ(ierr); 312 } 313 ierr = PetscFree(send_reqs);CHKERRQ(ierr); 314 ierr = PetscFree(recv_reqs);CHKERRQ(ierr); 315 ierr = PetscFree(send_buffer);CHKERRQ(ierr); 316 ierr = PetscFree(recv_buffer);CHKERRQ(ierr); 317 ierr = PetscFree(ptrs_buffer);CHKERRQ(ierr); 318 319 /* Compute B and B_delta (local actions) */ 320 ierr = PetscMalloc1(pcis->n_neigh,&aux_sums);CHKERRQ(ierr); 321 ierr = PetscMalloc1(n_local_lambda,&l2g_indices);CHKERRQ(ierr); 322 ierr = PetscMalloc1(n_local_lambda,&vals_B_delta);CHKERRQ(ierr); 323 ierr = PetscMalloc1(n_local_lambda,&cols_B_delta);CHKERRQ(ierr); 324 ierr = PetscMalloc1(n_local_lambda,&scaling_factors);CHKERRQ(ierr); 325 ierr = ISGetIndices(subset_n,&aux_global_numbering);CHKERRQ(ierr); 326 partial_sum=0; 327 cum = 0; 328 for (i=0;i<dual_size;i++) { 329 n_global_lambda = aux_global_numbering[cum]; 330 j = mat_graph->count[aux_local_numbering_1[i]]; 331 aux_sums[0]=0; 332 for (s=1;s<j;s++) { 333 aux_sums[s]=aux_sums[s-1]+j-s+1; 334 } 335 array = all_factors[aux_local_numbering_1[i]]; 336 n_neg_values = 0; 337 while(n_neg_values < j && mat_graph->neighbours_set[aux_local_numbering_1[i]][n_neg_values] < rank) {n_neg_values++;} 338 n_pos_values = j - n_neg_values; 339 if (fully_redundant) { 340 for (s=0;s<n_neg_values;s++) { 341 l2g_indices [partial_sum+s]=aux_sums[s]+n_neg_values-s-1+n_global_lambda; 342 cols_B_delta [partial_sum+s]=dual_dofs_boundary_indices[i]; 343 vals_B_delta [partial_sum+s]=-1.0; 344 scaling_factors[partial_sum+s]=array[s]; 345 } 346 for (s=0;s<n_pos_values;s++) { 347 l2g_indices [partial_sum+s+n_neg_values]=aux_sums[n_neg_values]+s+n_global_lambda; 348 cols_B_delta [partial_sum+s+n_neg_values]=dual_dofs_boundary_indices[i]; 349 vals_B_delta [partial_sum+s+n_neg_values]=1.0; 350 scaling_factors[partial_sum+s+n_neg_values]=array[s+n_neg_values]; 351 } 352 partial_sum += j; 353 } else { 354 /* l2g_indices and default cols and vals of B_delta */ 355 for (s=0;s<j;s++) { 356 l2g_indices [partial_sum+s]=n_global_lambda+s; 357 cols_B_delta [partial_sum+s]=dual_dofs_boundary_indices[i]; 358 vals_B_delta [partial_sum+s]=0.0; 359 } 360 /* B_delta */ 361 if ( n_neg_values > 0 ) { /* there's a rank next to me to the left */ 362 vals_B_delta [partial_sum+n_neg_values-1]=-1.0; 363 } 364 if ( n_neg_values < j ) { /* there's a rank next to me to the right */ 365 vals_B_delta [partial_sum+n_neg_values]=1.0; 366 } 367 /* scaling as in Klawonn-Widlund 1999*/ 368 for (s=0;s<n_neg_values;s++) { 369 scalar_value = 0.0; 370 for (k=0;k<s+1;k++) { 371 scalar_value += array[k]; 372 } 373 scaling_factors[partial_sum+s] = -scalar_value; 374 } 375 for (s=0;s<n_pos_values;s++) { 376 scalar_value = 0.0; 377 for (k=s+n_neg_values;k<j;k++) { 378 scalar_value += array[k]; 379 } 380 scaling_factors[partial_sum+s+n_neg_values] = scalar_value; 381 } 382 partial_sum += j; 383 } 384 cum += aux_local_numbering_2[i]; 385 } 386 ierr = ISRestoreIndices(subset_n,&aux_global_numbering);CHKERRQ(ierr); 387 ierr = ISDestroy(&subset_mult);CHKERRQ(ierr); 388 ierr = ISDestroy(&subset_n);CHKERRQ(ierr); 389 ierr = PetscFree(aux_sums);CHKERRQ(ierr); 390 ierr = PetscFree(aux_local_numbering_1);CHKERRQ(ierr); 391 ierr = PetscFree(dual_dofs_boundary_indices);CHKERRQ(ierr); 392 ierr = PetscFree(all_factors[0]);CHKERRQ(ierr); 393 ierr = PetscFree(all_factors);CHKERRQ(ierr); 394 395 /* Create local part of B_delta */ 396 ierr = MatCreate(PETSC_COMM_SELF,&fetidpmat_ctx->B_delta);CHKERRQ(ierr); 397 ierr = MatSetSizes(fetidpmat_ctx->B_delta,n_local_lambda,pcis->n_B,n_local_lambda,pcis->n_B);CHKERRQ(ierr); 398 ierr = MatSetType(fetidpmat_ctx->B_delta,MATSEQAIJ);CHKERRQ(ierr); 399 ierr = MatSeqAIJSetPreallocation(fetidpmat_ctx->B_delta,1,NULL);CHKERRQ(ierr); 400 ierr = MatSetOption(fetidpmat_ctx->B_delta,MAT_IGNORE_ZERO_ENTRIES,PETSC_TRUE);CHKERRQ(ierr); 401 for (i=0;i<n_local_lambda;i++) { 402 ierr = MatSetValue(fetidpmat_ctx->B_delta,i,cols_B_delta[i],vals_B_delta[i],INSERT_VALUES);CHKERRQ(ierr); 403 } 404 ierr = PetscFree(vals_B_delta);CHKERRQ(ierr); 405 ierr = MatAssemblyBegin(fetidpmat_ctx->B_delta,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 406 ierr = MatAssemblyEnd (fetidpmat_ctx->B_delta,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 407 408 if (fully_redundant) { 409 ierr = MatCreate(PETSC_COMM_SELF,&ScalingMat);CHKERRQ(ierr); 410 ierr = MatSetSizes(ScalingMat,n_local_lambda,n_local_lambda,n_local_lambda,n_local_lambda);CHKERRQ(ierr); 411 ierr = MatSetType(ScalingMat,MATSEQAIJ);CHKERRQ(ierr); 412 ierr = MatSeqAIJSetPreallocation(ScalingMat,1,NULL);CHKERRQ(ierr); 413 for (i=0;i<n_local_lambda;i++) { 414 ierr = MatSetValue(ScalingMat,i,i,scaling_factors[i],INSERT_VALUES);CHKERRQ(ierr); 415 } 416 ierr = MatAssemblyBegin(ScalingMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 417 ierr = MatAssemblyEnd (ScalingMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 418 ierr = MatMatMult(ScalingMat,fetidpmat_ctx->B_delta,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&fetidpmat_ctx->B_Ddelta);CHKERRQ(ierr); 419 ierr = MatDestroy(&ScalingMat);CHKERRQ(ierr); 420 } else { 421 ierr = MatCreate(PETSC_COMM_SELF,&fetidpmat_ctx->B_Ddelta);CHKERRQ(ierr); 422 ierr = MatSetSizes(fetidpmat_ctx->B_Ddelta,n_local_lambda,pcis->n_B,n_local_lambda,pcis->n_B);CHKERRQ(ierr); 423 ierr = MatSetType(fetidpmat_ctx->B_Ddelta,MATSEQAIJ);CHKERRQ(ierr); 424 ierr = MatSeqAIJSetPreallocation(fetidpmat_ctx->B_Ddelta,1,NULL);CHKERRQ(ierr); 425 for (i=0;i<n_local_lambda;i++) { 426 ierr = MatSetValue(fetidpmat_ctx->B_Ddelta,i,cols_B_delta[i],scaling_factors[i],INSERT_VALUES);CHKERRQ(ierr); 427 } 428 ierr = MatAssemblyBegin(fetidpmat_ctx->B_Ddelta,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 429 ierr = MatAssemblyEnd (fetidpmat_ctx->B_Ddelta,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 430 } 431 ierr = PetscFree(scaling_factors);CHKERRQ(ierr); 432 ierr = PetscFree(cols_B_delta);CHKERRQ(ierr); 433 434 /* Layout of multipliers */ 435 ierr = PetscLayoutCreate(comm,&llay);CHKERRQ(ierr); 436 ierr = PetscLayoutSetBlockSize(llay,1);CHKERRQ(ierr); 437 ierr = PetscLayoutSetSize(llay,fetidpmat_ctx->n_lambda);CHKERRQ(ierr); 438 ierr = PetscLayoutSetUp(llay);CHKERRQ(ierr); 439 ierr = PetscLayoutGetLocalSize(llay,&fetidpmat_ctx->n);CHKERRQ(ierr); 440 441 /* fetidpmat sizes */ 442 fetidpmat_ctx->n += nPgl; 443 fetidpmat_ctx->N = fetidpmat_ctx->n_lambda+nPg; 444 445 /* Local work vector of multipliers */ 446 ierr = VecCreate(PETSC_COMM_SELF,&fetidpmat_ctx->lambda_local);CHKERRQ(ierr); 447 ierr = VecSetSizes(fetidpmat_ctx->lambda_local,n_local_lambda,n_local_lambda);CHKERRQ(ierr); 448 ierr = VecSetType(fetidpmat_ctx->lambda_local,VECSEQ);CHKERRQ(ierr); 449 450 /* Global vector for FETI-DP linear system */ 451 ierr = VecCreate(comm,&fetidp_global);CHKERRQ(ierr); 452 ierr = VecSetSizes(fetidp_global,fetidpmat_ctx->n,fetidpmat_ctx->N);CHKERRQ(ierr); 453 ierr = VecSetType(fetidp_global,VECMPI);CHKERRQ(ierr); 454 ierr = VecSetUp(fetidp_global);CHKERRQ(ierr); 455 456 /* Decide layout for fetidp dofs: if it is a saddle point problem 457 pressure is ordered first in the local part of the global vector 458 of the FETI-DP linear system */ 459 if (nPg) { 460 IS IS_l2g_p,ais; 461 PetscLayout alay; 462 const PetscInt *idxs,*pranges,*aranges,*lranges; 463 PetscInt *l2g_indices_p,rst; 464 465 ierr = PetscMalloc1(nPl,&l2g_indices_p);CHKERRQ(ierr); 466 ierr = VecGetLayout(fetidp_global,&alay);CHKERRQ(ierr); 467 ierr = PetscLayoutGetRanges(alay,&aranges);CHKERRQ(ierr); 468 ierr = PetscLayoutGetRanges(play,&pranges);CHKERRQ(ierr); 469 ierr = PetscLayoutGetRanges(llay,&lranges);CHKERRQ(ierr); 470 ierr = ISLocalToGlobalMappingGetIndices(l2gmap_p,&idxs);CHKERRQ(ierr); 471 /* shift local to global indices for pressure */ 472 for (i=0;i<nPl;i++) { 473 PetscInt owner; 474 475 ierr = PetscLayoutFindOwner(play,idxs[i],&owner);CHKERRQ(ierr); 476 l2g_indices_p[i] = idxs[i]-pranges[owner]+aranges[owner]; 477 } 478 ierr = ISLocalToGlobalMappingRestoreIndices(l2gmap_p,&idxs);CHKERRQ(ierr); 479 ierr = ISCreateGeneral(comm,nPl,l2g_indices_p,PETSC_OWN_POINTER,&IS_l2g_p);CHKERRQ(ierr); 480 481 /* local to global scatter for interface pressure */ 482 ierr = VecScatterCreate(fetidpmat_ctx->vP,NULL,fetidp_global,IS_l2g_p,&fetidpmat_ctx->l2g_p);CHKERRQ(ierr); 483 ierr = ISDestroy(&IS_l2g_p);CHKERRQ(ierr); 484 485 /* shift local to global indices for multipliers */ 486 for (i=0;i<n_local_lambda;i++) { 487 PetscInt owner,ps; 488 489 ierr = PetscLayoutFindOwner(llay,l2g_indices[i],&owner);CHKERRQ(ierr); 490 ps = pranges[owner+1]-pranges[owner]; 491 l2g_indices[i] = l2g_indices[i]-lranges[owner]+aranges[owner]+ps; 492 } 493 494 /* scatter from alldofs to interface pressures global fetidp vector */ 495 ierr = PetscLayoutGetRange(alay,&rst,NULL);CHKERRQ(ierr); 496 ierr = ISCreateStride(comm,nPgl,rst,1,&ais);CHKERRQ(ierr); 497 ierr = VecScatterCreate(pcis->vec1_global,pP,fetidp_global,ais,&fetidpmat_ctx->g2g_p);CHKERRQ(ierr); 498 ierr = ISDestroy(&ais);CHKERRQ(ierr); 499 } 500 ierr = PetscLayoutDestroy(&llay);CHKERRQ(ierr); 501 ierr = PetscLayoutDestroy(&play);CHKERRQ(ierr); 502 ierr = ISCreateGeneral(comm,n_local_lambda,l2g_indices,PETSC_OWN_POINTER,&IS_l2g_lambda);CHKERRQ(ierr); 503 /* scatter from local to global multipliers */ 504 ierr = VecScatterCreate(fetidpmat_ctx->lambda_local,NULL,fetidp_global,IS_l2g_lambda,&fetidpmat_ctx->l2g_lambda);CHKERRQ(ierr); 505 ierr = ISDestroy(&IS_l2g_lambda);CHKERRQ(ierr); 506 ierr = ISLocalToGlobalMappingDestroy(&l2gmap_p);CHKERRQ(ierr); 507 508 /* Create some vectors needed by fetidp */ 509 ierr = VecDuplicate(pcis->vec1_B,&fetidpmat_ctx->temp_solution_B);CHKERRQ(ierr); 510 ierr = VecDuplicate(pcis->vec1_D,&fetidpmat_ctx->temp_solution_D);CHKERRQ(ierr); 511 512 test_fetidp = PETSC_FALSE; 513 ierr = PetscOptionsGetBool(NULL,NULL,"-fetidp_check",&test_fetidp,NULL);CHKERRQ(ierr); 514 515 if (test_fetidp && !pcbddc->use_deluxe_scaling) { 516 Vec test_vec,test_vec_p = NULL; 517 IS dirdofs; 518 PetscReal real_value; 519 520 ierr = PetscViewerASCIIGetStdout(comm,&viewer);CHKERRQ(ierr); 521 ierr = PetscViewerASCIIPushSynchronized(viewer);CHKERRQ(ierr); 522 ierr = PetscViewerASCIIPrintf(viewer,"----------FETI_DP SIZES--------------\n");CHKERRQ(ierr); 523 ierr = PetscViewerASCIISynchronizedPrintf(viewer,"[%04d]: Sizes %D %D %D %D %D\n",rank,fetidpmat_ctx->n,fetidpmat_ctx->N,nPgl,nPg,nPl);CHKERRQ(ierr); 524 ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 525 ierr = PetscViewerASCIIPrintf(viewer,"----------FETI_DP TESTS--------------\n");CHKERRQ(ierr); 526 ierr = PetscViewerASCIIPrintf(viewer,"All tests should return zero!\n");CHKERRQ(ierr); 527 ierr = PetscViewerASCIIPrintf(viewer,"FETIDP MAT context in the ");CHKERRQ(ierr); 528 if (fully_redundant) { 529 ierr = PetscViewerASCIIPrintf(viewer,"fully redundant case for lagrange multipliers.\n");CHKERRQ(ierr); 530 } else { 531 ierr = PetscViewerASCIIPrintf(viewer,"Non-fully redundant case for lagrange multiplier.\n");CHKERRQ(ierr); 532 } 533 ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 534 535 /******************************************************************/ 536 /* TEST A/B: Test numbering of global fetidp dofs */ 537 /******************************************************************/ 538 539 ierr = VecDuplicate(fetidpmat_ctx->lambda_local,&test_vec);CHKERRQ(ierr); 540 ierr = VecSet(fetidp_global,1.0);CHKERRQ(ierr); 541 ierr = VecSet(test_vec,1.);CHKERRQ(ierr); 542 ierr = VecScatterBegin(fetidpmat_ctx->l2g_lambda,fetidp_global,fetidpmat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 543 ierr = VecScatterEnd(fetidpmat_ctx->l2g_lambda,fetidp_global,fetidpmat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 544 if (fetidpmat_ctx->l2g_p) { 545 ierr = VecDuplicate(fetidpmat_ctx->vP,&test_vec_p);CHKERRQ(ierr); 546 ierr = VecSet(test_vec_p,1.);CHKERRQ(ierr); 547 ierr = VecScatterBegin(fetidpmat_ctx->l2g_p,fetidp_global,fetidpmat_ctx->vP,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 548 ierr = VecScatterEnd(fetidpmat_ctx->l2g_p,fetidp_global,fetidpmat_ctx->vP,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 549 } 550 scalar_value = -1.0; 551 ierr = VecAXPY(test_vec,scalar_value,fetidpmat_ctx->lambda_local);CHKERRQ(ierr); 552 ierr = VecNorm(test_vec,NORM_INFINITY,&real_value);CHKERRQ(ierr); 553 ierr = VecDestroy(&test_vec);CHKERRQ(ierr); 554 ierr = PetscViewerASCIISynchronizedPrintf(viewer,"A[%04d]: CHECK glob to loc: % 1.14e\n",rank,real_value);CHKERRQ(ierr); 555 if (fetidpmat_ctx->l2g_p) { 556 ierr = VecAXPY(test_vec_p,scalar_value,fetidpmat_ctx->vP);CHKERRQ(ierr); 557 ierr = VecNorm(test_vec_p,NORM_INFINITY,&real_value);CHKERRQ(ierr); 558 ierr = PetscViewerASCIISynchronizedPrintf(viewer,"A[%04d]: CHECK glob to loc (p): % 1.14e\n",rank,real_value);CHKERRQ(ierr); 559 } 560 ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 561 if (fully_redundant) { 562 ierr = VecSet(fetidp_global,0.0);CHKERRQ(ierr); 563 ierr = VecSet(fetidpmat_ctx->lambda_local,0.5);CHKERRQ(ierr); 564 ierr = VecScatterBegin(fetidpmat_ctx->l2g_lambda,fetidpmat_ctx->lambda_local,fetidp_global,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 565 ierr = VecScatterEnd(fetidpmat_ctx->l2g_lambda,fetidpmat_ctx->lambda_local,fetidp_global,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 566 ierr = VecSum(fetidp_global,&scalar_value);CHKERRQ(ierr); 567 ierr = PetscViewerASCIISynchronizedPrintf(viewer,"B[%04d]: CHECK loc to glob: % 1.14e\n",rank,PetscRealPart(scalar_value)-fetidpmat_ctx->n_lambda);CHKERRQ(ierr); 568 ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 569 } 570 if (fetidpmat_ctx->l2g_p) { 571 ierr = VecSet(pcis->vec1_N,1.0);CHKERRQ(ierr); 572 ierr = VecSet(pcis->vec1_global,0.0);CHKERRQ(ierr); 573 ierr = VecScatterBegin(matis->rctx,pcis->vec1_N,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 574 ierr = VecScatterEnd(matis->rctx,pcis->vec1_N,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 575 576 ierr = VecSet(fetidp_global,0.0);CHKERRQ(ierr); 577 ierr = VecSet(fetidpmat_ctx->vP,-1.0);CHKERRQ(ierr); 578 ierr = VecScatterBegin(fetidpmat_ctx->l2g_p,fetidpmat_ctx->vP,fetidp_global,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 579 ierr = VecScatterEnd(fetidpmat_ctx->l2g_p,fetidpmat_ctx->vP,fetidp_global,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 580 ierr = VecScatterBegin(fetidpmat_ctx->g2g_p,fetidp_global,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 581 ierr = VecScatterEnd(fetidpmat_ctx->g2g_p,fetidp_global,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 582 ierr = VecScatterBegin(fetidpmat_ctx->g2g_p,pcis->vec1_global,fetidp_global,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 583 ierr = VecScatterEnd(fetidpmat_ctx->g2g_p,pcis->vec1_global,fetidp_global,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 584 ierr = VecSum(fetidp_global,&scalar_value);CHKERRQ(ierr); 585 ierr = PetscViewerASCIISynchronizedPrintf(viewer,"B[%04d]: CHECK loc to glob (p): % 1.14e\n",rank,PetscRealPart(scalar_value));CHKERRQ(ierr); 586 ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 587 } 588 /******************************************************************/ 589 /* TEST C: It should holds B_delta*w=0, w\in\widehat{W} */ 590 /* This is the meaning of the B matrix */ 591 /******************************************************************/ 592 593 ierr = VecSetRandom(pcis->vec1_N,NULL);CHKERRQ(ierr); 594 ierr = VecSet(pcis->vec1_global,0.0);CHKERRQ(ierr); 595 ierr = VecScatterBegin(matis->rctx,pcis->vec1_N,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 596 ierr = VecScatterEnd(matis->rctx,pcis->vec1_N,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 597 ierr = VecScatterBegin(matis->rctx,pcis->vec1_global,pcis->vec1_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 598 ierr = VecScatterEnd(matis->rctx,pcis->vec1_global,pcis->vec1_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 599 ierr = VecScatterBegin(pcis->N_to_B,pcis->vec1_N,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 600 ierr = VecScatterEnd(pcis->N_to_B,pcis->vec1_N,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 601 /* Action of B_delta */ 602 ierr = MatMult(fetidpmat_ctx->B_delta,pcis->vec1_B,fetidpmat_ctx->lambda_local);CHKERRQ(ierr); 603 ierr = VecSet(fetidp_global,0.0);CHKERRQ(ierr); 604 ierr = VecScatterBegin(fetidpmat_ctx->l2g_lambda,fetidpmat_ctx->lambda_local,fetidp_global,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 605 ierr = VecScatterEnd(fetidpmat_ctx->l2g_lambda,fetidpmat_ctx->lambda_local,fetidp_global,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 606 ierr = VecNorm(fetidp_global,NORM_INFINITY,&real_value);CHKERRQ(ierr); 607 ierr = PetscViewerASCIIPrintf(viewer,"C[coll]: CHECK infty norm of B_delta*w (w continuous): % 1.14e\n",real_value);CHKERRQ(ierr); 608 ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 609 610 /******************************************************************/ 611 /* TEST D: It should holds E_Dw = w - P_Dw w\in\widetilde{W} */ 612 /* E_D = R_D^TR */ 613 /* P_D = B_{D,delta}^T B_{delta} */ 614 /* eq.44 Mandel Tezaur and Dohrmann 2005 */ 615 /******************************************************************/ 616 617 /* compute a random vector in \widetilde{W} */ 618 ierr = VecSetRandom(pcis->vec1_N,NULL);CHKERRQ(ierr); 619 /* set zero at vertices and essential dofs */ 620 scalar_value = 0.0; 621 ierr = VecGetArray(pcis->vec1_N,&array);CHKERRQ(ierr); 622 for (i=0;i<n_vertices;i++) { array[vertex_indices[i]]=scalar_value; } 623 ierr = PCBDDCGraphGetDirichletDofs(pcbddc->mat_graph,&dirdofs);CHKERRQ(ierr); 624 if (dirdofs) { 625 const PetscInt *idxs; 626 PetscInt ndir; 627 628 ierr = ISGetLocalSize(dirdofs,&ndir);CHKERRQ(ierr); 629 ierr = ISGetIndices(dirdofs,&idxs);CHKERRQ(ierr); 630 for (i=0;i<ndir;i++) { array[idxs[i]]=scalar_value; } 631 ierr = ISRestoreIndices(dirdofs,&idxs);CHKERRQ(ierr); 632 } 633 ierr = VecRestoreArray(pcis->vec1_N,&array);CHKERRQ(ierr); 634 /* store w for final comparison */ 635 ierr = VecDuplicate(pcis->vec1_B,&test_vec);CHKERRQ(ierr); 636 ierr = VecScatterBegin(pcis->N_to_B,pcis->vec1_N,test_vec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 637 ierr = VecScatterEnd(pcis->N_to_B,pcis->vec1_N,test_vec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 638 639 /* Jump operator P_D : results stored in pcis->vec1_B */ 640 641 ierr = VecScatterBegin(pcis->N_to_B,pcis->vec1_N,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 642 ierr = VecScatterEnd(pcis->N_to_B,pcis->vec1_N,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 643 /* Action of B_delta */ 644 ierr = MatMult(fetidpmat_ctx->B_delta,pcis->vec1_B,fetidpmat_ctx->lambda_local);CHKERRQ(ierr); 645 ierr = VecSet(fetidp_global,0.0);CHKERRQ(ierr); 646 ierr = VecScatterBegin(fetidpmat_ctx->l2g_lambda,fetidpmat_ctx->lambda_local,fetidp_global,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 647 ierr = VecScatterEnd(fetidpmat_ctx->l2g_lambda,fetidpmat_ctx->lambda_local,fetidp_global,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 648 /* Action of B_Ddelta^T */ 649 ierr = VecScatterBegin(fetidpmat_ctx->l2g_lambda,fetidp_global,fetidpmat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 650 ierr = VecScatterEnd(fetidpmat_ctx->l2g_lambda,fetidp_global,fetidpmat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 651 ierr = MatMultTranspose(fetidpmat_ctx->B_Ddelta,fetidpmat_ctx->lambda_local,pcis->vec1_B);CHKERRQ(ierr); 652 653 /* Average operator E_D : results stored in pcis->vec2_B */ 654 ierr = VecScatterBegin(pcis->N_to_B,pcis->vec1_N,pcis->vec2_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 655 ierr = VecScatterEnd(pcis->N_to_B,pcis->vec1_N,pcis->vec2_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 656 ierr = PCBDDCScalingExtension(fetidpmat_ctx->pc,pcis->vec2_B,pcis->vec1_global);CHKERRQ(ierr); 657 ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_global,pcis->vec2_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 658 ierr = VecScatterEnd(pcis->global_to_B,pcis->vec1_global,pcis->vec2_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 659 660 /* test E_D=I-P_D */ 661 scalar_value = 1.0; 662 ierr = VecAXPY(pcis->vec1_B,scalar_value,pcis->vec2_B);CHKERRQ(ierr); 663 scalar_value = -1.0; 664 ierr = VecAXPY(pcis->vec1_B,scalar_value,test_vec);CHKERRQ(ierr); 665 ierr = VecNorm(pcis->vec1_B,NORM_INFINITY,&real_value);CHKERRQ(ierr); 666 ierr = VecDestroy(&test_vec);CHKERRQ(ierr); 667 ierr = PetscViewerASCIISynchronizedPrintf(viewer,"D[%04d] CHECK infty norm of E_D + P_D - I: % 1.14e\n",rank,real_value);CHKERRQ(ierr); 668 ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 669 670 /******************************************************************/ 671 /* TEST E: It should hold R_D^TP_Dw=0 w\in\widetilde{W} */ 672 /* eq.48 Mandel Tezaur and Dohrmann 2005 */ 673 /******************************************************************/ 674 675 ierr = VecSetRandom(pcis->vec1_N,NULL);CHKERRQ(ierr); 676 /* set zero at vertices and essential dofs */ 677 scalar_value = 0.0; 678 ierr = VecGetArray(pcis->vec1_N,&array);CHKERRQ(ierr); 679 for (i=0;i<n_vertices;i++) { array[vertex_indices[i]]=scalar_value; } 680 if (dirdofs) { 681 const PetscInt *idxs; 682 PetscInt ndir; 683 684 ierr = ISGetLocalSize(dirdofs,&ndir);CHKERRQ(ierr); 685 ierr = ISGetIndices(dirdofs,&idxs);CHKERRQ(ierr); 686 for (i=0;i<ndir;i++) { array[idxs[i]]=scalar_value; } 687 ierr = ISRestoreIndices(dirdofs,&idxs);CHKERRQ(ierr); 688 } 689 ierr = VecRestoreArray(pcis->vec1_N,&array);CHKERRQ(ierr); 690 691 /* Jump operator P_D : results stored in pcis->vec1_B */ 692 693 ierr = VecScatterBegin(pcis->N_to_B,pcis->vec1_N,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 694 ierr = VecScatterEnd(pcis->N_to_B,pcis->vec1_N,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 695 /* Action of B_delta */ 696 ierr = MatMult(fetidpmat_ctx->B_delta,pcis->vec1_B,fetidpmat_ctx->lambda_local);CHKERRQ(ierr); 697 ierr = VecSet(fetidp_global,0.0);CHKERRQ(ierr); 698 ierr = VecScatterBegin(fetidpmat_ctx->l2g_lambda,fetidpmat_ctx->lambda_local,fetidp_global,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 699 ierr = VecScatterEnd(fetidpmat_ctx->l2g_lambda,fetidpmat_ctx->lambda_local,fetidp_global,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 700 /* Action of B_Ddelta^T */ 701 ierr = VecScatterBegin(fetidpmat_ctx->l2g_lambda,fetidp_global,fetidpmat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 702 ierr = VecScatterEnd(fetidpmat_ctx->l2g_lambda,fetidp_global,fetidpmat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 703 ierr = MatMultTranspose(fetidpmat_ctx->B_Ddelta,fetidpmat_ctx->lambda_local,pcis->vec1_B);CHKERRQ(ierr); 704 /* scaling */ 705 ierr = PCBDDCScalingExtension(fetidpmat_ctx->pc,pcis->vec1_B,pcis->vec1_global);CHKERRQ(ierr); 706 ierr = VecNorm(pcis->vec1_global,NORM_INFINITY,&real_value);CHKERRQ(ierr); 707 ierr = PetscViewerASCIIPrintf(viewer,"E[coll]: CHECK infty norm of R^T_D P_D: % 1.14e\n",real_value);CHKERRQ(ierr); 708 ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 709 710 if (!fully_redundant) { 711 /******************************************************************/ 712 /* TEST F: It should holds B_{delta}B^T_{D,delta}=I */ 713 /* Corollary thm 14 Mandel Tezaur and Dohrmann 2005 */ 714 /******************************************************************/ 715 ierr = VecDuplicate(fetidp_global,&test_vec);CHKERRQ(ierr); 716 ierr = VecSetRandom(fetidp_global,NULL);CHKERRQ(ierr); 717 if (fetidpmat_ctx->l2g_p) { 718 ierr = VecSet(fetidpmat_ctx->vP,0.);CHKERRQ(ierr); 719 ierr = VecScatterBegin(fetidpmat_ctx->l2g_p,fetidpmat_ctx->vP,fetidp_global,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 720 ierr = VecScatterEnd(fetidpmat_ctx->l2g_p,fetidpmat_ctx->vP,fetidp_global,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 721 } 722 /* Action of B_Ddelta^T */ 723 ierr = VecScatterBegin(fetidpmat_ctx->l2g_lambda,fetidp_global,fetidpmat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 724 ierr = VecScatterEnd(fetidpmat_ctx->l2g_lambda,fetidp_global,fetidpmat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 725 ierr = MatMultTranspose(fetidpmat_ctx->B_Ddelta,fetidpmat_ctx->lambda_local,pcis->vec1_B);CHKERRQ(ierr); 726 /* Action of B_delta */ 727 ierr = MatMult(fetidpmat_ctx->B_delta,pcis->vec1_B,fetidpmat_ctx->lambda_local);CHKERRQ(ierr); 728 ierr = VecSet(test_vec,0.0);CHKERRQ(ierr); 729 ierr = VecScatterBegin(fetidpmat_ctx->l2g_lambda,fetidpmat_ctx->lambda_local,test_vec,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 730 ierr = VecScatterEnd(fetidpmat_ctx->l2g_lambda,fetidpmat_ctx->lambda_local,test_vec,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 731 scalar_value = -1.0; 732 ierr = VecAXPY(fetidp_global,scalar_value,test_vec);CHKERRQ(ierr); 733 ierr = VecNorm(fetidp_global,NORM_INFINITY,&real_value);CHKERRQ(ierr); 734 ierr = PetscViewerASCIIPrintf(viewer,"E[coll]: CHECK infty norm of P^T_D - I: % 1.14e\n",real_value);CHKERRQ(ierr); 735 ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 736 ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 737 ierr = VecDestroy(&test_vec);CHKERRQ(ierr); 738 } 739 ierr = VecDestroy(&test_vec_p);CHKERRQ(ierr); 740 ierr = ISDestroy(&dirdofs);CHKERRQ(ierr); 741 } 742 /* final cleanup */ 743 ierr = VecDestroy(&fetidp_global);CHKERRQ(ierr); 744 PetscFunctionReturn(0); 745 } 746 747 #undef __FUNCT__ 748 #define __FUNCT__ "PCBDDCSetupFETIDPPCContext" 749 PetscErrorCode PCBDDCSetupFETIDPPCContext(Mat fetimat, FETIDPPC_ctx fetidppc_ctx) 750 { 751 FETIDPMat_ctx mat_ctx; 752 PC_IS *pcis; 753 PetscBool lumped = PETSC_FALSE; 754 PetscErrorCode ierr; 755 756 PetscFunctionBegin; 757 ierr = MatShellGetContext(fetimat,(void**)&mat_ctx);CHKERRQ(ierr); 758 /* get references from objects created when setting up feti mat context */ 759 ierr = PetscObjectReference((PetscObject)mat_ctx->lambda_local);CHKERRQ(ierr); 760 fetidppc_ctx->lambda_local = mat_ctx->lambda_local; 761 ierr = PetscObjectReference((PetscObject)mat_ctx->B_Ddelta);CHKERRQ(ierr); 762 fetidppc_ctx->B_Ddelta = mat_ctx->B_Ddelta; 763 ierr = PetscObjectReference((PetscObject)mat_ctx->l2g_lambda);CHKERRQ(ierr); 764 fetidppc_ctx->l2g_lambda = mat_ctx->l2g_lambda; 765 /* create preconditioner */ 766 pcis = (PC_IS*)fetidppc_ctx->pc->data; 767 /* Dirichlet preconditioner */ 768 ierr = PetscOptionsGetBool(NULL,NULL,"-fetidp_lumped",&lumped,NULL);CHKERRQ(ierr); 769 if (!lumped) { 770 ierr = MatCreateSchurComplement(pcis->A_II,pcis->A_II,pcis->A_IB,pcis->A_BI,pcis->A_BB,&fetidppc_ctx->S_j);CHKERRQ(ierr); 771 ierr = MatSchurComplementSetKSP(fetidppc_ctx->S_j,pcis->ksp_D);CHKERRQ(ierr); 772 } else { 773 ierr = PetscObjectReference((PetscObject)pcis->A_BB);CHKERRQ(ierr); 774 fetidppc_ctx->S_j = pcis->A_BB; 775 } 776 /* saddle-point */ 777 if (mat_ctx->xPg) { 778 ierr = PetscObjectReference((PetscObject)mat_ctx->xPg);CHKERRQ(ierr); 779 fetidppc_ctx->xPg = mat_ctx->xPg; 780 ierr = PetscObjectReference((PetscObject)mat_ctx->yPg);CHKERRQ(ierr); 781 fetidppc_ctx->yPg = mat_ctx->yPg; 782 ierr = PetscObjectQuery((PetscObject)fetidppc_ctx->pc,"__KSPFETIDP_PKSP",(PetscObject*)&fetidppc_ctx->kP);CHKERRQ(ierr); 783 ierr = PetscObjectReference((PetscObject)fetidppc_ctx->kP);CHKERRQ(ierr); 784 } 785 PetscFunctionReturn(0); 786 } 787 788 #undef __FUNCT__ 789 #define __FUNCT__ "FETIDPMatMult" 790 PetscErrorCode FETIDPMatMult(Mat fetimat, Vec x, Vec y) 791 { 792 FETIDPMat_ctx mat_ctx; 793 PC_BDDC *pcbddc; 794 PC_IS *pcis; 795 PetscErrorCode ierr; 796 797 PetscFunctionBegin; 798 ierr = MatShellGetContext(fetimat,(void**)&mat_ctx);CHKERRQ(ierr); 799 pcis = (PC_IS*)mat_ctx->pc->data; 800 pcbddc = (PC_BDDC*)mat_ctx->pc->data; 801 /* Application of B_delta^T */ 802 ierr = VecSet(pcis->vec1_B,0.);CHKERRQ(ierr); 803 ierr = VecScatterBegin(mat_ctx->l2g_lambda,x,mat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 804 ierr = VecScatterEnd(mat_ctx->l2g_lambda,x,mat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 805 ierr = MatMultTranspose(mat_ctx->B_delta,mat_ctx->lambda_local,pcis->vec1_B);CHKERRQ(ierr); 806 807 /* Add contribution from saddle point */ 808 if (mat_ctx->l2g_p) { 809 ierr = VecScatterBegin(mat_ctx->l2g_p,x,mat_ctx->vP,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 810 ierr = VecScatterEnd(mat_ctx->l2g_p,x,mat_ctx->vP,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 811 if (pcbddc->switch_static) { 812 ierr = MatMult(mat_ctx->Bt_BI,mat_ctx->vP,pcis->vec1_D);CHKERRQ(ierr); 813 } 814 ierr = MatMultAdd(mat_ctx->Bt_BB,mat_ctx->vP,pcis->vec1_B,pcis->vec1_B);CHKERRQ(ierr); 815 } else { 816 if (pcbddc->switch_static) { 817 ierr = VecSet(pcis->vec1_D,0.0);CHKERRQ(ierr); 818 } 819 } 820 /* Application of \widetilde{S}^-1 */ 821 ierr = PetscMemzero(pcbddc->benign_p0,pcbddc->benign_n*sizeof(PetscScalar));CHKERRQ(ierr); 822 ierr = PCBDDCApplyInterfacePreconditioner(mat_ctx->pc,PETSC_FALSE);CHKERRQ(ierr); 823 ierr = PetscMemzero(pcbddc->benign_p0,pcbddc->benign_n*sizeof(PetscScalar));CHKERRQ(ierr); 824 ierr = VecSet(y,0.0);CHKERRQ(ierr); 825 /* Application of B_delta */ 826 ierr = MatMult(mat_ctx->B_delta,pcis->vec1_B,mat_ctx->lambda_local);CHKERRQ(ierr); 827 /* Contribution from boundary pressures */ 828 if (mat_ctx->C) { 829 const PetscScalar *lx; 830 PetscScalar *ly; 831 832 /* pressures ordered first in x and y */ 833 ierr = VecGetArrayRead(x,&lx);CHKERRQ(ierr); 834 ierr = VecGetArray(y,&ly);CHKERRQ(ierr); 835 ierr = VecPlaceArray(mat_ctx->xPg,lx);CHKERRQ(ierr); 836 ierr = VecPlaceArray(mat_ctx->yPg,ly);CHKERRQ(ierr); 837 ierr = MatMult(mat_ctx->C,mat_ctx->xPg,mat_ctx->yPg);CHKERRQ(ierr); 838 ierr = VecResetArray(mat_ctx->xPg);CHKERRQ(ierr); 839 ierr = VecResetArray(mat_ctx->yPg);CHKERRQ(ierr); 840 ierr = VecRestoreArrayRead(x,&lx);CHKERRQ(ierr); 841 ierr = VecRestoreArray(y,&ly);CHKERRQ(ierr); 842 } 843 /* Add contribution from saddle point */ 844 if (mat_ctx->l2g_p) { 845 ierr = MatMult(mat_ctx->B_BB,pcis->vec1_B,mat_ctx->vP);CHKERRQ(ierr); 846 if (pcbddc->switch_static) { 847 ierr = MatMultAdd(mat_ctx->B_BI,pcis->vec1_D,mat_ctx->vP,mat_ctx->vP);CHKERRQ(ierr); 848 } 849 ierr = VecScatterBegin(mat_ctx->l2g_p,mat_ctx->vP,y,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 850 ierr = VecScatterEnd(mat_ctx->l2g_p,mat_ctx->vP,y,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 851 } 852 ierr = VecScatterBegin(mat_ctx->l2g_lambda,mat_ctx->lambda_local,y,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 853 ierr = VecScatterEnd(mat_ctx->l2g_lambda,mat_ctx->lambda_local,y,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 854 PetscFunctionReturn(0); 855 } 856 857 #undef __FUNCT__ 858 #define __FUNCT__ "FETIDPMatMultTranspose" 859 PetscErrorCode FETIDPMatMultTranspose(Mat fetimat, Vec x, Vec y) 860 { 861 FETIDPMat_ctx mat_ctx; 862 PC_IS *pcis; 863 PetscErrorCode ierr; 864 865 PetscFunctionBegin; 866 ierr = MatShellGetContext(fetimat,(void**)&mat_ctx);CHKERRQ(ierr); 867 pcis = (PC_IS*)mat_ctx->pc->data; 868 /* Application of B_delta^T */ 869 ierr = VecScatterBegin(mat_ctx->l2g_lambda,x,mat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 870 ierr = VecScatterEnd(mat_ctx->l2g_lambda,x,mat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 871 ierr = MatMultTranspose(mat_ctx->B_delta,mat_ctx->lambda_local,pcis->vec1_B);CHKERRQ(ierr); 872 /* Application of \widetilde{S}^-1 */ 873 ierr = VecSet(pcis->vec1_D,0.0);CHKERRQ(ierr); 874 ierr = PCBDDCApplyInterfacePreconditioner(mat_ctx->pc,PETSC_TRUE);CHKERRQ(ierr); 875 /* Application of B_delta */ 876 ierr = MatMult(mat_ctx->B_delta,pcis->vec1_B,mat_ctx->lambda_local);CHKERRQ(ierr); 877 ierr = VecSet(y,0.0);CHKERRQ(ierr); 878 ierr = VecScatterBegin(mat_ctx->l2g_lambda,mat_ctx->lambda_local,y,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 879 ierr = VecScatterEnd(mat_ctx->l2g_lambda,mat_ctx->lambda_local,y,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 880 PetscFunctionReturn(0); 881 } 882 883 #undef __FUNCT__ 884 #define __FUNCT__ "FETIDPPCApply" 885 PetscErrorCode FETIDPPCApply(PC fetipc, Vec x, Vec y) 886 { 887 FETIDPPC_ctx pc_ctx; 888 PC_IS *pcis; 889 PetscErrorCode ierr; 890 891 PetscFunctionBegin; 892 ierr = PCShellGetContext(fetipc,(void**)&pc_ctx);CHKERRQ(ierr); 893 pcis = (PC_IS*)pc_ctx->pc->data; 894 /* Application of B_Ddelta^T */ 895 ierr = VecScatterBegin(pc_ctx->l2g_lambda,x,pc_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 896 ierr = VecScatterEnd(pc_ctx->l2g_lambda,x,pc_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 897 ierr = VecSet(pcis->vec2_B,0.0);CHKERRQ(ierr); 898 ierr = MatMultTranspose(pc_ctx->B_Ddelta,pc_ctx->lambda_local,pcis->vec2_B);CHKERRQ(ierr); 899 /* Application of local Schur complement */ 900 ierr = MatMult(pc_ctx->S_j,pcis->vec2_B,pcis->vec1_B);CHKERRQ(ierr); 901 /* Application of B_Ddelta */ 902 ierr = MatMult(pc_ctx->B_Ddelta,pcis->vec1_B,pc_ctx->lambda_local);CHKERRQ(ierr); 903 ierr = VecSet(y,0.0);CHKERRQ(ierr); 904 ierr = VecScatterBegin(pc_ctx->l2g_lambda,pc_ctx->lambda_local,y,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 905 ierr = VecScatterEnd(pc_ctx->l2g_lambda,pc_ctx->lambda_local,y,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 906 /* interface pressure preconditioner */ 907 if (pc_ctx->kP) { 908 const PetscScalar *lx; 909 PetscScalar *ly; 910 911 /* pressures ordered first in x and y */ 912 ierr = VecGetArrayRead(x,&lx);CHKERRQ(ierr); 913 ierr = VecGetArray(y,&ly);CHKERRQ(ierr); 914 ierr = VecPlaceArray(pc_ctx->xPg,lx);CHKERRQ(ierr); 915 ierr = VecPlaceArray(pc_ctx->yPg,ly);CHKERRQ(ierr); 916 ierr = KSPSolve(pc_ctx->kP,pc_ctx->xPg,pc_ctx->yPg);CHKERRQ(ierr); 917 ierr = VecResetArray(pc_ctx->xPg);CHKERRQ(ierr); 918 ierr = VecResetArray(pc_ctx->yPg);CHKERRQ(ierr); 919 ierr = VecRestoreArrayRead(x,&lx);CHKERRQ(ierr); 920 ierr = VecRestoreArray(y,&ly);CHKERRQ(ierr); 921 } 922 PetscFunctionReturn(0); 923 } 924 925 #undef __FUNCT__ 926 #define __FUNCT__ "FETIDPPCApplyTranspose" 927 PetscErrorCode FETIDPPCApplyTranspose(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 = MatMultTranspose(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 PetscFunctionReturn(0); 949 } 950