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