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