1674ae819SStefano Zampini #include "bddc.h" 2674ae819SStefano Zampini #include "bddcprivate.h" 3674ae819SStefano Zampini 4674ae819SStefano Zampini #undef __FUNCT__ 5674ae819SStefano Zampini #define __FUNCT__ "PCBDDCCreateFETIDPMatContext" 6674ae819SStefano Zampini PetscErrorCode PCBDDCCreateFETIDPMatContext(PC pc, FETIDPMat_ctx *fetidpmat_ctx) 7674ae819SStefano Zampini { 8674ae819SStefano Zampini FETIDPMat_ctx newctx; 9674ae819SStefano Zampini PetscErrorCode ierr; 10674ae819SStefano Zampini 11674ae819SStefano Zampini PetscFunctionBegin; 12674ae819SStefano Zampini ierr = PetscMalloc(sizeof(*newctx),&newctx);CHKERRQ(ierr); 13674ae819SStefano Zampini newctx->lambda_local = 0; 14674ae819SStefano Zampini newctx->temp_solution_B = 0; 15674ae819SStefano Zampini newctx->temp_solution_D = 0; 16674ae819SStefano Zampini newctx->B_delta = 0; 17674ae819SStefano Zampini newctx->B_Ddelta = 0; /* theoretically belongs to the FETIDP preconditioner */ 18674ae819SStefano Zampini newctx->l2g_lambda = 0; 19674ae819SStefano Zampini /* increase the reference count for BDDC preconditioner */ 20674ae819SStefano Zampini ierr = PetscObjectReference((PetscObject)pc);CHKERRQ(ierr); 21674ae819SStefano Zampini newctx->pc = pc; 22674ae819SStefano Zampini *fetidpmat_ctx = newctx; 23674ae819SStefano Zampini PetscFunctionReturn(0); 24674ae819SStefano Zampini } 25674ae819SStefano Zampini 26674ae819SStefano Zampini #undef __FUNCT__ 27674ae819SStefano Zampini #define __FUNCT__ "PCBDDCCreateFETIDPPCContext" 28674ae819SStefano Zampini PetscErrorCode PCBDDCCreateFETIDPPCContext(PC pc, FETIDPPC_ctx *fetidppc_ctx) 29674ae819SStefano Zampini { 30674ae819SStefano Zampini FETIDPPC_ctx newctx; 31674ae819SStefano Zampini PetscErrorCode ierr; 32674ae819SStefano Zampini 33674ae819SStefano Zampini PetscFunctionBegin; 34674ae819SStefano Zampini ierr = PetscMalloc(sizeof(*newctx),&newctx);CHKERRQ(ierr); 35674ae819SStefano Zampini newctx->lambda_local = 0; 36674ae819SStefano Zampini newctx->B_Ddelta = 0; 37674ae819SStefano Zampini newctx->l2g_lambda = 0; 38674ae819SStefano Zampini /* increase the reference count for BDDC preconditioner */ 39674ae819SStefano Zampini ierr = PetscObjectReference((PetscObject)pc);CHKERRQ(ierr); 40674ae819SStefano Zampini newctx->pc = pc; 41674ae819SStefano Zampini *fetidppc_ctx = newctx; 42674ae819SStefano Zampini PetscFunctionReturn(0); 43674ae819SStefano Zampini } 44674ae819SStefano Zampini 45674ae819SStefano Zampini #undef __FUNCT__ 46674ae819SStefano Zampini #define __FUNCT__ "PCBDDCDestroyFETIDPMat" 47674ae819SStefano Zampini PetscErrorCode PCBDDCDestroyFETIDPMat(Mat A) 48674ae819SStefano Zampini { 49674ae819SStefano Zampini FETIDPMat_ctx mat_ctx; 50674ae819SStefano Zampini PetscErrorCode ierr; 51674ae819SStefano Zampini 52674ae819SStefano Zampini PetscFunctionBegin; 53674ae819SStefano Zampini ierr = MatShellGetContext(A,(void**)&mat_ctx);CHKERRQ(ierr); 54674ae819SStefano Zampini ierr = VecDestroy(&mat_ctx->lambda_local);CHKERRQ(ierr); 55674ae819SStefano Zampini ierr = VecDestroy(&mat_ctx->temp_solution_D);CHKERRQ(ierr); 56674ae819SStefano Zampini ierr = VecDestroy(&mat_ctx->temp_solution_B);CHKERRQ(ierr); 57674ae819SStefano Zampini ierr = MatDestroy(&mat_ctx->B_delta);CHKERRQ(ierr); 58674ae819SStefano Zampini ierr = MatDestroy(&mat_ctx->B_Ddelta);CHKERRQ(ierr); 59674ae819SStefano Zampini ierr = VecScatterDestroy(&mat_ctx->l2g_lambda);CHKERRQ(ierr); 60674ae819SStefano Zampini ierr = PCDestroy(&mat_ctx->pc);CHKERRQ(ierr); /* decrease PCBDDC reference count */ 61674ae819SStefano Zampini ierr = PetscFree(mat_ctx);CHKERRQ(ierr); 62674ae819SStefano Zampini PetscFunctionReturn(0); 63674ae819SStefano Zampini } 64674ae819SStefano Zampini 65674ae819SStefano Zampini #undef __FUNCT__ 66674ae819SStefano Zampini #define __FUNCT__ "PCBDDCDestroyFETIDPPC" 67674ae819SStefano Zampini PetscErrorCode PCBDDCDestroyFETIDPPC(PC pc) 68674ae819SStefano Zampini { 69674ae819SStefano Zampini FETIDPPC_ctx pc_ctx; 70674ae819SStefano Zampini PetscErrorCode ierr; 71674ae819SStefano Zampini 72674ae819SStefano Zampini PetscFunctionBegin; 73674ae819SStefano Zampini ierr = PCShellGetContext(pc,(void**)&pc_ctx);CHKERRQ(ierr); 74674ae819SStefano Zampini ierr = VecDestroy(&pc_ctx->lambda_local);CHKERRQ(ierr); 75674ae819SStefano Zampini ierr = MatDestroy(&pc_ctx->B_Ddelta);CHKERRQ(ierr); 76674ae819SStefano Zampini ierr = VecScatterDestroy(&pc_ctx->l2g_lambda);CHKERRQ(ierr); 77674ae819SStefano Zampini ierr = PCDestroy(&pc_ctx->pc);CHKERRQ(ierr); /* decrease PCBDDC reference count */ 78674ae819SStefano Zampini ierr = PetscFree(pc_ctx);CHKERRQ(ierr); 79674ae819SStefano Zampini PetscFunctionReturn(0); 80674ae819SStefano Zampini } 81674ae819SStefano Zampini 82674ae819SStefano Zampini #undef __FUNCT__ 83674ae819SStefano Zampini #define __FUNCT__ "PCBDDCSetupFETIDPMatContext" 84674ae819SStefano Zampini PetscErrorCode PCBDDCSetupFETIDPMatContext(FETIDPMat_ctx fetidpmat_ctx ) 85674ae819SStefano Zampini { 86674ae819SStefano Zampini PetscErrorCode ierr; 87674ae819SStefano Zampini PC_IS *pcis=(PC_IS*)fetidpmat_ctx->pc->data; 88674ae819SStefano Zampini PC_BDDC *pcbddc=(PC_BDDC*)fetidpmat_ctx->pc->data; 89674ae819SStefano Zampini PCBDDCGraph mat_graph=pcbddc->mat_graph; 90674ae819SStefano Zampini Mat_IS *matis = (Mat_IS*)fetidpmat_ctx->pc->pmat->data; 91674ae819SStefano Zampini MPI_Comm comm; 92674ae819SStefano Zampini Mat ScalingMat; 93674ae819SStefano Zampini Vec lambda_global; 94674ae819SStefano Zampini IS IS_l2g_lambda; 95674ae819SStefano Zampini PetscBool skip_node,fully_redundant; 96674ae819SStefano Zampini PetscInt i,j,k,s,n_boundary_dofs,n_global_lambda,n_vertices,partial_sum; 97674ae819SStefano Zampini PetscInt n_local_lambda,n_lambda_for_dof,dual_size,n_neg_values,n_pos_values; 98674ae819SStefano Zampini PetscMPIInt rank,nprocs,buf_size,neigh; 99674ae819SStefano Zampini PetscScalar scalar_value; 100674ae819SStefano Zampini PetscInt *vertex_indices; 101674ae819SStefano Zampini PetscInt *dual_dofs_boundary_indices,*aux_local_numbering_1,*aux_global_numbering; 102674ae819SStefano Zampini PetscInt *aux_sums,*cols_B_delta,*l2g_indices; 103674ae819SStefano Zampini PetscScalar *array,*scaling_factors,*vals_B_delta; 104674ae819SStefano Zampini PetscInt *aux_local_numbering_2; 105674ae819SStefano Zampini /* For communication of scaling factors */ 106674ae819SStefano Zampini PetscInt *ptrs_buffer,neigh_position; 107674ae819SStefano Zampini PetscScalar **all_factors,*send_buffer,*recv_buffer; 108674ae819SStefano Zampini MPI_Request *send_reqs,*recv_reqs; 109674ae819SStefano Zampini /* tests */ 110674ae819SStefano Zampini Vec test_vec; 111674ae819SStefano Zampini PetscBool test_fetidp; 112674ae819SStefano Zampini PetscViewer viewer; 113674ae819SStefano Zampini 114674ae819SStefano Zampini PetscFunctionBegin; 115674ae819SStefano Zampini ierr = PetscObjectGetComm((PetscObject)(fetidpmat_ctx->pc),&comm);CHKERRQ(ierr); 116674ae819SStefano Zampini ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 117674ae819SStefano Zampini ierr = MPI_Comm_size(comm,&nprocs);CHKERRQ(ierr); 118674ae819SStefano Zampini 119674ae819SStefano Zampini /* Default type of lagrange multipliers is non-redundant */ 120674ae819SStefano Zampini fully_redundant = PETSC_FALSE; 121674ae819SStefano Zampini ierr = PetscOptionsGetBool(NULL,"-fetidp_fullyredundant",&fully_redundant,NULL);CHKERRQ(ierr); 122674ae819SStefano Zampini 123674ae819SStefano Zampini /* Evaluate local and global number of lagrange multipliers */ 124674ae819SStefano Zampini ierr = VecSet(pcis->vec1_N,0.0);CHKERRQ(ierr); 125674ae819SStefano Zampini n_local_lambda = 0; 126674ae819SStefano Zampini partial_sum = 0; 127674ae819SStefano Zampini n_boundary_dofs = 0; 128674ae819SStefano Zampini s = 0; 129674ae819SStefano Zampini /* Get Vertices used to define the BDDC */ 130674ae819SStefano Zampini ierr = PCBDDCGetPrimalVerticesLocalIdx(fetidpmat_ctx->pc,&n_vertices,&vertex_indices);CHKERRQ(ierr); 131674ae819SStefano Zampini dual_size = pcis->n_B-n_vertices; 132674ae819SStefano Zampini ierr = PetscSortInt(n_vertices,vertex_indices);CHKERRQ(ierr); 133*785e854fSJed Brown ierr = PetscMalloc1(dual_size,&dual_dofs_boundary_indices);CHKERRQ(ierr); 134*785e854fSJed Brown ierr = PetscMalloc1(dual_size,&aux_local_numbering_1);CHKERRQ(ierr); 135*785e854fSJed Brown ierr = PetscMalloc1(dual_size,&aux_local_numbering_2);CHKERRQ(ierr); 136674ae819SStefano Zampini 137674ae819SStefano Zampini ierr = VecGetArray(pcis->vec1_N,&array);CHKERRQ(ierr); 138674ae819SStefano Zampini for (i=0;i<pcis->n;i++){ 139674ae819SStefano Zampini j = mat_graph->count[i]; /* RECALL: mat_graph->count[i] does not count myself */ 140674ae819SStefano Zampini if ( j > 0 ) { 141674ae819SStefano Zampini n_boundary_dofs++; 142674ae819SStefano Zampini } 143674ae819SStefano Zampini skip_node = PETSC_FALSE; 144674ae819SStefano Zampini if ( s < n_vertices && vertex_indices[s]==i) { /* it works for a sorted set of vertices */ 145674ae819SStefano Zampini skip_node = PETSC_TRUE; 146674ae819SStefano Zampini s++; 147674ae819SStefano Zampini } 148674ae819SStefano Zampini if (j < 1) { 149674ae819SStefano Zampini skip_node = PETSC_TRUE; 150674ae819SStefano Zampini } 151674ae819SStefano Zampini if ( !skip_node ) { 152674ae819SStefano Zampini if (fully_redundant) { 153674ae819SStefano Zampini /* fully redundant set of lagrange multipliers */ 154674ae819SStefano Zampini n_lambda_for_dof = (j*(j+1))/2; 155674ae819SStefano Zampini } else { 156674ae819SStefano Zampini n_lambda_for_dof = j; 157674ae819SStefano Zampini } 158674ae819SStefano Zampini n_local_lambda += j; 159674ae819SStefano Zampini /* needed to evaluate global number of lagrange multipliers */ 160674ae819SStefano Zampini array[i]=(1.0*n_lambda_for_dof)/(j+1.0); /* already scaled for the next global sum */ 161674ae819SStefano Zampini /* store some data needed */ 162674ae819SStefano Zampini dual_dofs_boundary_indices[partial_sum] = n_boundary_dofs-1; 163674ae819SStefano Zampini aux_local_numbering_1[partial_sum] = i; 164674ae819SStefano Zampini aux_local_numbering_2[partial_sum] = n_lambda_for_dof; 165674ae819SStefano Zampini partial_sum++; 166674ae819SStefano Zampini } 167674ae819SStefano Zampini } 168674ae819SStefano Zampini ierr = VecRestoreArray(pcis->vec1_N,&array);CHKERRQ(ierr); 169674ae819SStefano Zampini 170674ae819SStefano Zampini ierr = VecSet(pcis->vec1_global,0.0);CHKERRQ(ierr); 171674ae819SStefano Zampini ierr = VecScatterBegin(matis->ctx,pcis->vec1_N,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 172674ae819SStefano Zampini ierr = VecScatterEnd (matis->ctx,pcis->vec1_N,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 173674ae819SStefano Zampini ierr = VecSum(pcis->vec1_global,&scalar_value);CHKERRQ(ierr); 1745b08dc53SStefano Zampini fetidpmat_ctx->n_lambda = (PetscInt)PetscRealPart(scalar_value); 175674ae819SStefano Zampini 176674ae819SStefano Zampini /* compute global ordering of lagrange multipliers and associate l2g map */ 177674ae819SStefano Zampini ierr = PCBDDCSubsetNumbering(comm,matis->mapping,partial_sum,aux_local_numbering_1,aux_local_numbering_2,&i,&aux_global_numbering);CHKERRQ(ierr); 178674ae819SStefano Zampini if (i != fetidpmat_ctx->n_lambda) { 179674ae819SStefano Zampini SETERRQ3(PETSC_COMM_WORLD,PETSC_ERR_PLIB,"Error in %s: global number of multipliers mismatch! (%d!=%d)\n",__FUNCT__,fetidpmat_ctx->n_lambda,i); 180674ae819SStefano Zampini } 181674ae819SStefano Zampini ierr = PetscFree(aux_local_numbering_2);CHKERRQ(ierr); 182674ae819SStefano Zampini 183674ae819SStefano Zampini /* init data for scaling factors exchange */ 184674ae819SStefano Zampini partial_sum = 0; 185674ae819SStefano Zampini j = 0; 186*785e854fSJed Brown ierr = PetscMalloc1(pcis->n_neigh,&ptrs_buffer);CHKERRQ(ierr); 187*785e854fSJed Brown ierr = PetscMalloc1((pcis->n_neigh-1),&send_reqs);CHKERRQ(ierr); 188*785e854fSJed Brown ierr = PetscMalloc1((pcis->n_neigh-1),&recv_reqs);CHKERRQ(ierr); 189*785e854fSJed Brown ierr = PetscMalloc1(pcis->n,&all_factors);CHKERRQ(ierr); 190674ae819SStefano Zampini ptrs_buffer[0]=0; 191674ae819SStefano Zampini for (i=1;i<pcis->n_neigh;i++) { 192674ae819SStefano Zampini partial_sum += pcis->n_shared[i]; 193674ae819SStefano Zampini ptrs_buffer[i] = ptrs_buffer[i-1]+pcis->n_shared[i]; 194674ae819SStefano Zampini } 195*785e854fSJed Brown ierr = PetscMalloc1(partial_sum,&send_buffer);CHKERRQ(ierr); 196*785e854fSJed Brown ierr = PetscMalloc1(partial_sum,&recv_buffer);CHKERRQ(ierr); 197*785e854fSJed Brown ierr = PetscMalloc1(partial_sum,&all_factors[0]);CHKERRQ(ierr); 198674ae819SStefano Zampini for (i=0;i<pcis->n-1;i++) { 199674ae819SStefano Zampini j = mat_graph->count[i]; 200674ae819SStefano Zampini all_factors[i+1]=all_factors[i]+j; 201674ae819SStefano Zampini } 202674ae819SStefano Zampini /* scatter B scaling to N vec */ 203674ae819SStefano Zampini ierr = VecScatterBegin(pcis->N_to_B,pcis->D,pcis->vec1_N,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 204674ae819SStefano Zampini ierr = VecScatterEnd(pcis->N_to_B,pcis->D,pcis->vec1_N,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 205674ae819SStefano Zampini /* communications */ 206674ae819SStefano Zampini ierr = VecGetArray(pcis->vec1_N,&array);CHKERRQ(ierr); 207674ae819SStefano Zampini for (i=1;i<pcis->n_neigh;i++) { 208674ae819SStefano Zampini for (j=0;j<pcis->n_shared[i];j++) { 209674ae819SStefano Zampini send_buffer[ptrs_buffer[i-1]+j]=array[pcis->shared[i][j]]; 210674ae819SStefano Zampini } 211674ae819SStefano Zampini ierr = PetscMPIIntCast(ptrs_buffer[i]-ptrs_buffer[i-1],&buf_size);CHKERRQ(ierr); 212674ae819SStefano Zampini ierr = PetscMPIIntCast(pcis->neigh[i],&neigh);CHKERRQ(ierr); 213674ae819SStefano Zampini ierr = MPI_Isend(&send_buffer[ptrs_buffer[i-1]],buf_size,MPIU_SCALAR,neigh,0,comm,&send_reqs[i-1]);CHKERRQ(ierr); 214674ae819SStefano Zampini ierr = MPI_Irecv(&recv_buffer[ptrs_buffer[i-1]],buf_size,MPIU_SCALAR,neigh,0,comm,&recv_reqs[i-1]);CHKERRQ(ierr); 215674ae819SStefano Zampini } 216674ae819SStefano Zampini ierr = VecRestoreArray(pcis->vec1_N,&array);CHKERRQ(ierr); 217674ae819SStefano Zampini ierr = MPI_Waitall((pcis->n_neigh-1),recv_reqs,MPI_STATUSES_IGNORE);CHKERRQ(ierr); 218674ae819SStefano Zampini /* put values in correct places */ 219674ae819SStefano Zampini for (i=1;i<pcis->n_neigh;i++) { 220674ae819SStefano Zampini for (j=0;j<pcis->n_shared[i];j++) { 221674ae819SStefano Zampini k = pcis->shared[i][j]; 222674ae819SStefano Zampini neigh_position = 0; 223674ae819SStefano Zampini while(mat_graph->neighbours_set[k][neigh_position] != pcis->neigh[i]) {neigh_position++;} 224674ae819SStefano Zampini all_factors[k][neigh_position]=recv_buffer[ptrs_buffer[i-1]+j]; 225674ae819SStefano Zampini } 226674ae819SStefano Zampini } 227674ae819SStefano Zampini ierr = MPI_Waitall((pcis->n_neigh-1),send_reqs,MPI_STATUSES_IGNORE);CHKERRQ(ierr); 228674ae819SStefano Zampini ierr = PetscFree(send_reqs);CHKERRQ(ierr); 229674ae819SStefano Zampini ierr = PetscFree(recv_reqs);CHKERRQ(ierr); 230674ae819SStefano Zampini ierr = PetscFree(send_buffer);CHKERRQ(ierr); 231674ae819SStefano Zampini ierr = PetscFree(recv_buffer);CHKERRQ(ierr); 232674ae819SStefano Zampini ierr = PetscFree(ptrs_buffer);CHKERRQ(ierr); 233674ae819SStefano Zampini 234674ae819SStefano Zampini /* Compute B and B_delta (local actions) */ 235*785e854fSJed Brown ierr = PetscMalloc1(pcis->n_neigh,&aux_sums);CHKERRQ(ierr); 236*785e854fSJed Brown ierr = PetscMalloc1(n_local_lambda,&l2g_indices);CHKERRQ(ierr); 237*785e854fSJed Brown ierr = PetscMalloc1(n_local_lambda,&vals_B_delta);CHKERRQ(ierr); 238*785e854fSJed Brown ierr = PetscMalloc1(n_local_lambda,&cols_B_delta);CHKERRQ(ierr); 239*785e854fSJed Brown ierr = PetscMalloc1(n_local_lambda,&scaling_factors);CHKERRQ(ierr); 240674ae819SStefano Zampini n_global_lambda=0; 241674ae819SStefano Zampini partial_sum=0; 242674ae819SStefano Zampini for (i=0;i<dual_size;i++) { 243674ae819SStefano Zampini n_global_lambda = aux_global_numbering[i]; 244674ae819SStefano Zampini j = mat_graph->count[aux_local_numbering_1[i]]; 245674ae819SStefano Zampini aux_sums[0]=0; 246674ae819SStefano Zampini for (s=1;s<j;s++) { 247674ae819SStefano Zampini aux_sums[s]=aux_sums[s-1]+j-s+1; 248674ae819SStefano Zampini } 249674ae819SStefano Zampini array = all_factors[aux_local_numbering_1[i]]; 250674ae819SStefano Zampini n_neg_values = 0; 2512a7da448SStefano Zampini while(n_neg_values < j && mat_graph->neighbours_set[aux_local_numbering_1[i]][n_neg_values] < rank) {n_neg_values++;} 252674ae819SStefano Zampini n_pos_values = j - n_neg_values; 253674ae819SStefano Zampini if (fully_redundant) { 254674ae819SStefano Zampini for (s=0;s<n_neg_values;s++) { 255674ae819SStefano Zampini l2g_indices [partial_sum+s]=aux_sums[s]+n_neg_values-s-1+n_global_lambda; 256674ae819SStefano Zampini cols_B_delta [partial_sum+s]=dual_dofs_boundary_indices[i]; 257674ae819SStefano Zampini vals_B_delta [partial_sum+s]=-1.0; 258674ae819SStefano Zampini scaling_factors[partial_sum+s]=array[s]; 259674ae819SStefano Zampini } 260674ae819SStefano Zampini for (s=0;s<n_pos_values;s++) { 261674ae819SStefano Zampini l2g_indices [partial_sum+s+n_neg_values]=aux_sums[n_neg_values]+s+n_global_lambda; 262674ae819SStefano Zampini cols_B_delta [partial_sum+s+n_neg_values]=dual_dofs_boundary_indices[i]; 263674ae819SStefano Zampini vals_B_delta [partial_sum+s+n_neg_values]=1.0; 264674ae819SStefano Zampini scaling_factors[partial_sum+s+n_neg_values]=array[s+n_neg_values]; 265674ae819SStefano Zampini } 266674ae819SStefano Zampini partial_sum += j; 267674ae819SStefano Zampini } else { 268674ae819SStefano Zampini /* l2g_indices and default cols and vals of B_delta */ 269674ae819SStefano Zampini for (s=0;s<j;s++) { 270674ae819SStefano Zampini l2g_indices [partial_sum+s]=n_global_lambda+s; 271674ae819SStefano Zampini cols_B_delta [partial_sum+s]=dual_dofs_boundary_indices[i]; 272674ae819SStefano Zampini vals_B_delta [partial_sum+s]=0.0; 273674ae819SStefano Zampini } 274674ae819SStefano Zampini /* B_delta */ 275674ae819SStefano Zampini if ( n_neg_values > 0 ) { /* there's a rank next to me to the left */ 276674ae819SStefano Zampini vals_B_delta [partial_sum+n_neg_values-1]=-1.0; 277674ae819SStefano Zampini } 278674ae819SStefano Zampini if ( n_neg_values < j ) { /* there's a rank next to me to the right */ 279674ae819SStefano Zampini vals_B_delta [partial_sum+n_neg_values]=1.0; 280674ae819SStefano Zampini } 281674ae819SStefano Zampini /* scaling as in Klawonn-Widlund 1999*/ 282674ae819SStefano Zampini for (s=0;s<n_neg_values;s++) { 283674ae819SStefano Zampini scalar_value = 0.0; 284674ae819SStefano Zampini for (k=0;k<s+1;k++) { 285674ae819SStefano Zampini scalar_value += array[k]; 286674ae819SStefano Zampini } 287674ae819SStefano Zampini scaling_factors[partial_sum+s] = -scalar_value; 288674ae819SStefano Zampini } 289674ae819SStefano Zampini for (s=0;s<n_pos_values;s++) { 290674ae819SStefano Zampini scalar_value = 0.0; 291674ae819SStefano Zampini for (k=s+n_neg_values;k<j;k++) { 292674ae819SStefano Zampini scalar_value += array[k]; 293674ae819SStefano Zampini } 294674ae819SStefano Zampini scaling_factors[partial_sum+s+n_neg_values] = scalar_value; 295674ae819SStefano Zampini } 296674ae819SStefano Zampini partial_sum += j; 297674ae819SStefano Zampini } 298674ae819SStefano Zampini } 299674ae819SStefano Zampini ierr = PetscFree(aux_global_numbering);CHKERRQ(ierr); 300674ae819SStefano Zampini ierr = PetscFree(aux_sums);CHKERRQ(ierr); 301674ae819SStefano Zampini ierr = PetscFree(aux_local_numbering_1);CHKERRQ(ierr); 302674ae819SStefano Zampini ierr = PetscFree(dual_dofs_boundary_indices);CHKERRQ(ierr); 303674ae819SStefano Zampini ierr = PetscFree(all_factors[0]);CHKERRQ(ierr); 304674ae819SStefano Zampini ierr = PetscFree(all_factors);CHKERRQ(ierr); 305674ae819SStefano Zampini 306674ae819SStefano Zampini /* Local to global mapping of fetidpmat */ 307674ae819SStefano Zampini ierr = VecCreate(PETSC_COMM_SELF,&fetidpmat_ctx->lambda_local);CHKERRQ(ierr); 308674ae819SStefano Zampini ierr = VecSetSizes(fetidpmat_ctx->lambda_local,n_local_lambda,n_local_lambda);CHKERRQ(ierr); 309674ae819SStefano Zampini ierr = VecSetType(fetidpmat_ctx->lambda_local,VECSEQ);CHKERRQ(ierr); 310674ae819SStefano Zampini ierr = VecCreate(comm,&lambda_global);CHKERRQ(ierr); 311674ae819SStefano Zampini ierr = VecSetSizes(lambda_global,PETSC_DECIDE,fetidpmat_ctx->n_lambda);CHKERRQ(ierr); 312674ae819SStefano Zampini ierr = VecSetType(lambda_global,VECMPI);CHKERRQ(ierr); 313674ae819SStefano Zampini ierr = ISCreateGeneral(comm,n_local_lambda,l2g_indices,PETSC_OWN_POINTER,&IS_l2g_lambda);CHKERRQ(ierr); 314674ae819SStefano Zampini ierr = VecScatterCreate(fetidpmat_ctx->lambda_local,(IS)0,lambda_global,IS_l2g_lambda,&fetidpmat_ctx->l2g_lambda);CHKERRQ(ierr); 315674ae819SStefano Zampini ierr = ISDestroy(&IS_l2g_lambda);CHKERRQ(ierr); 316674ae819SStefano Zampini 317674ae819SStefano Zampini /* Create local part of B_delta */ 318674ae819SStefano Zampini ierr = MatCreate(PETSC_COMM_SELF,&fetidpmat_ctx->B_delta); 319674ae819SStefano Zampini ierr = MatSetSizes(fetidpmat_ctx->B_delta,n_local_lambda,pcis->n_B,n_local_lambda,pcis->n_B);CHKERRQ(ierr); 320674ae819SStefano Zampini ierr = MatSetType(fetidpmat_ctx->B_delta,MATSEQAIJ);CHKERRQ(ierr); 321674ae819SStefano Zampini ierr = MatSeqAIJSetPreallocation(fetidpmat_ctx->B_delta,1,NULL);CHKERRQ(ierr); 322674ae819SStefano Zampini ierr = MatSetOption(fetidpmat_ctx->B_delta,MAT_IGNORE_ZERO_ENTRIES,PETSC_TRUE);CHKERRQ(ierr); 323674ae819SStefano Zampini for (i=0;i<n_local_lambda;i++) { 324674ae819SStefano Zampini ierr = MatSetValue(fetidpmat_ctx->B_delta,i,cols_B_delta[i],vals_B_delta[i],INSERT_VALUES);CHKERRQ(ierr); 325674ae819SStefano Zampini } 326674ae819SStefano Zampini ierr = PetscFree(vals_B_delta);CHKERRQ(ierr); 327674ae819SStefano Zampini ierr = MatAssemblyBegin(fetidpmat_ctx->B_delta,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 328674ae819SStefano Zampini ierr = MatAssemblyEnd (fetidpmat_ctx->B_delta,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 329674ae819SStefano Zampini 330674ae819SStefano Zampini if (fully_redundant) { 331674ae819SStefano Zampini ierr = MatCreate(PETSC_COMM_SELF,&ScalingMat); 332674ae819SStefano Zampini ierr = MatSetSizes(ScalingMat,n_local_lambda,n_local_lambda,n_local_lambda,n_local_lambda);CHKERRQ(ierr); 333674ae819SStefano Zampini ierr = MatSetType(ScalingMat,MATSEQAIJ);CHKERRQ(ierr); 334674ae819SStefano Zampini ierr = MatSeqAIJSetPreallocation(ScalingMat,1,NULL);CHKERRQ(ierr); 335674ae819SStefano Zampini for (i=0;i<n_local_lambda;i++) { 336674ae819SStefano Zampini ierr = MatSetValue(ScalingMat,i,i,scaling_factors[i],INSERT_VALUES);CHKERRQ(ierr); 337674ae819SStefano Zampini } 338674ae819SStefano Zampini ierr = MatAssemblyBegin(ScalingMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 339674ae819SStefano Zampini ierr = MatAssemblyEnd (ScalingMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 340674ae819SStefano Zampini ierr = MatMatMult(ScalingMat,fetidpmat_ctx->B_delta,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&fetidpmat_ctx->B_Ddelta);CHKERRQ(ierr); 341674ae819SStefano Zampini ierr = MatDestroy(&ScalingMat);CHKERRQ(ierr); 342674ae819SStefano Zampini } else { 343674ae819SStefano Zampini ierr = MatCreate(PETSC_COMM_SELF,&fetidpmat_ctx->B_Ddelta); 344674ae819SStefano Zampini ierr = MatSetSizes(fetidpmat_ctx->B_Ddelta,n_local_lambda,pcis->n_B,n_local_lambda,pcis->n_B);CHKERRQ(ierr); 345674ae819SStefano Zampini ierr = MatSetType(fetidpmat_ctx->B_Ddelta,MATSEQAIJ);CHKERRQ(ierr); 346674ae819SStefano Zampini ierr = MatSeqAIJSetPreallocation(fetidpmat_ctx->B_Ddelta,1,NULL);CHKERRQ(ierr); 347674ae819SStefano Zampini for (i=0;i<n_local_lambda;i++) { 348674ae819SStefano Zampini ierr = MatSetValue(fetidpmat_ctx->B_Ddelta,i,cols_B_delta[i],scaling_factors[i],INSERT_VALUES);CHKERRQ(ierr); 349674ae819SStefano Zampini } 350674ae819SStefano Zampini ierr = MatAssemblyBegin(fetidpmat_ctx->B_Ddelta,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 351674ae819SStefano Zampini ierr = MatAssemblyEnd (fetidpmat_ctx->B_Ddelta,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 352674ae819SStefano Zampini } 353674ae819SStefano Zampini ierr = PetscFree(scaling_factors);CHKERRQ(ierr); 354674ae819SStefano Zampini ierr = PetscFree(cols_B_delta);CHKERRQ(ierr); 355674ae819SStefano Zampini 356674ae819SStefano Zampini /* Create some vectors needed by fetidp */ 357674ae819SStefano Zampini ierr = VecDuplicate(pcis->vec1_B,&fetidpmat_ctx->temp_solution_B);CHKERRQ(ierr); 358674ae819SStefano Zampini ierr = VecDuplicate(pcis->vec1_D,&fetidpmat_ctx->temp_solution_D);CHKERRQ(ierr); 359674ae819SStefano Zampini 360674ae819SStefano Zampini test_fetidp = PETSC_FALSE; 361674ae819SStefano Zampini ierr = PetscOptionsGetBool(NULL,"-fetidp_check",&test_fetidp,NULL);CHKERRQ(ierr); 362674ae819SStefano Zampini 363674ae819SStefano Zampini if (test_fetidp && !pcbddc->use_deluxe_scaling) { 364674ae819SStefano Zampini 3655b08dc53SStefano Zampini PetscReal real_value; 3665b08dc53SStefano Zampini 367674ae819SStefano Zampini ierr = PetscViewerASCIIGetStdout(comm,&viewer);CHKERRQ(ierr); 368674ae819SStefano Zampini ierr = PetscViewerASCIISynchronizedAllow(viewer,PETSC_TRUE);CHKERRQ(ierr); 369674ae819SStefano Zampini ierr = PetscViewerASCIIPrintf(viewer,"----------FETI_DP TESTS--------------\n");CHKERRQ(ierr); 370674ae819SStefano Zampini ierr = PetscViewerASCIIPrintf(viewer,"All tests should return zero!\n");CHKERRQ(ierr); 371674ae819SStefano Zampini ierr = PetscViewerASCIIPrintf(viewer,"FETIDP MAT context in the ");CHKERRQ(ierr); 372674ae819SStefano Zampini if (fully_redundant) { 373674ae819SStefano Zampini ierr = PetscViewerASCIIPrintf(viewer,"fully redundant case for lagrange multipliers.\n");CHKERRQ(ierr); 374674ae819SStefano Zampini } else { 375674ae819SStefano Zampini ierr = PetscViewerASCIIPrintf(viewer,"Non-fully redundant case for lagrange multiplier.\n");CHKERRQ(ierr); 376674ae819SStefano Zampini } 377674ae819SStefano Zampini ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 378674ae819SStefano Zampini 379674ae819SStefano Zampini /******************************************************************/ 380674ae819SStefano Zampini /* TEST A/B: Test numbering of global lambda dofs */ 381674ae819SStefano Zampini /******************************************************************/ 382674ae819SStefano Zampini 383674ae819SStefano Zampini ierr = VecDuplicate(fetidpmat_ctx->lambda_local,&test_vec);CHKERRQ(ierr); 384674ae819SStefano Zampini ierr = VecSet(lambda_global,1.0);CHKERRQ(ierr); 385674ae819SStefano Zampini ierr = VecSet(test_vec,1.0);CHKERRQ(ierr); 386674ae819SStefano Zampini ierr = VecScatterBegin(fetidpmat_ctx->l2g_lambda,lambda_global,fetidpmat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 387674ae819SStefano Zampini ierr = VecScatterEnd (fetidpmat_ctx->l2g_lambda,lambda_global,fetidpmat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 388674ae819SStefano Zampini scalar_value = -1.0; 389674ae819SStefano Zampini ierr = VecAXPY(test_vec,scalar_value,fetidpmat_ctx->lambda_local);CHKERRQ(ierr); 3905b08dc53SStefano Zampini ierr = VecNorm(test_vec,NORM_INFINITY,&real_value);CHKERRQ(ierr); 391674ae819SStefano Zampini ierr = VecDestroy(&test_vec);CHKERRQ(ierr); 3925b08dc53SStefano Zampini ierr = PetscViewerASCIISynchronizedPrintf(viewer,"A[%04d]: CHECK glob to loc: % 1.14e\n",rank,real_value);CHKERRQ(ierr); 393674ae819SStefano Zampini ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 394674ae819SStefano Zampini if (fully_redundant) { 395674ae819SStefano Zampini ierr = VecSet(lambda_global,0.0);CHKERRQ(ierr); 396674ae819SStefano Zampini ierr = VecSet(fetidpmat_ctx->lambda_local,0.5);CHKERRQ(ierr); 397674ae819SStefano Zampini ierr = VecScatterBegin(fetidpmat_ctx->l2g_lambda,fetidpmat_ctx->lambda_local,lambda_global,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 398674ae819SStefano Zampini ierr = VecScatterEnd (fetidpmat_ctx->l2g_lambda,fetidpmat_ctx->lambda_local,lambda_global,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 399674ae819SStefano Zampini ierr = VecSum(lambda_global,&scalar_value);CHKERRQ(ierr); 4005b08dc53SStefano Zampini ierr = PetscViewerASCIISynchronizedPrintf(viewer,"B[%04d]: CHECK loc to glob: % 1.14e\n",rank,PetscRealPart(scalar_value)-fetidpmat_ctx->n_lambda);CHKERRQ(ierr); 401674ae819SStefano Zampini ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 402674ae819SStefano Zampini } 403674ae819SStefano Zampini 404674ae819SStefano Zampini /******************************************************************/ 405674ae819SStefano Zampini /* TEST C: It should holds B_delta*w=0, w\in\widehat{W} */ 406674ae819SStefano Zampini /* This is the meaning of the B matrix */ 407674ae819SStefano Zampini /******************************************************************/ 408674ae819SStefano Zampini 409674ae819SStefano Zampini ierr = VecSetRandom(pcis->vec1_N,NULL);CHKERRQ(ierr); 410674ae819SStefano Zampini ierr = VecSet(pcis->vec1_global,0.0);CHKERRQ(ierr); 411674ae819SStefano Zampini ierr = VecScatterBegin(matis->ctx,pcis->vec1_N,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 412674ae819SStefano Zampini ierr = VecScatterEnd (matis->ctx,pcis->vec1_N,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 413674ae819SStefano Zampini ierr = VecScatterBegin(matis->ctx,pcis->vec1_global,pcis->vec1_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 414674ae819SStefano Zampini ierr = VecScatterEnd (matis->ctx,pcis->vec1_global,pcis->vec1_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 415674ae819SStefano Zampini ierr = VecScatterBegin(pcis->N_to_B,pcis->vec1_N,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 416674ae819SStefano Zampini ierr = VecScatterEnd (pcis->N_to_B,pcis->vec1_N,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 417674ae819SStefano Zampini /* Action of B_delta */ 418674ae819SStefano Zampini ierr = MatMult(fetidpmat_ctx->B_delta,pcis->vec1_B,fetidpmat_ctx->lambda_local);CHKERRQ(ierr); 419674ae819SStefano Zampini ierr = VecSet(lambda_global,0.0);CHKERRQ(ierr); 420674ae819SStefano Zampini ierr = VecScatterBegin(fetidpmat_ctx->l2g_lambda,fetidpmat_ctx->lambda_local,lambda_global,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 421674ae819SStefano Zampini ierr = VecScatterEnd (fetidpmat_ctx->l2g_lambda,fetidpmat_ctx->lambda_local,lambda_global,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 4225b08dc53SStefano Zampini ierr = VecNorm(lambda_global,NORM_INFINITY,&real_value);CHKERRQ(ierr); 4235b08dc53SStefano Zampini ierr = PetscViewerASCIIPrintf(viewer,"C[coll]: CHECK infty norm of B_delta*w (w continuous): % 1.14e\n",real_value);CHKERRQ(ierr); 424674ae819SStefano Zampini ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 425674ae819SStefano Zampini 426674ae819SStefano Zampini /******************************************************************/ 427674ae819SStefano Zampini /* TEST D: It should holds E_Dw = w - P_Dw w\in\widetilde{W} */ 428674ae819SStefano Zampini /* E_D = R_D^TR */ 429674ae819SStefano Zampini /* P_D = B_{D,delta}^T B_{delta} */ 430674ae819SStefano Zampini /* eq.44 Mandel Tezaur and Dohrmann 2005 */ 431674ae819SStefano Zampini /******************************************************************/ 432674ae819SStefano Zampini 433674ae819SStefano Zampini /* compute a random vector in \widetilde{W} */ 434674ae819SStefano Zampini ierr = VecSetRandom(pcis->vec1_N,NULL);CHKERRQ(ierr); 435674ae819SStefano Zampini scalar_value = 0.0; /* set zero at vertices */ 436674ae819SStefano Zampini ierr = VecGetArray(pcis->vec1_N,&array);CHKERRQ(ierr); 437674ae819SStefano Zampini for (i=0;i<n_vertices;i++) { array[vertex_indices[i]]=scalar_value; } 438674ae819SStefano Zampini ierr = VecRestoreArray(pcis->vec1_N,&array);CHKERRQ(ierr); 439674ae819SStefano Zampini /* store w for final comparison */ 440674ae819SStefano Zampini ierr = VecDuplicate(pcis->vec1_B,&test_vec);CHKERRQ(ierr); 441674ae819SStefano Zampini ierr = VecScatterBegin(pcis->N_to_B,pcis->vec1_N,test_vec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 442674ae819SStefano Zampini ierr = VecScatterEnd (pcis->N_to_B,pcis->vec1_N,test_vec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 443674ae819SStefano Zampini 444674ae819SStefano Zampini /* Jump operator P_D : results stored in pcis->vec1_B */ 445674ae819SStefano Zampini 446674ae819SStefano Zampini ierr = VecScatterBegin(pcis->N_to_B,pcis->vec1_N,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 447674ae819SStefano Zampini ierr = VecScatterEnd (pcis->N_to_B,pcis->vec1_N,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 448674ae819SStefano Zampini /* Action of B_delta */ 449674ae819SStefano Zampini ierr = MatMult(fetidpmat_ctx->B_delta,pcis->vec1_B,fetidpmat_ctx->lambda_local);CHKERRQ(ierr); 450674ae819SStefano Zampini ierr = VecSet(lambda_global,0.0);CHKERRQ(ierr); 451674ae819SStefano Zampini ierr = VecScatterBegin(fetidpmat_ctx->l2g_lambda,fetidpmat_ctx->lambda_local,lambda_global,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 452674ae819SStefano Zampini ierr = VecScatterEnd (fetidpmat_ctx->l2g_lambda,fetidpmat_ctx->lambda_local,lambda_global,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 453674ae819SStefano Zampini /* Action of B_Ddelta^T */ 454674ae819SStefano Zampini ierr = VecScatterBegin(fetidpmat_ctx->l2g_lambda,lambda_global,fetidpmat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 455674ae819SStefano Zampini ierr = VecScatterEnd (fetidpmat_ctx->l2g_lambda,lambda_global,fetidpmat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 456674ae819SStefano Zampini ierr = MatMultTranspose(fetidpmat_ctx->B_Ddelta,fetidpmat_ctx->lambda_local,pcis->vec1_B);CHKERRQ(ierr); 457674ae819SStefano Zampini 458674ae819SStefano Zampini /* Average operator E_D : results stored in pcis->vec2_B */ 459674ae819SStefano Zampini ierr = VecScatterBegin(pcis->N_to_B,pcis->vec1_N,pcis->vec2_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 460674ae819SStefano Zampini ierr = VecScatterEnd (pcis->N_to_B,pcis->vec1_N,pcis->vec2_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 461674ae819SStefano Zampini ierr = PCBDDCScalingExtension(fetidpmat_ctx->pc,pcis->vec2_B,pcis->vec1_global);CHKERRQ(ierr); 462674ae819SStefano Zampini ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_global,pcis->vec2_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 463674ae819SStefano Zampini ierr = VecScatterEnd (pcis->global_to_B,pcis->vec1_global,pcis->vec2_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 464674ae819SStefano Zampini 465674ae819SStefano Zampini /* test E_D=I-P_D */ 466674ae819SStefano Zampini scalar_value = 1.0; 467674ae819SStefano Zampini ierr = VecAXPY(pcis->vec1_B,scalar_value,pcis->vec2_B);CHKERRQ(ierr); 468674ae819SStefano Zampini scalar_value = -1.0; 469674ae819SStefano Zampini ierr = VecAXPY(pcis->vec1_B,scalar_value,test_vec);CHKERRQ(ierr); 4705b08dc53SStefano Zampini ierr = VecNorm(pcis->vec1_B,NORM_INFINITY,&real_value);CHKERRQ(ierr); 471674ae819SStefano Zampini ierr = VecDestroy(&test_vec);CHKERRQ(ierr); 4725b08dc53SStefano Zampini ierr = PetscViewerASCIISynchronizedPrintf(viewer,"D[%04d] CHECK infty norm of E_D + P_D - I: % 1.14e\n",rank,real_value);CHKERRQ(ierr); 473674ae819SStefano Zampini ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 474674ae819SStefano Zampini 475674ae819SStefano Zampini /******************************************************************/ 476674ae819SStefano Zampini /* TEST E: It should holds R_D^TP_Dw=0 w\in\widetilde{W} */ 477674ae819SStefano Zampini /* eq.48 Mandel Tezaur and Dohrmann 2005 */ 478674ae819SStefano Zampini /******************************************************************/ 479674ae819SStefano Zampini 480674ae819SStefano Zampini ierr = VecSetRandom(pcis->vec1_N,NULL);CHKERRQ(ierr); 481674ae819SStefano Zampini ierr = VecGetArray(pcis->vec1_N,&array);CHKERRQ(ierr); 482674ae819SStefano Zampini scalar_value = 0.0; /* set zero at vertices */ 483674ae819SStefano Zampini for (i=0;i<n_vertices;i++) { array[vertex_indices[i]]=scalar_value; } 484674ae819SStefano Zampini ierr = VecRestoreArray(pcis->vec1_N,&array);CHKERRQ(ierr); 485674ae819SStefano Zampini 486674ae819SStefano Zampini /* Jump operator P_D : results stored in pcis->vec1_B */ 487674ae819SStefano Zampini 488674ae819SStefano Zampini ierr = VecScatterBegin(pcis->N_to_B,pcis->vec1_N,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 489674ae819SStefano Zampini ierr = VecScatterEnd (pcis->N_to_B,pcis->vec1_N,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 490674ae819SStefano Zampini /* Action of B_delta */ 491674ae819SStefano Zampini ierr = MatMult(fetidpmat_ctx->B_delta,pcis->vec1_B,fetidpmat_ctx->lambda_local);CHKERRQ(ierr); 492674ae819SStefano Zampini ierr = VecSet(lambda_global,0.0);CHKERRQ(ierr); 493674ae819SStefano Zampini ierr = VecScatterBegin(fetidpmat_ctx->l2g_lambda,fetidpmat_ctx->lambda_local,lambda_global,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 494674ae819SStefano Zampini ierr = VecScatterEnd (fetidpmat_ctx->l2g_lambda,fetidpmat_ctx->lambda_local,lambda_global,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 495674ae819SStefano Zampini /* Action of B_Ddelta^T */ 496674ae819SStefano Zampini ierr = VecScatterBegin(fetidpmat_ctx->l2g_lambda,lambda_global,fetidpmat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 497674ae819SStefano Zampini ierr = VecScatterEnd (fetidpmat_ctx->l2g_lambda,lambda_global,fetidpmat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 498674ae819SStefano Zampini ierr = MatMultTranspose(fetidpmat_ctx->B_Ddelta,fetidpmat_ctx->lambda_local,pcis->vec1_B);CHKERRQ(ierr); 499674ae819SStefano Zampini /* scaling */ 500674ae819SStefano Zampini ierr = PCBDDCScalingExtension(fetidpmat_ctx->pc,pcis->vec1_B,pcis->vec1_global);CHKERRQ(ierr); 5015b08dc53SStefano Zampini ierr = VecNorm(pcis->vec1_global,NORM_INFINITY,&real_value);CHKERRQ(ierr); 5025b08dc53SStefano Zampini ierr = PetscViewerASCIIPrintf(viewer,"E[coll]: CHECK infty norm of R^T_D P_D: % 1.14e\n",real_value);CHKERRQ(ierr); 503674ae819SStefano Zampini ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 504674ae819SStefano Zampini 505674ae819SStefano Zampini if (!fully_redundant) { 506674ae819SStefano Zampini /******************************************************************/ 507674ae819SStefano Zampini /* TEST F: It should holds B_{delta}B^T_{D,delta}=I */ 508674ae819SStefano Zampini /* Corollary thm 14 Mandel Tezaur and Dohrmann 2005 */ 509674ae819SStefano Zampini /******************************************************************/ 510674ae819SStefano Zampini ierr = VecDuplicate(lambda_global,&test_vec);CHKERRQ(ierr); 511674ae819SStefano Zampini ierr = VecSetRandom(lambda_global,NULL);CHKERRQ(ierr); 512674ae819SStefano Zampini /* Action of B_Ddelta^T */ 513674ae819SStefano Zampini ierr = VecScatterBegin(fetidpmat_ctx->l2g_lambda,lambda_global,fetidpmat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 514674ae819SStefano Zampini ierr = VecScatterEnd (fetidpmat_ctx->l2g_lambda,lambda_global,fetidpmat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 515674ae819SStefano Zampini ierr = MatMultTranspose(fetidpmat_ctx->B_Ddelta,fetidpmat_ctx->lambda_local,pcis->vec1_B);CHKERRQ(ierr); 516674ae819SStefano Zampini /* Action of B_delta */ 517674ae819SStefano Zampini ierr = MatMult(fetidpmat_ctx->B_delta,pcis->vec1_B,fetidpmat_ctx->lambda_local);CHKERRQ(ierr); 518674ae819SStefano Zampini ierr = VecSet(test_vec,0.0);CHKERRQ(ierr); 519674ae819SStefano Zampini ierr = VecScatterBegin(fetidpmat_ctx->l2g_lambda,fetidpmat_ctx->lambda_local,test_vec,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 520674ae819SStefano Zampini ierr = VecScatterEnd (fetidpmat_ctx->l2g_lambda,fetidpmat_ctx->lambda_local,test_vec,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 521674ae819SStefano Zampini scalar_value = -1.0; 522674ae819SStefano Zampini ierr = VecAXPY(lambda_global,scalar_value,test_vec);CHKERRQ(ierr); 5235b08dc53SStefano Zampini ierr = VecNorm(lambda_global,NORM_INFINITY,&real_value);CHKERRQ(ierr); 5245b08dc53SStefano Zampini ierr = PetscViewerASCIIPrintf(viewer,"E[coll]: CHECK infty norm of P^T_D - I: % 1.14e\n",real_value);CHKERRQ(ierr); 525674ae819SStefano Zampini ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 526674ae819SStefano Zampini ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 527674ae819SStefano Zampini ierr = VecDestroy(&test_vec);CHKERRQ(ierr); 528674ae819SStefano Zampini } 529674ae819SStefano Zampini } 530674ae819SStefano Zampini /* final cleanup */ 531674ae819SStefano Zampini ierr = PetscFree(vertex_indices);CHKERRQ(ierr); 532674ae819SStefano Zampini ierr = VecDestroy(&lambda_global);CHKERRQ(ierr); 533674ae819SStefano Zampini 534674ae819SStefano Zampini PetscFunctionReturn(0); 535674ae819SStefano Zampini } 536674ae819SStefano Zampini 537674ae819SStefano Zampini #undef __FUNCT__ 538674ae819SStefano Zampini #define __FUNCT__ "PCBDDCSetupFETIDPPCContext" 539674ae819SStefano Zampini PetscErrorCode PCBDDCSetupFETIDPPCContext(Mat fetimat, FETIDPPC_ctx fetidppc_ctx) 540674ae819SStefano Zampini { 541674ae819SStefano Zampini FETIDPMat_ctx mat_ctx; 542674ae819SStefano Zampini PetscErrorCode ierr; 543674ae819SStefano Zampini 544674ae819SStefano Zampini PetscFunctionBegin; 545674ae819SStefano Zampini ierr = MatShellGetContext(fetimat,(void**)&mat_ctx);CHKERRQ(ierr); 546674ae819SStefano Zampini /* get references from objects created when setting up feti mat context */ 547674ae819SStefano Zampini ierr = PetscObjectReference((PetscObject)mat_ctx->lambda_local);CHKERRQ(ierr); 548674ae819SStefano Zampini fetidppc_ctx->lambda_local = mat_ctx->lambda_local; 549674ae819SStefano Zampini ierr = PetscObjectReference((PetscObject)mat_ctx->B_Ddelta);CHKERRQ(ierr); 550674ae819SStefano Zampini fetidppc_ctx->B_Ddelta = mat_ctx->B_Ddelta; 551674ae819SStefano Zampini ierr = PetscObjectReference((PetscObject)mat_ctx->l2g_lambda);CHKERRQ(ierr); 552674ae819SStefano Zampini fetidppc_ctx->l2g_lambda = mat_ctx->l2g_lambda; 553674ae819SStefano Zampini PetscFunctionReturn(0); 554674ae819SStefano Zampini } 555674ae819SStefano Zampini 556674ae819SStefano Zampini #undef __FUNCT__ 557674ae819SStefano Zampini #define __FUNCT__ "FETIDPMatMult" 558674ae819SStefano Zampini PetscErrorCode FETIDPMatMult(Mat fetimat, Vec x, Vec y) 559674ae819SStefano Zampini { 560674ae819SStefano Zampini FETIDPMat_ctx mat_ctx; 561674ae819SStefano Zampini PC_IS *pcis; 562674ae819SStefano Zampini PetscErrorCode ierr; 563674ae819SStefano Zampini 564674ae819SStefano Zampini PetscFunctionBegin; 565674ae819SStefano Zampini ierr = MatShellGetContext(fetimat,(void**)&mat_ctx);CHKERRQ(ierr); 566674ae819SStefano Zampini pcis = (PC_IS*)mat_ctx->pc->data; 567674ae819SStefano Zampini /* Application of B_delta^T */ 568674ae819SStefano Zampini ierr = VecScatterBegin(mat_ctx->l2g_lambda,x,mat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 569674ae819SStefano Zampini ierr = VecScatterEnd(mat_ctx->l2g_lambda,x,mat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 570674ae819SStefano Zampini ierr = MatMultTranspose(mat_ctx->B_delta,mat_ctx->lambda_local,pcis->vec1_B);CHKERRQ(ierr); 571674ae819SStefano Zampini /* Application of \widetilde{S}^-1 */ 572674ae819SStefano Zampini ierr = VecSet(pcis->vec1_D,0.0);CHKERRQ(ierr); 573674ae819SStefano Zampini ierr = PCBDDCApplyInterfacePreconditioner(mat_ctx->pc);CHKERRQ(ierr); 574674ae819SStefano Zampini /* Application of B_delta */ 575674ae819SStefano Zampini ierr = MatMult(mat_ctx->B_delta,pcis->vec1_B,mat_ctx->lambda_local);CHKERRQ(ierr); 576674ae819SStefano Zampini ierr = VecSet(y,0.0);CHKERRQ(ierr); 577674ae819SStefano Zampini ierr = VecScatterBegin(mat_ctx->l2g_lambda,mat_ctx->lambda_local,y,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 578674ae819SStefano Zampini ierr = VecScatterEnd(mat_ctx->l2g_lambda,mat_ctx->lambda_local,y,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 579674ae819SStefano Zampini PetscFunctionReturn(0); 580674ae819SStefano Zampini } 581674ae819SStefano Zampini 582674ae819SStefano Zampini #undef __FUNCT__ 583674ae819SStefano Zampini #define __FUNCT__ "FETIDPPCApply" 584674ae819SStefano Zampini PetscErrorCode FETIDPPCApply(PC fetipc, Vec x, Vec y) 585674ae819SStefano Zampini { 586674ae819SStefano Zampini FETIDPPC_ctx pc_ctx; 587674ae819SStefano Zampini PC_IS *pcis; 588674ae819SStefano Zampini PetscErrorCode ierr; 589674ae819SStefano Zampini 590674ae819SStefano Zampini PetscFunctionBegin; 591674ae819SStefano Zampini ierr = PCShellGetContext(fetipc,(void**)&pc_ctx); 592674ae819SStefano Zampini pcis = (PC_IS*)pc_ctx->pc->data; 593674ae819SStefano Zampini /* Application of B_Ddelta^T */ 594674ae819SStefano Zampini ierr = VecScatterBegin(pc_ctx->l2g_lambda,x,pc_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 595674ae819SStefano Zampini ierr = VecScatterEnd(pc_ctx->l2g_lambda,x,pc_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 596674ae819SStefano Zampini ierr = VecSet(pcis->vec2_B,0.0);CHKERRQ(ierr); 597674ae819SStefano Zampini ierr = MatMultTranspose(pc_ctx->B_Ddelta,pc_ctx->lambda_local,pcis->vec2_B);CHKERRQ(ierr); 598674ae819SStefano Zampini /* Application of S */ 599674ae819SStefano Zampini ierr = PCISApplySchur(pc_ctx->pc,pcis->vec2_B,pcis->vec1_B,(Vec)0,pcis->vec1_D,pcis->vec2_D);CHKERRQ(ierr); 600674ae819SStefano Zampini /* Application of B_Ddelta */ 601674ae819SStefano Zampini ierr = MatMult(pc_ctx->B_Ddelta,pcis->vec1_B,pc_ctx->lambda_local);CHKERRQ(ierr); 602674ae819SStefano Zampini ierr = VecSet(y,0.0);CHKERRQ(ierr); 603674ae819SStefano Zampini ierr = VecScatterBegin(pc_ctx->l2g_lambda,pc_ctx->lambda_local,y,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 604674ae819SStefano Zampini ierr = VecScatterEnd(pc_ctx->l2g_lambda,pc_ctx->lambda_local,y,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 605674ae819SStefano Zampini PetscFunctionReturn(0); 606674ae819SStefano Zampini } 607