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