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