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 = PetscMalloc(sizeof(*newctx),&newctx);CHKERRQ(ierr); 13 newctx->lambda_local = 0; 14 newctx->temp_solution_B = 0; 15 newctx->temp_solution_D = 0; 16 newctx->B_delta = 0; 17 newctx->B_Ddelta = 0; /* theoretically belongs to the FETIDP preconditioner */ 18 newctx->l2g_lambda = 0; 19 /* increase the reference count for BDDC preconditioner */ 20 ierr = PetscObjectReference((PetscObject)pc);CHKERRQ(ierr); 21 newctx->pc = pc; 22 *fetidpmat_ctx = newctx; 23 PetscFunctionReturn(0); 24 } 25 26 #undef __FUNCT__ 27 #define __FUNCT__ "PCBDDCCreateFETIDPPCContext" 28 PetscErrorCode PCBDDCCreateFETIDPPCContext(PC pc, FETIDPPC_ctx *fetidppc_ctx) 29 { 30 FETIDPPC_ctx newctx; 31 PetscErrorCode ierr; 32 33 PetscFunctionBegin; 34 ierr = PetscMalloc(sizeof(*newctx),&newctx);CHKERRQ(ierr); 35 newctx->lambda_local = 0; 36 newctx->B_Ddelta = 0; 37 newctx->l2g_lambda = 0; 38 /* increase the reference count for BDDC preconditioner */ 39 ierr = PetscObjectReference((PetscObject)pc);CHKERRQ(ierr); 40 newctx->pc = pc; 41 *fetidppc_ctx = newctx; 42 PetscFunctionReturn(0); 43 } 44 45 #undef __FUNCT__ 46 #define __FUNCT__ "PCBDDCDestroyFETIDPMat" 47 PetscErrorCode PCBDDCDestroyFETIDPMat(Mat A) 48 { 49 FETIDPMat_ctx mat_ctx; 50 PetscErrorCode ierr; 51 52 PetscFunctionBegin; 53 ierr = MatShellGetContext(A,(void**)&mat_ctx);CHKERRQ(ierr); 54 ierr = VecDestroy(&mat_ctx->lambda_local);CHKERRQ(ierr); 55 ierr = VecDestroy(&mat_ctx->temp_solution_D);CHKERRQ(ierr); 56 ierr = VecDestroy(&mat_ctx->temp_solution_B);CHKERRQ(ierr); 57 ierr = MatDestroy(&mat_ctx->B_delta);CHKERRQ(ierr); 58 ierr = MatDestroy(&mat_ctx->B_Ddelta);CHKERRQ(ierr); 59 ierr = VecScatterDestroy(&mat_ctx->l2g_lambda);CHKERRQ(ierr); 60 ierr = PCDestroy(&mat_ctx->pc);CHKERRQ(ierr); /* decrease PCBDDC reference count */ 61 ierr = PetscFree(mat_ctx);CHKERRQ(ierr); 62 PetscFunctionReturn(0); 63 } 64 65 #undef __FUNCT__ 66 #define __FUNCT__ "PCBDDCDestroyFETIDPPC" 67 PetscErrorCode PCBDDCDestroyFETIDPPC(PC pc) 68 { 69 FETIDPPC_ctx pc_ctx; 70 PetscErrorCode ierr; 71 72 PetscFunctionBegin; 73 ierr = PCShellGetContext(pc,(void**)&pc_ctx);CHKERRQ(ierr); 74 ierr = VecDestroy(&pc_ctx->lambda_local);CHKERRQ(ierr); 75 ierr = MatDestroy(&pc_ctx->B_Ddelta);CHKERRQ(ierr); 76 ierr = VecScatterDestroy(&pc_ctx->l2g_lambda);CHKERRQ(ierr); 77 ierr = PCDestroy(&pc_ctx->pc);CHKERRQ(ierr); /* decrease PCBDDC reference count */ 78 ierr = PetscFree(pc_ctx);CHKERRQ(ierr); 79 PetscFunctionReturn(0); 80 } 81 82 #undef __FUNCT__ 83 #define __FUNCT__ "PCBDDCSetupFETIDPMatContext" 84 PetscErrorCode PCBDDCSetupFETIDPMatContext(FETIDPMat_ctx fetidpmat_ctx ) 85 { 86 PetscErrorCode ierr; 87 PC_IS *pcis=(PC_IS*)fetidpmat_ctx->pc->data; 88 PC_BDDC *pcbddc=(PC_BDDC*)fetidpmat_ctx->pc->data; 89 PCBDDCGraph mat_graph=pcbddc->mat_graph; 90 Mat_IS *matis = (Mat_IS*)fetidpmat_ctx->pc->pmat->data; 91 MPI_Comm comm; 92 Mat ScalingMat; 93 Vec lambda_global; 94 IS IS_l2g_lambda; 95 PetscBool skip_node,fully_redundant; 96 PetscInt i,j,k,s,n_boundary_dofs,n_global_lambda,n_vertices,partial_sum; 97 PetscInt n_local_lambda,n_lambda_for_dof,dual_size,n_neg_values,n_pos_values; 98 PetscMPIInt rank,size,buf_size,neigh; 99 PetscScalar scalar_value; 100 PetscInt *vertex_indices; 101 PetscInt *dual_dofs_boundary_indices,*aux_local_numbering_1,*aux_global_numbering; 102 PetscInt *aux_sums,*cols_B_delta,*l2g_indices; 103 PetscScalar *array,*scaling_factors,*vals_B_delta; 104 PetscInt *aux_local_numbering_2; 105 /* For communication of scaling factors */ 106 PetscInt *ptrs_buffer,neigh_position; 107 PetscScalar **all_factors,*send_buffer,*recv_buffer; 108 MPI_Request *send_reqs,*recv_reqs; 109 /* tests */ 110 Vec test_vec; 111 PetscBool test_fetidp; 112 PetscViewer viewer; 113 114 PetscFunctionBegin; 115 ierr = PetscObjectGetComm((PetscObject)(fetidpmat_ctx->pc),&comm);CHKERRQ(ierr); 116 ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 117 ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 118 119 /* Default type of lagrange multipliers is non-redundant */ 120 fully_redundant = PETSC_FALSE; 121 ierr = PetscOptionsGetBool(NULL,"-fetidp_fullyredundant",&fully_redundant,NULL);CHKERRQ(ierr); 122 123 /* Evaluate local and global number of lagrange multipliers */ 124 ierr = VecSet(pcis->vec1_N,0.0);CHKERRQ(ierr); 125 n_local_lambda = 0; 126 partial_sum = 0; 127 n_boundary_dofs = 0; 128 s = 0; 129 /* Get Vertices used to define the BDDC */ 130 ierr = PCBDDCGetPrimalVerticesLocalIdx(fetidpmat_ctx->pc,&n_vertices,&vertex_indices);CHKERRQ(ierr); 131 dual_size = pcis->n_B-n_vertices; 132 ierr = PetscSortInt(n_vertices,vertex_indices);CHKERRQ(ierr); 133 ierr = PetscMalloc1(dual_size,&dual_dofs_boundary_indices);CHKERRQ(ierr); 134 ierr = PetscMalloc1(dual_size,&aux_local_numbering_1);CHKERRQ(ierr); 135 ierr = PetscMalloc1(dual_size,&aux_local_numbering_2);CHKERRQ(ierr); 136 137 ierr = VecGetArray(pcis->vec1_N,&array);CHKERRQ(ierr); 138 for (i=0;i<pcis->n;i++){ 139 j = mat_graph->count[i]; /* RECALL: mat_graph->count[i] does not count myself */ 140 if ( j > 0 ) { 141 n_boundary_dofs++; 142 } 143 skip_node = PETSC_FALSE; 144 if ( s < n_vertices && vertex_indices[s]==i) { /* it works for a sorted set of vertices */ 145 skip_node = PETSC_TRUE; 146 s++; 147 } 148 if (j < 1) { 149 skip_node = PETSC_TRUE; 150 } 151 if ( !skip_node ) { 152 if (fully_redundant) { 153 /* fully redundant set of lagrange multipliers */ 154 n_lambda_for_dof = (j*(j+1))/2; 155 } else { 156 n_lambda_for_dof = j; 157 } 158 n_local_lambda += j; 159 /* needed to evaluate global number of lagrange multipliers */ 160 array[i]=(1.0*n_lambda_for_dof)/(j+1.0); /* already scaled for the next global sum */ 161 /* store some data needed */ 162 dual_dofs_boundary_indices[partial_sum] = n_boundary_dofs-1; 163 aux_local_numbering_1[partial_sum] = i; 164 aux_local_numbering_2[partial_sum] = n_lambda_for_dof; 165 partial_sum++; 166 } 167 } 168 ierr = VecRestoreArray(pcis->vec1_N,&array);CHKERRQ(ierr); 169 170 ierr = VecSet(pcis->vec1_global,0.0);CHKERRQ(ierr); 171 ierr = VecScatterBegin(matis->ctx,pcis->vec1_N,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 172 ierr = VecScatterEnd (matis->ctx,pcis->vec1_N,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 173 ierr = VecSum(pcis->vec1_global,&scalar_value);CHKERRQ(ierr); 174 fetidpmat_ctx->n_lambda = (PetscInt)PetscRealPart(scalar_value); 175 176 /* compute global ordering of lagrange multipliers and associate l2g map */ 177 ierr = PCBDDCSubsetNumbering(comm,matis->mapping,partial_sum,aux_local_numbering_1,aux_local_numbering_2,&i,&aux_global_numbering);CHKERRQ(ierr); 178 if (i != fetidpmat_ctx->n_lambda) { 179 SETERRQ3(PETSC_COMM_WORLD,PETSC_ERR_PLIB,"Error in %s: global number of multipliers mismatch! (%d!=%d)\n",__FUNCT__,fetidpmat_ctx->n_lambda,i); 180 } 181 ierr = PetscFree(aux_local_numbering_2);CHKERRQ(ierr); 182 183 /* init data for scaling factors exchange */ 184 partial_sum = 0; 185 j = 0; 186 ierr = PetscMalloc1(pcis->n_neigh,&ptrs_buffer);CHKERRQ(ierr); 187 ierr = PetscMalloc1((pcis->n_neigh-1),&send_reqs);CHKERRQ(ierr); 188 ierr = PetscMalloc1((pcis->n_neigh-1),&recv_reqs);CHKERRQ(ierr); 189 ierr = PetscMalloc1(pcis->n,&all_factors);CHKERRQ(ierr); 190 ptrs_buffer[0]=0; 191 for (i=1;i<pcis->n_neigh;i++) { 192 partial_sum += pcis->n_shared[i]; 193 ptrs_buffer[i] = ptrs_buffer[i-1]+pcis->n_shared[i]; 194 } 195 ierr = PetscMalloc1(partial_sum,&send_buffer);CHKERRQ(ierr); 196 ierr = PetscMalloc1(partial_sum,&recv_buffer);CHKERRQ(ierr); 197 ierr = PetscMalloc1(partial_sum,&all_factors[0]);CHKERRQ(ierr); 198 for (i=0;i<pcis->n-1;i++) { 199 j = mat_graph->count[i]; 200 all_factors[i+1]=all_factors[i]+j; 201 } 202 /* scatter B scaling to N vec */ 203 ierr = VecScatterBegin(pcis->N_to_B,pcis->D,pcis->vec1_N,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 204 ierr = VecScatterEnd(pcis->N_to_B,pcis->D,pcis->vec1_N,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 205 /* communications */ 206 ierr = VecGetArray(pcis->vec1_N,&array);CHKERRQ(ierr); 207 for (i=1;i<pcis->n_neigh;i++) { 208 for (j=0;j<pcis->n_shared[i];j++) { 209 send_buffer[ptrs_buffer[i-1]+j]=array[pcis->shared[i][j]]; 210 } 211 ierr = PetscMPIIntCast(ptrs_buffer[i]-ptrs_buffer[i-1],&buf_size);CHKERRQ(ierr); 212 ierr = PetscMPIIntCast(pcis->neigh[i],&neigh);CHKERRQ(ierr); 213 ierr = MPI_Isend(&send_buffer[ptrs_buffer[i-1]],buf_size,MPIU_SCALAR,neigh,0,comm,&send_reqs[i-1]);CHKERRQ(ierr); 214 ierr = MPI_Irecv(&recv_buffer[ptrs_buffer[i-1]],buf_size,MPIU_SCALAR,neigh,0,comm,&recv_reqs[i-1]);CHKERRQ(ierr); 215 } 216 ierr = VecRestoreArray(pcis->vec1_N,&array);CHKERRQ(ierr); 217 ierr = MPI_Waitall((pcis->n_neigh-1),recv_reqs,MPI_STATUSES_IGNORE);CHKERRQ(ierr); 218 /* put values in correct places */ 219 for (i=1;i<pcis->n_neigh;i++) { 220 for (j=0;j<pcis->n_shared[i];j++) { 221 k = pcis->shared[i][j]; 222 neigh_position = 0; 223 while(mat_graph->neighbours_set[k][neigh_position] != pcis->neigh[i]) {neigh_position++;} 224 all_factors[k][neigh_position]=recv_buffer[ptrs_buffer[i-1]+j]; 225 } 226 } 227 ierr = MPI_Waitall((pcis->n_neigh-1),send_reqs,MPI_STATUSES_IGNORE);CHKERRQ(ierr); 228 ierr = PetscFree(send_reqs);CHKERRQ(ierr); 229 ierr = PetscFree(recv_reqs);CHKERRQ(ierr); 230 ierr = PetscFree(send_buffer);CHKERRQ(ierr); 231 ierr = PetscFree(recv_buffer);CHKERRQ(ierr); 232 ierr = PetscFree(ptrs_buffer);CHKERRQ(ierr); 233 234 /* Compute B and B_delta (local actions) */ 235 ierr = PetscMalloc1(pcis->n_neigh,&aux_sums);CHKERRQ(ierr); 236 ierr = PetscMalloc1(n_local_lambda,&l2g_indices);CHKERRQ(ierr); 237 ierr = PetscMalloc1(n_local_lambda,&vals_B_delta);CHKERRQ(ierr); 238 ierr = PetscMalloc1(n_local_lambda,&cols_B_delta);CHKERRQ(ierr); 239 ierr = PetscMalloc1(n_local_lambda,&scaling_factors);CHKERRQ(ierr); 240 n_global_lambda=0; 241 partial_sum=0; 242 for (i=0;i<dual_size;i++) { 243 n_global_lambda = aux_global_numbering[i]; 244 j = mat_graph->count[aux_local_numbering_1[i]]; 245 aux_sums[0]=0; 246 for (s=1;s<j;s++) { 247 aux_sums[s]=aux_sums[s-1]+j-s+1; 248 } 249 array = all_factors[aux_local_numbering_1[i]]; 250 n_neg_values = 0; 251 while(n_neg_values < j && mat_graph->neighbours_set[aux_local_numbering_1[i]][n_neg_values] < rank) {n_neg_values++;} 252 n_pos_values = j - n_neg_values; 253 if (fully_redundant) { 254 for (s=0;s<n_neg_values;s++) { 255 l2g_indices [partial_sum+s]=aux_sums[s]+n_neg_values-s-1+n_global_lambda; 256 cols_B_delta [partial_sum+s]=dual_dofs_boundary_indices[i]; 257 vals_B_delta [partial_sum+s]=-1.0; 258 scaling_factors[partial_sum+s]=array[s]; 259 } 260 for (s=0;s<n_pos_values;s++) { 261 l2g_indices [partial_sum+s+n_neg_values]=aux_sums[n_neg_values]+s+n_global_lambda; 262 cols_B_delta [partial_sum+s+n_neg_values]=dual_dofs_boundary_indices[i]; 263 vals_B_delta [partial_sum+s+n_neg_values]=1.0; 264 scaling_factors[partial_sum+s+n_neg_values]=array[s+n_neg_values]; 265 } 266 partial_sum += j; 267 } else { 268 /* l2g_indices and default cols and vals of B_delta */ 269 for (s=0;s<j;s++) { 270 l2g_indices [partial_sum+s]=n_global_lambda+s; 271 cols_B_delta [partial_sum+s]=dual_dofs_boundary_indices[i]; 272 vals_B_delta [partial_sum+s]=0.0; 273 } 274 /* B_delta */ 275 if ( n_neg_values > 0 ) { /* there's a rank next to me to the left */ 276 vals_B_delta [partial_sum+n_neg_values-1]=-1.0; 277 } 278 if ( n_neg_values < j ) { /* there's a rank next to me to the right */ 279 vals_B_delta [partial_sum+n_neg_values]=1.0; 280 } 281 /* scaling as in Klawonn-Widlund 1999*/ 282 for (s=0;s<n_neg_values;s++) { 283 scalar_value = 0.0; 284 for (k=0;k<s+1;k++) { 285 scalar_value += array[k]; 286 } 287 scaling_factors[partial_sum+s] = -scalar_value; 288 } 289 for (s=0;s<n_pos_values;s++) { 290 scalar_value = 0.0; 291 for (k=s+n_neg_values;k<j;k++) { 292 scalar_value += array[k]; 293 } 294 scaling_factors[partial_sum+s+n_neg_values] = scalar_value; 295 } 296 partial_sum += j; 297 } 298 } 299 ierr = PetscFree(aux_global_numbering);CHKERRQ(ierr); 300 ierr = PetscFree(aux_sums);CHKERRQ(ierr); 301 ierr = PetscFree(aux_local_numbering_1);CHKERRQ(ierr); 302 ierr = PetscFree(dual_dofs_boundary_indices);CHKERRQ(ierr); 303 ierr = PetscFree(all_factors[0]);CHKERRQ(ierr); 304 ierr = PetscFree(all_factors);CHKERRQ(ierr); 305 306 /* Local to global mapping of fetidpmat */ 307 ierr = VecCreate(PETSC_COMM_SELF,&fetidpmat_ctx->lambda_local);CHKERRQ(ierr); 308 ierr = VecSetSizes(fetidpmat_ctx->lambda_local,n_local_lambda,n_local_lambda);CHKERRQ(ierr); 309 ierr = VecSetType(fetidpmat_ctx->lambda_local,VECSEQ);CHKERRQ(ierr); 310 ierr = VecCreate(comm,&lambda_global);CHKERRQ(ierr); 311 ierr = VecSetSizes(lambda_global,PETSC_DECIDE,fetidpmat_ctx->n_lambda);CHKERRQ(ierr); 312 ierr = VecSetType(lambda_global,VECMPI);CHKERRQ(ierr); 313 ierr = ISCreateGeneral(comm,n_local_lambda,l2g_indices,PETSC_OWN_POINTER,&IS_l2g_lambda);CHKERRQ(ierr); 314 ierr = VecScatterCreate(fetidpmat_ctx->lambda_local,(IS)0,lambda_global,IS_l2g_lambda,&fetidpmat_ctx->l2g_lambda);CHKERRQ(ierr); 315 ierr = ISDestroy(&IS_l2g_lambda);CHKERRQ(ierr); 316 317 /* Create local part of B_delta */ 318 ierr = MatCreate(PETSC_COMM_SELF,&fetidpmat_ctx->B_delta); 319 ierr = MatSetSizes(fetidpmat_ctx->B_delta,n_local_lambda,pcis->n_B,n_local_lambda,pcis->n_B);CHKERRQ(ierr); 320 ierr = MatSetType(fetidpmat_ctx->B_delta,MATSEQAIJ);CHKERRQ(ierr); 321 ierr = MatSeqAIJSetPreallocation(fetidpmat_ctx->B_delta,1,NULL);CHKERRQ(ierr); 322 ierr = MatSetOption(fetidpmat_ctx->B_delta,MAT_IGNORE_ZERO_ENTRIES,PETSC_TRUE);CHKERRQ(ierr); 323 for (i=0;i<n_local_lambda;i++) { 324 ierr = MatSetValue(fetidpmat_ctx->B_delta,i,cols_B_delta[i],vals_B_delta[i],INSERT_VALUES);CHKERRQ(ierr); 325 } 326 ierr = PetscFree(vals_B_delta);CHKERRQ(ierr); 327 ierr = MatAssemblyBegin(fetidpmat_ctx->B_delta,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 328 ierr = MatAssemblyEnd (fetidpmat_ctx->B_delta,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 329 330 if (fully_redundant) { 331 ierr = MatCreate(PETSC_COMM_SELF,&ScalingMat); 332 ierr = MatSetSizes(ScalingMat,n_local_lambda,n_local_lambda,n_local_lambda,n_local_lambda);CHKERRQ(ierr); 333 ierr = MatSetType(ScalingMat,MATSEQAIJ);CHKERRQ(ierr); 334 ierr = MatSeqAIJSetPreallocation(ScalingMat,1,NULL);CHKERRQ(ierr); 335 for (i=0;i<n_local_lambda;i++) { 336 ierr = MatSetValue(ScalingMat,i,i,scaling_factors[i],INSERT_VALUES);CHKERRQ(ierr); 337 } 338 ierr = MatAssemblyBegin(ScalingMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 339 ierr = MatAssemblyEnd (ScalingMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 340 ierr = MatMatMult(ScalingMat,fetidpmat_ctx->B_delta,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&fetidpmat_ctx->B_Ddelta);CHKERRQ(ierr); 341 ierr = MatDestroy(&ScalingMat);CHKERRQ(ierr); 342 } else { 343 ierr = MatCreate(PETSC_COMM_SELF,&fetidpmat_ctx->B_Ddelta); 344 ierr = MatSetSizes(fetidpmat_ctx->B_Ddelta,n_local_lambda,pcis->n_B,n_local_lambda,pcis->n_B);CHKERRQ(ierr); 345 ierr = MatSetType(fetidpmat_ctx->B_Ddelta,MATSEQAIJ);CHKERRQ(ierr); 346 ierr = MatSeqAIJSetPreallocation(fetidpmat_ctx->B_Ddelta,1,NULL);CHKERRQ(ierr); 347 for (i=0;i<n_local_lambda;i++) { 348 ierr = MatSetValue(fetidpmat_ctx->B_Ddelta,i,cols_B_delta[i],scaling_factors[i],INSERT_VALUES);CHKERRQ(ierr); 349 } 350 ierr = MatAssemblyBegin(fetidpmat_ctx->B_Ddelta,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 351 ierr = MatAssemblyEnd (fetidpmat_ctx->B_Ddelta,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 352 } 353 ierr = PetscFree(scaling_factors);CHKERRQ(ierr); 354 ierr = PetscFree(cols_B_delta);CHKERRQ(ierr); 355 356 /* Create some vectors needed by fetidp */ 357 ierr = VecDuplicate(pcis->vec1_B,&fetidpmat_ctx->temp_solution_B);CHKERRQ(ierr); 358 ierr = VecDuplicate(pcis->vec1_D,&fetidpmat_ctx->temp_solution_D);CHKERRQ(ierr); 359 360 test_fetidp = PETSC_FALSE; 361 ierr = PetscOptionsGetBool(NULL,"-fetidp_check",&test_fetidp,NULL);CHKERRQ(ierr); 362 363 if (test_fetidp && !pcbddc->use_deluxe_scaling) { 364 365 PetscReal real_value; 366 367 ierr = PetscViewerASCIIGetStdout(comm,&viewer);CHKERRQ(ierr); 368 ierr = PetscViewerASCIISynchronizedAllow(viewer,PETSC_TRUE);CHKERRQ(ierr); 369 ierr = PetscViewerASCIIPrintf(viewer,"----------FETI_DP TESTS--------------\n");CHKERRQ(ierr); 370 ierr = PetscViewerASCIIPrintf(viewer,"All tests should return zero!\n");CHKERRQ(ierr); 371 ierr = PetscViewerASCIIPrintf(viewer,"FETIDP MAT context in the ");CHKERRQ(ierr); 372 if (fully_redundant) { 373 ierr = PetscViewerASCIIPrintf(viewer,"fully redundant case for lagrange multipliers.\n");CHKERRQ(ierr); 374 } else { 375 ierr = PetscViewerASCIIPrintf(viewer,"Non-fully redundant case for lagrange multiplier.\n");CHKERRQ(ierr); 376 } 377 ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 378 379 /******************************************************************/ 380 /* TEST A/B: Test numbering of global lambda dofs */ 381 /******************************************************************/ 382 383 ierr = VecDuplicate(fetidpmat_ctx->lambda_local,&test_vec);CHKERRQ(ierr); 384 ierr = VecSet(lambda_global,1.0);CHKERRQ(ierr); 385 ierr = VecSet(test_vec,1.0);CHKERRQ(ierr); 386 ierr = VecScatterBegin(fetidpmat_ctx->l2g_lambda,lambda_global,fetidpmat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 387 ierr = VecScatterEnd (fetidpmat_ctx->l2g_lambda,lambda_global,fetidpmat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 388 scalar_value = -1.0; 389 ierr = VecAXPY(test_vec,scalar_value,fetidpmat_ctx->lambda_local);CHKERRQ(ierr); 390 ierr = VecNorm(test_vec,NORM_INFINITY,&real_value);CHKERRQ(ierr); 391 ierr = VecDestroy(&test_vec);CHKERRQ(ierr); 392 ierr = PetscViewerASCIISynchronizedPrintf(viewer,"A[%04d]: CHECK glob to loc: % 1.14e\n",rank,real_value);CHKERRQ(ierr); 393 ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 394 if (fully_redundant) { 395 ierr = VecSet(lambda_global,0.0);CHKERRQ(ierr); 396 ierr = VecSet(fetidpmat_ctx->lambda_local,0.5);CHKERRQ(ierr); 397 ierr = VecScatterBegin(fetidpmat_ctx->l2g_lambda,fetidpmat_ctx->lambda_local,lambda_global,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 398 ierr = VecScatterEnd (fetidpmat_ctx->l2g_lambda,fetidpmat_ctx->lambda_local,lambda_global,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 399 ierr = VecSum(lambda_global,&scalar_value);CHKERRQ(ierr); 400 ierr = PetscViewerASCIISynchronizedPrintf(viewer,"B[%04d]: CHECK loc to glob: % 1.14e\n",rank,PetscRealPart(scalar_value)-fetidpmat_ctx->n_lambda);CHKERRQ(ierr); 401 ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 402 } 403 404 /******************************************************************/ 405 /* TEST C: It should holds B_delta*w=0, w\in\widehat{W} */ 406 /* This is the meaning of the B matrix */ 407 /******************************************************************/ 408 409 ierr = VecSetRandom(pcis->vec1_N,NULL);CHKERRQ(ierr); 410 ierr = VecSet(pcis->vec1_global,0.0);CHKERRQ(ierr); 411 ierr = VecScatterBegin(matis->ctx,pcis->vec1_N,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 412 ierr = VecScatterEnd (matis->ctx,pcis->vec1_N,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 413 ierr = VecScatterBegin(matis->ctx,pcis->vec1_global,pcis->vec1_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 414 ierr = VecScatterEnd (matis->ctx,pcis->vec1_global,pcis->vec1_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 415 ierr = VecScatterBegin(pcis->N_to_B,pcis->vec1_N,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 416 ierr = VecScatterEnd (pcis->N_to_B,pcis->vec1_N,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 417 /* Action of B_delta */ 418 ierr = MatMult(fetidpmat_ctx->B_delta,pcis->vec1_B,fetidpmat_ctx->lambda_local);CHKERRQ(ierr); 419 ierr = VecSet(lambda_global,0.0);CHKERRQ(ierr); 420 ierr = VecScatterBegin(fetidpmat_ctx->l2g_lambda,fetidpmat_ctx->lambda_local,lambda_global,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 421 ierr = VecScatterEnd (fetidpmat_ctx->l2g_lambda,fetidpmat_ctx->lambda_local,lambda_global,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 422 ierr = VecNorm(lambda_global,NORM_INFINITY,&real_value);CHKERRQ(ierr); 423 ierr = PetscViewerASCIIPrintf(viewer,"C[coll]: CHECK infty norm of B_delta*w (w continuous): % 1.14e\n",real_value);CHKERRQ(ierr); 424 ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 425 426 /******************************************************************/ 427 /* TEST D: It should holds E_Dw = w - P_Dw w\in\widetilde{W} */ 428 /* E_D = R_D^TR */ 429 /* P_D = B_{D,delta}^T B_{delta} */ 430 /* eq.44 Mandel Tezaur and Dohrmann 2005 */ 431 /******************************************************************/ 432 433 /* compute a random vector in \widetilde{W} */ 434 ierr = VecSetRandom(pcis->vec1_N,NULL);CHKERRQ(ierr); 435 scalar_value = 0.0; /* set zero at vertices */ 436 ierr = VecGetArray(pcis->vec1_N,&array);CHKERRQ(ierr); 437 for (i=0;i<n_vertices;i++) { array[vertex_indices[i]]=scalar_value; } 438 ierr = VecRestoreArray(pcis->vec1_N,&array);CHKERRQ(ierr); 439 /* store w for final comparison */ 440 ierr = VecDuplicate(pcis->vec1_B,&test_vec);CHKERRQ(ierr); 441 ierr = VecScatterBegin(pcis->N_to_B,pcis->vec1_N,test_vec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 442 ierr = VecScatterEnd (pcis->N_to_B,pcis->vec1_N,test_vec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 443 444 /* Jump operator P_D : results stored in pcis->vec1_B */ 445 446 ierr = VecScatterBegin(pcis->N_to_B,pcis->vec1_N,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 447 ierr = VecScatterEnd (pcis->N_to_B,pcis->vec1_N,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 448 /* Action of B_delta */ 449 ierr = MatMult(fetidpmat_ctx->B_delta,pcis->vec1_B,fetidpmat_ctx->lambda_local);CHKERRQ(ierr); 450 ierr = VecSet(lambda_global,0.0);CHKERRQ(ierr); 451 ierr = VecScatterBegin(fetidpmat_ctx->l2g_lambda,fetidpmat_ctx->lambda_local,lambda_global,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 452 ierr = VecScatterEnd (fetidpmat_ctx->l2g_lambda,fetidpmat_ctx->lambda_local,lambda_global,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 453 /* Action of B_Ddelta^T */ 454 ierr = VecScatterBegin(fetidpmat_ctx->l2g_lambda,lambda_global,fetidpmat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 455 ierr = VecScatterEnd (fetidpmat_ctx->l2g_lambda,lambda_global,fetidpmat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 456 ierr = MatMultTranspose(fetidpmat_ctx->B_Ddelta,fetidpmat_ctx->lambda_local,pcis->vec1_B);CHKERRQ(ierr); 457 458 /* Average operator E_D : results stored in pcis->vec2_B */ 459 ierr = VecScatterBegin(pcis->N_to_B,pcis->vec1_N,pcis->vec2_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 460 ierr = VecScatterEnd (pcis->N_to_B,pcis->vec1_N,pcis->vec2_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 461 ierr = PCBDDCScalingExtension(fetidpmat_ctx->pc,pcis->vec2_B,pcis->vec1_global);CHKERRQ(ierr); 462 ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_global,pcis->vec2_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 463 ierr = VecScatterEnd (pcis->global_to_B,pcis->vec1_global,pcis->vec2_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 464 465 /* test E_D=I-P_D */ 466 scalar_value = 1.0; 467 ierr = VecAXPY(pcis->vec1_B,scalar_value,pcis->vec2_B);CHKERRQ(ierr); 468 scalar_value = -1.0; 469 ierr = VecAXPY(pcis->vec1_B,scalar_value,test_vec);CHKERRQ(ierr); 470 ierr = VecNorm(pcis->vec1_B,NORM_INFINITY,&real_value);CHKERRQ(ierr); 471 ierr = VecDestroy(&test_vec);CHKERRQ(ierr); 472 ierr = PetscViewerASCIISynchronizedPrintf(viewer,"D[%04d] CHECK infty norm of E_D + P_D - I: % 1.14e\n",rank,real_value);CHKERRQ(ierr); 473 ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 474 475 /******************************************************************/ 476 /* TEST E: It should holds R_D^TP_Dw=0 w\in\widetilde{W} */ 477 /* eq.48 Mandel Tezaur and Dohrmann 2005 */ 478 /******************************************************************/ 479 480 ierr = VecSetRandom(pcis->vec1_N,NULL);CHKERRQ(ierr); 481 ierr = VecGetArray(pcis->vec1_N,&array);CHKERRQ(ierr); 482 scalar_value = 0.0; /* set zero at vertices */ 483 for (i=0;i<n_vertices;i++) { array[vertex_indices[i]]=scalar_value; } 484 ierr = VecRestoreArray(pcis->vec1_N,&array);CHKERRQ(ierr); 485 486 /* Jump operator P_D : results stored in pcis->vec1_B */ 487 488 ierr = VecScatterBegin(pcis->N_to_B,pcis->vec1_N,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 489 ierr = VecScatterEnd (pcis->N_to_B,pcis->vec1_N,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 490 /* Action of B_delta */ 491 ierr = MatMult(fetidpmat_ctx->B_delta,pcis->vec1_B,fetidpmat_ctx->lambda_local);CHKERRQ(ierr); 492 ierr = VecSet(lambda_global,0.0);CHKERRQ(ierr); 493 ierr = VecScatterBegin(fetidpmat_ctx->l2g_lambda,fetidpmat_ctx->lambda_local,lambda_global,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 494 ierr = VecScatterEnd (fetidpmat_ctx->l2g_lambda,fetidpmat_ctx->lambda_local,lambda_global,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 495 /* Action of B_Ddelta^T */ 496 ierr = VecScatterBegin(fetidpmat_ctx->l2g_lambda,lambda_global,fetidpmat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 497 ierr = VecScatterEnd (fetidpmat_ctx->l2g_lambda,lambda_global,fetidpmat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 498 ierr = MatMultTranspose(fetidpmat_ctx->B_Ddelta,fetidpmat_ctx->lambda_local,pcis->vec1_B);CHKERRQ(ierr); 499 /* scaling */ 500 ierr = PCBDDCScalingExtension(fetidpmat_ctx->pc,pcis->vec1_B,pcis->vec1_global);CHKERRQ(ierr); 501 ierr = VecNorm(pcis->vec1_global,NORM_INFINITY,&real_value);CHKERRQ(ierr); 502 ierr = PetscViewerASCIIPrintf(viewer,"E[coll]: CHECK infty norm of R^T_D P_D: % 1.14e\n",real_value);CHKERRQ(ierr); 503 ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 504 505 if (!fully_redundant) { 506 /******************************************************************/ 507 /* TEST F: It should holds B_{delta}B^T_{D,delta}=I */ 508 /* Corollary thm 14 Mandel Tezaur and Dohrmann 2005 */ 509 /******************************************************************/ 510 ierr = VecDuplicate(lambda_global,&test_vec);CHKERRQ(ierr); 511 ierr = VecSetRandom(lambda_global,NULL);CHKERRQ(ierr); 512 /* Action of B_Ddelta^T */ 513 ierr = VecScatterBegin(fetidpmat_ctx->l2g_lambda,lambda_global,fetidpmat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 514 ierr = VecScatterEnd (fetidpmat_ctx->l2g_lambda,lambda_global,fetidpmat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 515 ierr = MatMultTranspose(fetidpmat_ctx->B_Ddelta,fetidpmat_ctx->lambda_local,pcis->vec1_B);CHKERRQ(ierr); 516 /* Action of B_delta */ 517 ierr = MatMult(fetidpmat_ctx->B_delta,pcis->vec1_B,fetidpmat_ctx->lambda_local);CHKERRQ(ierr); 518 ierr = VecSet(test_vec,0.0);CHKERRQ(ierr); 519 ierr = VecScatterBegin(fetidpmat_ctx->l2g_lambda,fetidpmat_ctx->lambda_local,test_vec,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 520 ierr = VecScatterEnd (fetidpmat_ctx->l2g_lambda,fetidpmat_ctx->lambda_local,test_vec,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 521 scalar_value = -1.0; 522 ierr = VecAXPY(lambda_global,scalar_value,test_vec);CHKERRQ(ierr); 523 ierr = VecNorm(lambda_global,NORM_INFINITY,&real_value);CHKERRQ(ierr); 524 ierr = PetscViewerASCIIPrintf(viewer,"E[coll]: CHECK infty norm of P^T_D - I: % 1.14e\n",real_value);CHKERRQ(ierr); 525 ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 526 ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 527 ierr = VecDestroy(&test_vec);CHKERRQ(ierr); 528 } 529 } 530 /* final cleanup */ 531 ierr = PetscFree(vertex_indices);CHKERRQ(ierr); 532 ierr = VecDestroy(&lambda_global);CHKERRQ(ierr); 533 534 PetscFunctionReturn(0); 535 } 536 537 #undef __FUNCT__ 538 #define __FUNCT__ "PCBDDCSetupFETIDPPCContext" 539 PetscErrorCode PCBDDCSetupFETIDPPCContext(Mat fetimat, FETIDPPC_ctx fetidppc_ctx) 540 { 541 FETIDPMat_ctx mat_ctx; 542 PetscErrorCode ierr; 543 544 PetscFunctionBegin; 545 ierr = MatShellGetContext(fetimat,(void**)&mat_ctx);CHKERRQ(ierr); 546 /* get references from objects created when setting up feti mat context */ 547 ierr = PetscObjectReference((PetscObject)mat_ctx->lambda_local);CHKERRQ(ierr); 548 fetidppc_ctx->lambda_local = mat_ctx->lambda_local; 549 ierr = PetscObjectReference((PetscObject)mat_ctx->B_Ddelta);CHKERRQ(ierr); 550 fetidppc_ctx->B_Ddelta = mat_ctx->B_Ddelta; 551 ierr = PetscObjectReference((PetscObject)mat_ctx->l2g_lambda);CHKERRQ(ierr); 552 fetidppc_ctx->l2g_lambda = mat_ctx->l2g_lambda; 553 PetscFunctionReturn(0); 554 } 555 556 #undef __FUNCT__ 557 #define __FUNCT__ "FETIDPMatMult" 558 PetscErrorCode FETIDPMatMult(Mat fetimat, Vec x, Vec y) 559 { 560 FETIDPMat_ctx mat_ctx; 561 PC_IS *pcis; 562 PetscErrorCode ierr; 563 564 PetscFunctionBegin; 565 ierr = MatShellGetContext(fetimat,(void**)&mat_ctx);CHKERRQ(ierr); 566 pcis = (PC_IS*)mat_ctx->pc->data; 567 /* Application of B_delta^T */ 568 ierr = VecScatterBegin(mat_ctx->l2g_lambda,x,mat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 569 ierr = VecScatterEnd(mat_ctx->l2g_lambda,x,mat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 570 ierr = MatMultTranspose(mat_ctx->B_delta,mat_ctx->lambda_local,pcis->vec1_B);CHKERRQ(ierr); 571 /* Application of \widetilde{S}^-1 */ 572 ierr = VecSet(pcis->vec1_D,0.0);CHKERRQ(ierr); 573 ierr = PCBDDCApplyInterfacePreconditioner(mat_ctx->pc,PETSC_FALSE);CHKERRQ(ierr); 574 /* Application of B_delta */ 575 ierr = MatMult(mat_ctx->B_delta,pcis->vec1_B,mat_ctx->lambda_local);CHKERRQ(ierr); 576 ierr = VecSet(y,0.0);CHKERRQ(ierr); 577 ierr = VecScatterBegin(mat_ctx->l2g_lambda,mat_ctx->lambda_local,y,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 578 ierr = VecScatterEnd(mat_ctx->l2g_lambda,mat_ctx->lambda_local,y,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 579 PetscFunctionReturn(0); 580 } 581 582 #undef __FUNCT__ 583 #define __FUNCT__ "FETIDPMatMultTranspose" 584 PetscErrorCode FETIDPMatMultTranspose(Mat fetimat, Vec x, Vec y) 585 { 586 FETIDPMat_ctx mat_ctx; 587 PC_IS *pcis; 588 PetscErrorCode ierr; 589 590 PetscFunctionBegin; 591 ierr = MatShellGetContext(fetimat,(void**)&mat_ctx);CHKERRQ(ierr); 592 pcis = (PC_IS*)mat_ctx->pc->data; 593 /* Application of B_delta^T */ 594 ierr = VecScatterBegin(mat_ctx->l2g_lambda,x,mat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 595 ierr = VecScatterEnd(mat_ctx->l2g_lambda,x,mat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 596 ierr = MatMultTranspose(mat_ctx->B_delta,mat_ctx->lambda_local,pcis->vec1_B);CHKERRQ(ierr); 597 /* Application of \widetilde{S}^-1 */ 598 ierr = VecSet(pcis->vec1_D,0.0);CHKERRQ(ierr); 599 ierr = PCBDDCApplyInterfacePreconditioner(mat_ctx->pc,PETSC_TRUE);CHKERRQ(ierr); 600 /* Application of B_delta */ 601 ierr = MatMult(mat_ctx->B_delta,pcis->vec1_B,mat_ctx->lambda_local);CHKERRQ(ierr); 602 ierr = VecSet(y,0.0);CHKERRQ(ierr); 603 ierr = VecScatterBegin(mat_ctx->l2g_lambda,mat_ctx->lambda_local,y,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 604 ierr = VecScatterEnd(mat_ctx->l2g_lambda,mat_ctx->lambda_local,y,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 605 PetscFunctionReturn(0); 606 } 607 608 #undef __FUNCT__ 609 #define __FUNCT__ "FETIDPPCApply" 610 PetscErrorCode FETIDPPCApply(PC fetipc, Vec x, Vec y) 611 { 612 FETIDPPC_ctx pc_ctx; 613 PC_IS *pcis; 614 PetscErrorCode ierr; 615 616 PetscFunctionBegin; 617 ierr = PCShellGetContext(fetipc,(void**)&pc_ctx); 618 pcis = (PC_IS*)pc_ctx->pc->data; 619 /* Application of B_Ddelta^T */ 620 ierr = VecScatterBegin(pc_ctx->l2g_lambda,x,pc_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 621 ierr = VecScatterEnd(pc_ctx->l2g_lambda,x,pc_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 622 ierr = VecSet(pcis->vec2_B,0.0);CHKERRQ(ierr); 623 ierr = MatMultTranspose(pc_ctx->B_Ddelta,pc_ctx->lambda_local,pcis->vec2_B);CHKERRQ(ierr); 624 /* Application of S */ 625 ierr = PCBDDCApplySchur(pc_ctx->pc,pcis->vec2_B,pcis->vec1_B,(Vec)0,pcis->vec1_D,pcis->vec2_D);CHKERRQ(ierr); 626 /* Application of B_Ddelta */ 627 ierr = MatMult(pc_ctx->B_Ddelta,pcis->vec1_B,pc_ctx->lambda_local);CHKERRQ(ierr); 628 ierr = VecSet(y,0.0);CHKERRQ(ierr); 629 ierr = VecScatterBegin(pc_ctx->l2g_lambda,pc_ctx->lambda_local,y,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 630 ierr = VecScatterEnd(pc_ctx->l2g_lambda,pc_ctx->lambda_local,y,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 631 PetscFunctionReturn(0); 632 } 633 634 #undef __FUNCT__ 635 #define __FUNCT__ "FETIDPPCApplyTranspose" 636 PetscErrorCode FETIDPPCApplyTranspose(PC fetipc, Vec x, Vec y) 637 { 638 FETIDPPC_ctx pc_ctx; 639 PC_IS *pcis; 640 PetscErrorCode ierr; 641 642 PetscFunctionBegin; 643 ierr = PCShellGetContext(fetipc,(void**)&pc_ctx); 644 pcis = (PC_IS*)pc_ctx->pc->data; 645 /* Application of B_Ddelta^T */ 646 ierr = VecScatterBegin(pc_ctx->l2g_lambda,x,pc_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 647 ierr = VecScatterEnd(pc_ctx->l2g_lambda,x,pc_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 648 ierr = VecSet(pcis->vec2_B,0.0);CHKERRQ(ierr); 649 ierr = MatMultTranspose(pc_ctx->B_Ddelta,pc_ctx->lambda_local,pcis->vec2_B);CHKERRQ(ierr); 650 /* Application of S */ 651 ierr = PCBDDCApplySchurTranspose(pc_ctx->pc,pcis->vec2_B,pcis->vec1_B,(Vec)0,pcis->vec1_D,pcis->vec2_D);CHKERRQ(ierr); 652 /* Application of B_Ddelta */ 653 ierr = MatMult(pc_ctx->B_Ddelta,pcis->vec1_B,pc_ctx->lambda_local);CHKERRQ(ierr); 654 ierr = VecSet(y,0.0);CHKERRQ(ierr); 655 ierr = VecScatterBegin(pc_ctx->l2g_lambda,pc_ctx->lambda_local,y,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 656 ierr = VecScatterEnd(pc_ctx->l2g_lambda,pc_ctx->lambda_local,y,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 657 PetscFunctionReturn(0); 658 } 659