1ab5c6b0cSJed Brown #include <../src/ksp/pc/impls/bddc/bddc.h> 2ab5c6b0cSJed Brown #include <../src/ksp/pc/impls/bddc/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; 12854ce69bSBarry Smith ierr = PetscNew(&newctx);CHKERRQ(ierr); 13674ae819SStefano Zampini /* increase the reference count for BDDC preconditioner */ 14674ae819SStefano Zampini ierr = PetscObjectReference((PetscObject)pc);CHKERRQ(ierr); 15674ae819SStefano Zampini newctx->pc = pc; 16674ae819SStefano Zampini *fetidpmat_ctx = newctx; 17674ae819SStefano Zampini PetscFunctionReturn(0); 18674ae819SStefano Zampini } 19674ae819SStefano Zampini 20674ae819SStefano Zampini #undef __FUNCT__ 21674ae819SStefano Zampini #define __FUNCT__ "PCBDDCCreateFETIDPPCContext" 22674ae819SStefano Zampini PetscErrorCode PCBDDCCreateFETIDPPCContext(PC pc, FETIDPPC_ctx *fetidppc_ctx) 23674ae819SStefano Zampini { 24674ae819SStefano Zampini FETIDPPC_ctx newctx; 25674ae819SStefano Zampini PetscErrorCode ierr; 26674ae819SStefano Zampini 27674ae819SStefano Zampini PetscFunctionBegin; 28854ce69bSBarry Smith ierr = PetscNew(&newctx);CHKERRQ(ierr); 29674ae819SStefano Zampini /* increase the reference count for BDDC preconditioner */ 30674ae819SStefano Zampini ierr = PetscObjectReference((PetscObject)pc);CHKERRQ(ierr); 31674ae819SStefano Zampini newctx->pc = pc; 32674ae819SStefano Zampini *fetidppc_ctx = newctx; 33674ae819SStefano Zampini PetscFunctionReturn(0); 34674ae819SStefano Zampini } 35674ae819SStefano Zampini 36674ae819SStefano Zampini #undef __FUNCT__ 37674ae819SStefano Zampini #define __FUNCT__ "PCBDDCDestroyFETIDPMat" 38674ae819SStefano Zampini PetscErrorCode PCBDDCDestroyFETIDPMat(Mat A) 39674ae819SStefano Zampini { 40674ae819SStefano Zampini FETIDPMat_ctx mat_ctx; 41674ae819SStefano Zampini PetscErrorCode ierr; 42674ae819SStefano Zampini 43674ae819SStefano Zampini PetscFunctionBegin; 44674ae819SStefano Zampini ierr = MatShellGetContext(A,(void**)&mat_ctx);CHKERRQ(ierr); 45674ae819SStefano Zampini ierr = VecDestroy(&mat_ctx->lambda_local);CHKERRQ(ierr); 46674ae819SStefano Zampini ierr = VecDestroy(&mat_ctx->temp_solution_D);CHKERRQ(ierr); 47674ae819SStefano Zampini ierr = VecDestroy(&mat_ctx->temp_solution_B);CHKERRQ(ierr); 48674ae819SStefano Zampini ierr = MatDestroy(&mat_ctx->B_delta);CHKERRQ(ierr); 49674ae819SStefano Zampini ierr = MatDestroy(&mat_ctx->B_Ddelta);CHKERRQ(ierr); 50*457d4a33Sstefano_zampini ierr = MatDestroy(&mat_ctx->B_BB);CHKERRQ(ierr); 51*457d4a33Sstefano_zampini ierr = MatDestroy(&mat_ctx->B_BI);CHKERRQ(ierr); 52*457d4a33Sstefano_zampini ierr = MatDestroy(&mat_ctx->Bt_BB);CHKERRQ(ierr); 53*457d4a33Sstefano_zampini ierr = MatDestroy(&mat_ctx->Bt_BI);CHKERRQ(ierr); 54*457d4a33Sstefano_zampini ierr = MatDestroy(&mat_ctx->C);CHKERRQ(ierr); 55*457d4a33Sstefano_zampini ierr = VecDestroy(&mat_ctx->vP);CHKERRQ(ierr); 56*457d4a33Sstefano_zampini ierr = VecDestroy(&mat_ctx->xPg);CHKERRQ(ierr); 57*457d4a33Sstefano_zampini ierr = VecDestroy(&mat_ctx->yPg);CHKERRQ(ierr); 58*457d4a33Sstefano_zampini ierr = VecDestroy(&mat_ctx->work);CHKERRQ(ierr); 59674ae819SStefano Zampini ierr = VecScatterDestroy(&mat_ctx->l2g_lambda);CHKERRQ(ierr); 60*457d4a33Sstefano_zampini ierr = VecScatterDestroy(&mat_ctx->l2g_p);CHKERRQ(ierr); 61*457d4a33Sstefano_zampini ierr = VecScatterDestroy(&mat_ctx->g2g);CHKERRQ(ierr); 62*457d4a33Sstefano_zampini ierr = VecScatterDestroy(&mat_ctx->g2l_p);CHKERRQ(ierr); 63674ae819SStefano Zampini ierr = PCDestroy(&mat_ctx->pc);CHKERRQ(ierr); /* decrease PCBDDC reference count */ 64674ae819SStefano Zampini ierr = PetscFree(mat_ctx);CHKERRQ(ierr); 65674ae819SStefano Zampini PetscFunctionReturn(0); 66674ae819SStefano Zampini } 67674ae819SStefano Zampini 68674ae819SStefano Zampini #undef __FUNCT__ 69674ae819SStefano Zampini #define __FUNCT__ "PCBDDCDestroyFETIDPPC" 70674ae819SStefano Zampini PetscErrorCode PCBDDCDestroyFETIDPPC(PC pc) 71674ae819SStefano Zampini { 72674ae819SStefano Zampini FETIDPPC_ctx pc_ctx; 73674ae819SStefano Zampini PetscErrorCode ierr; 74674ae819SStefano Zampini 75674ae819SStefano Zampini PetscFunctionBegin; 76674ae819SStefano Zampini ierr = PCShellGetContext(pc,(void**)&pc_ctx);CHKERRQ(ierr); 77674ae819SStefano Zampini ierr = VecDestroy(&pc_ctx->lambda_local);CHKERRQ(ierr); 78674ae819SStefano Zampini ierr = MatDestroy(&pc_ctx->B_Ddelta);CHKERRQ(ierr); 79674ae819SStefano Zampini ierr = VecScatterDestroy(&pc_ctx->l2g_lambda);CHKERRQ(ierr); 8068668abeSStefano Zampini ierr = MatDestroy(&pc_ctx->S_j);CHKERRQ(ierr); 81674ae819SStefano Zampini ierr = PCDestroy(&pc_ctx->pc);CHKERRQ(ierr); /* decrease PCBDDC reference count */ 82*457d4a33Sstefano_zampini ierr = KSPDestroy(&pc_ctx->kP);CHKERRQ(ierr); 83*457d4a33Sstefano_zampini ierr = VecDestroy(&pc_ctx->xPg);CHKERRQ(ierr); 84*457d4a33Sstefano_zampini ierr = VecDestroy(&pc_ctx->yPg);CHKERRQ(ierr); 85674ae819SStefano Zampini ierr = PetscFree(pc_ctx);CHKERRQ(ierr); 86674ae819SStefano Zampini PetscFunctionReturn(0); 87674ae819SStefano Zampini } 88674ae819SStefano Zampini 89674ae819SStefano Zampini #undef __FUNCT__ 90674ae819SStefano Zampini #define __FUNCT__ "PCBDDCSetupFETIDPMatContext" 91674ae819SStefano Zampini PetscErrorCode PCBDDCSetupFETIDPMatContext(FETIDPMat_ctx fetidpmat_ctx ) 92674ae819SStefano Zampini { 93674ae819SStefano Zampini PetscErrorCode ierr; 94674ae819SStefano Zampini PC_IS *pcis=(PC_IS*)fetidpmat_ctx->pc->data; 95674ae819SStefano Zampini PC_BDDC *pcbddc=(PC_BDDC*)fetidpmat_ctx->pc->data; 96674ae819SStefano Zampini PCBDDCGraph mat_graph=pcbddc->mat_graph; 97674ae819SStefano Zampini Mat_IS *matis = (Mat_IS*)fetidpmat_ctx->pc->pmat->data; 98674ae819SStefano Zampini MPI_Comm comm; 99674ae819SStefano Zampini Mat ScalingMat; 100*457d4a33Sstefano_zampini Vec fetidp_global; 101674ae819SStefano Zampini IS IS_l2g_lambda; 102dc456d91SStefano Zampini IS subset,subset_mult,subset_n; 103674ae819SStefano Zampini PetscBool skip_node,fully_redundant; 104674ae819SStefano Zampini PetscInt i,j,k,s,n_boundary_dofs,n_global_lambda,n_vertices,partial_sum; 105dc456d91SStefano Zampini PetscInt cum,n_local_lambda,n_lambda_for_dof,dual_size,n_neg_values,n_pos_values; 10676ec1555SBarry Smith PetscMPIInt rank,size,buf_size,neigh; 107674ae819SStefano Zampini PetscScalar scalar_value; 108674ae819SStefano Zampini PetscInt *vertex_indices; 109dc456d91SStefano Zampini PetscInt *dual_dofs_boundary_indices,*aux_local_numbering_1; 110dc456d91SStefano Zampini const PetscInt *aux_global_numbering; 111674ae819SStefano Zampini PetscInt *aux_sums,*cols_B_delta,*l2g_indices; 112674ae819SStefano Zampini PetscScalar *array,*scaling_factors,*vals_B_delta; 113674ae819SStefano Zampini PetscInt *aux_local_numbering_2; 114*457d4a33Sstefano_zampini PetscLayout llay; 115674ae819SStefano Zampini /* For communication of scaling factors */ 116674ae819SStefano Zampini PetscInt *ptrs_buffer,neigh_position; 117674ae819SStefano Zampini PetscScalar **all_factors,*send_buffer,*recv_buffer; 118674ae819SStefano Zampini MPI_Request *send_reqs,*recv_reqs; 119674ae819SStefano Zampini /* tests */ 120674ae819SStefano Zampini Vec test_vec; 121674ae819SStefano Zampini PetscBool test_fetidp; 122674ae819SStefano Zampini PetscViewer viewer; 123*457d4a33Sstefano_zampini /* saddlepoint */ 124*457d4a33Sstefano_zampini Mat oA; 125*457d4a33Sstefano_zampini ISLocalToGlobalMapping l2gmap_p; 126*457d4a33Sstefano_zampini PetscLayout play; 127*457d4a33Sstefano_zampini IS gP,pP; 128*457d4a33Sstefano_zampini PetscInt nPl,nPg,nPgl; 129674ae819SStefano Zampini 130674ae819SStefano Zampini PetscFunctionBegin; 131674ae819SStefano Zampini ierr = PetscObjectGetComm((PetscObject)(fetidpmat_ctx->pc),&comm);CHKERRQ(ierr); 132674ae819SStefano Zampini ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 13376ec1555SBarry Smith ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 134674ae819SStefano Zampini 135*457d4a33Sstefano_zampini /* saddlepoint */ 136*457d4a33Sstefano_zampini nPl = 0; 137*457d4a33Sstefano_zampini nPg = 0; 138*457d4a33Sstefano_zampini nPgl = 0; 139*457d4a33Sstefano_zampini gP = NULL; 140*457d4a33Sstefano_zampini pP = NULL; 141*457d4a33Sstefano_zampini l2gmap_p = NULL; 142*457d4a33Sstefano_zampini play = NULL; 143*457d4a33Sstefano_zampini oA = NULL; 144*457d4a33Sstefano_zampini ierr = PetscObjectQuery((PetscObject)fetidpmat_ctx->pc,"__KSPFETIDP_A",(PetscObject*)&oA);CHKERRQ(ierr); 145*457d4a33Sstefano_zampini if (oA) { /* oA is the original matrix */ 146*457d4a33Sstefano_zampini /* pressures in parallel layout of oA */ 147*457d4a33Sstefano_zampini ierr = PetscObjectQuery((PetscObject)fetidpmat_ctx->pc,"__KSPFETIDP_pP",(PetscObject*)&pP);CHKERRQ(ierr); 148*457d4a33Sstefano_zampini if (!pP) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"pP not present"); 149*457d4a33Sstefano_zampini 150*457d4a33Sstefano_zampini /* subdomain pressures in global numbering */ 151*457d4a33Sstefano_zampini ierr = PetscObjectQuery((PetscObject)fetidpmat_ctx->pc,"__KSPFETIDP_gP",(PetscObject*)&gP);CHKERRQ(ierr); 152*457d4a33Sstefano_zampini if (!gP) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"gP not present"); 153*457d4a33Sstefano_zampini ierr = ISGetLocalSize(gP,&nPl);CHKERRQ(ierr); 154*457d4a33Sstefano_zampini ierr = VecCreate(PETSC_COMM_SELF,&fetidpmat_ctx->vP);CHKERRQ(ierr); 155*457d4a33Sstefano_zampini ierr = VecSetSizes(fetidpmat_ctx->vP,nPl,nPl);CHKERRQ(ierr); 156*457d4a33Sstefano_zampini ierr = VecSetType(fetidpmat_ctx->vP,VECSTANDARD);CHKERRQ(ierr); 157*457d4a33Sstefano_zampini ierr = VecSetUp(fetidpmat_ctx->vP);CHKERRQ(ierr); 158*457d4a33Sstefano_zampini 159*457d4a33Sstefano_zampini /* interface pressure matrix */ 160*457d4a33Sstefano_zampini ierr = PetscObjectQuery((PetscObject)fetidpmat_ctx->pc,"__KSPFETIDP_C",(PetscObject*)&fetidpmat_ctx->C);CHKERRQ(ierr); 161*457d4a33Sstefano_zampini if (!fetidpmat_ctx->C) { /* null pressure block, compute layout and global numbering for interface pressures */ 162*457d4a33Sstefano_zampini IS Pg; 163*457d4a33Sstefano_zampini 164*457d4a33Sstefano_zampini ierr = ISRenumber(gP,NULL,&nPg,&Pg);CHKERRQ(ierr); 165*457d4a33Sstefano_zampini ierr = ISLocalToGlobalMappingCreateIS(Pg,&l2gmap_p);CHKERRQ(ierr); 166*457d4a33Sstefano_zampini ierr = ISDestroy(&Pg);CHKERRQ(ierr); 167*457d4a33Sstefano_zampini ierr = PetscLayoutCreate(comm,&play);CHKERRQ(ierr); 168*457d4a33Sstefano_zampini ierr = PetscLayoutSetBlockSize(play,1);CHKERRQ(ierr); 169*457d4a33Sstefano_zampini ierr = PetscLayoutSetSize(play,nPg);CHKERRQ(ierr); 170*457d4a33Sstefano_zampini ierr = ISGetLocalSize(pP,&nPgl);CHKERRQ(ierr); 171*457d4a33Sstefano_zampini ierr = PetscLayoutSetLocalSize(play,nPgl);CHKERRQ(ierr); 172*457d4a33Sstefano_zampini ierr = PetscLayoutSetUp(play);CHKERRQ(ierr); 173*457d4a33Sstefano_zampini } else { 174*457d4a33Sstefano_zampini ierr = PetscObjectReference((PetscObject)fetidpmat_ctx->C);CHKERRQ(ierr); 175*457d4a33Sstefano_zampini ierr = MatGetLocalToGlobalMapping(fetidpmat_ctx->C,&l2gmap_p,NULL);CHKERRQ(ierr); 176*457d4a33Sstefano_zampini ierr = PetscObjectReference((PetscObject)l2gmap_p);CHKERRQ(ierr); 177*457d4a33Sstefano_zampini ierr = MatGetSize(fetidpmat_ctx->C,&nPg,NULL);CHKERRQ(ierr); 178*457d4a33Sstefano_zampini ierr = MatGetLocalSize(fetidpmat_ctx->C,NULL,&nPgl);CHKERRQ(ierr); 179*457d4a33Sstefano_zampini ierr = MatGetLayouts(fetidpmat_ctx->C,NULL,&llay);CHKERRQ(ierr); 180*457d4a33Sstefano_zampini ierr = PetscLayoutReference(llay,&play);CHKERRQ(ierr); 181*457d4a33Sstefano_zampini } 182*457d4a33Sstefano_zampini ierr = VecCreateMPIWithArray(comm,1,nPgl,nPg,NULL,&fetidpmat_ctx->xPg);CHKERRQ(ierr); 183*457d4a33Sstefano_zampini ierr = VecCreateMPIWithArray(comm,1,nPgl,nPg,NULL,&fetidpmat_ctx->yPg);CHKERRQ(ierr); 184*457d4a33Sstefano_zampini 185*457d4a33Sstefano_zampini /* import matrices for interface pressures coupling */ 186*457d4a33Sstefano_zampini ierr = PetscObjectQuery((PetscObject)fetidpmat_ctx->pc,"__KSPFETIDP_B_BI",(PetscObject*)&fetidpmat_ctx->B_BI);CHKERRQ(ierr); 187*457d4a33Sstefano_zampini if (!fetidpmat_ctx->B_BI) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"B_BI not present"); 188*457d4a33Sstefano_zampini ierr = PetscObjectReference((PetscObject)fetidpmat_ctx->B_BI);CHKERRQ(ierr); 189*457d4a33Sstefano_zampini 190*457d4a33Sstefano_zampini ierr = PetscObjectQuery((PetscObject)fetidpmat_ctx->pc,"__KSPFETIDP_B_BB",(PetscObject*)&fetidpmat_ctx->B_BB);CHKERRQ(ierr); 191*457d4a33Sstefano_zampini if (!fetidpmat_ctx->B_BB) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"B_BB not present"); 192*457d4a33Sstefano_zampini ierr = PetscObjectReference((PetscObject)fetidpmat_ctx->B_BB);CHKERRQ(ierr); 193*457d4a33Sstefano_zampini 194*457d4a33Sstefano_zampini ierr = PetscObjectQuery((PetscObject)fetidpmat_ctx->pc,"__KSPFETIDP_Bt_BI",(PetscObject*)&fetidpmat_ctx->Bt_BI);CHKERRQ(ierr); 195*457d4a33Sstefano_zampini if (!fetidpmat_ctx->Bt_BI) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Bt_BI not present"); 196*457d4a33Sstefano_zampini ierr = PetscObjectReference((PetscObject)fetidpmat_ctx->Bt_BI);CHKERRQ(ierr); 197*457d4a33Sstefano_zampini 198*457d4a33Sstefano_zampini ierr = PetscObjectQuery((PetscObject)fetidpmat_ctx->pc,"__KSPFETIDP_Bt_BB",(PetscObject*)&fetidpmat_ctx->Bt_BB);CHKERRQ(ierr); 199*457d4a33Sstefano_zampini if (!fetidpmat_ctx->Bt_BB) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Bt_BB not present"); 200*457d4a33Sstefano_zampini ierr = PetscObjectReference((PetscObject)fetidpmat_ctx->Bt_BB);CHKERRQ(ierr); 201*457d4a33Sstefano_zampini } 202*457d4a33Sstefano_zampini 203674ae819SStefano Zampini /* Default type of lagrange multipliers is non-redundant */ 204329cd39dSStefano Zampini fully_redundant = fetidpmat_ctx->fully_redundant; 205674ae819SStefano Zampini 206674ae819SStefano Zampini /* Evaluate local and global number of lagrange multipliers */ 207674ae819SStefano Zampini ierr = VecSet(pcis->vec1_N,0.0);CHKERRQ(ierr); 208674ae819SStefano Zampini n_local_lambda = 0; 209674ae819SStefano Zampini partial_sum = 0; 210674ae819SStefano Zampini n_boundary_dofs = 0; 211674ae819SStefano Zampini s = 0; 212674ae819SStefano Zampini /* Get Vertices used to define the BDDC */ 213b371cd4fSStefano Zampini n_vertices = pcbddc->n_vertices; 2140e6343abSStefano Zampini vertex_indices = pcbddc->local_primal_ref_node; 215674ae819SStefano Zampini dual_size = pcis->n_B-n_vertices; 216785e854fSJed Brown ierr = PetscMalloc1(dual_size,&dual_dofs_boundary_indices);CHKERRQ(ierr); 217785e854fSJed Brown ierr = PetscMalloc1(dual_size,&aux_local_numbering_1);CHKERRQ(ierr); 218785e854fSJed Brown ierr = PetscMalloc1(dual_size,&aux_local_numbering_2);CHKERRQ(ierr); 219674ae819SStefano Zampini 220674ae819SStefano Zampini ierr = VecGetArray(pcis->vec1_N,&array);CHKERRQ(ierr); 221674ae819SStefano Zampini for (i=0;i<pcis->n;i++){ 222674ae819SStefano Zampini j = mat_graph->count[i]; /* RECALL: mat_graph->count[i] does not count myself */ 223674ae819SStefano Zampini if ( j > 0 ) { 224674ae819SStefano Zampini n_boundary_dofs++; 225674ae819SStefano Zampini } 226674ae819SStefano Zampini skip_node = PETSC_FALSE; 227674ae819SStefano Zampini if ( s < n_vertices && vertex_indices[s]==i) { /* it works for a sorted set of vertices */ 228674ae819SStefano Zampini skip_node = PETSC_TRUE; 229674ae819SStefano Zampini s++; 230674ae819SStefano Zampini } 231674ae819SStefano Zampini if (j < 1) { 232674ae819SStefano Zampini skip_node = PETSC_TRUE; 233674ae819SStefano Zampini } 234674ae819SStefano Zampini if ( !skip_node ) { 235674ae819SStefano Zampini if (fully_redundant) { 236674ae819SStefano Zampini /* fully redundant set of lagrange multipliers */ 237674ae819SStefano Zampini n_lambda_for_dof = (j*(j+1))/2; 238674ae819SStefano Zampini } else { 239674ae819SStefano Zampini n_lambda_for_dof = j; 240674ae819SStefano Zampini } 241674ae819SStefano Zampini n_local_lambda += j; 242674ae819SStefano Zampini /* needed to evaluate global number of lagrange multipliers */ 243674ae819SStefano Zampini array[i]=(1.0*n_lambda_for_dof)/(j+1.0); /* already scaled for the next global sum */ 244674ae819SStefano Zampini /* store some data needed */ 245674ae819SStefano Zampini dual_dofs_boundary_indices[partial_sum] = n_boundary_dofs-1; 246674ae819SStefano Zampini aux_local_numbering_1[partial_sum] = i; 247674ae819SStefano Zampini aux_local_numbering_2[partial_sum] = n_lambda_for_dof; 248674ae819SStefano Zampini partial_sum++; 249674ae819SStefano Zampini } 250674ae819SStefano Zampini } 251674ae819SStefano Zampini ierr = VecRestoreArray(pcis->vec1_N,&array);CHKERRQ(ierr); 252674ae819SStefano Zampini 253674ae819SStefano Zampini /* compute global ordering of lagrange multipliers and associate l2g map */ 254dc456d91SStefano Zampini ierr = ISCreateGeneral(comm,partial_sum,aux_local_numbering_1,PETSC_COPY_VALUES,&subset_n);CHKERRQ(ierr); 2553bbff08aSStefano Zampini ierr = ISLocalToGlobalMappingApplyIS(pcis->mapping,subset_n,&subset);CHKERRQ(ierr); 256dc456d91SStefano Zampini ierr = ISDestroy(&subset_n);CHKERRQ(ierr); 257dc456d91SStefano Zampini ierr = ISCreateGeneral(comm,partial_sum,aux_local_numbering_2,PETSC_OWN_POINTER,&subset_mult);CHKERRQ(ierr); 2583d996552SStefano Zampini ierr = ISRenumber(subset,subset_mult,&fetidpmat_ctx->n_lambda,&subset_n);CHKERRQ(ierr); 259dc456d91SStefano Zampini ierr = ISDestroy(&subset);CHKERRQ(ierr); 2603d996552SStefano Zampini 2614fe826edSStefano Zampini #if defined(PETSC_USE_DEBUG) 2624fe826edSStefano Zampini ierr = VecSet(pcis->vec1_global,0.0);CHKERRQ(ierr); 2634fe826edSStefano Zampini ierr = VecScatterBegin(matis->rctx,pcis->vec1_N,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 2644fe826edSStefano Zampini ierr = VecScatterEnd(matis->rctx,pcis->vec1_N,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 2654fe826edSStefano Zampini ierr = VecSum(pcis->vec1_global,&scalar_value);CHKERRQ(ierr); 2664fe826edSStefano Zampini i = (PetscInt)PetscRealPart(scalar_value); 2676c4ed002SBarry Smith 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); 2684fe826edSStefano Zampini #endif 269674ae819SStefano Zampini 270674ae819SStefano Zampini /* init data for scaling factors exchange */ 271674ae819SStefano Zampini partial_sum = 0; 272785e854fSJed Brown ierr = PetscMalloc1(pcis->n_neigh,&ptrs_buffer);CHKERRQ(ierr); 2734b2aedd3SStefano Zampini ierr = PetscMalloc1(PetscMax(pcis->n_neigh-1,0),&send_reqs);CHKERRQ(ierr); 2744b2aedd3SStefano Zampini ierr = PetscMalloc1(PetscMax(pcis->n_neigh-1,0),&recv_reqs);CHKERRQ(ierr); 275785e854fSJed Brown ierr = PetscMalloc1(pcis->n,&all_factors);CHKERRQ(ierr); 2764b2aedd3SStefano Zampini if (pcis->n_neigh > 0) ptrs_buffer[0]=0; 277674ae819SStefano Zampini for (i=1;i<pcis->n_neigh;i++) { 278674ae819SStefano Zampini partial_sum += pcis->n_shared[i]; 279674ae819SStefano Zampini ptrs_buffer[i] = ptrs_buffer[i-1]+pcis->n_shared[i]; 280674ae819SStefano Zampini } 281785e854fSJed Brown ierr = PetscMalloc1(partial_sum,&send_buffer);CHKERRQ(ierr); 282785e854fSJed Brown ierr = PetscMalloc1(partial_sum,&recv_buffer);CHKERRQ(ierr); 283785e854fSJed Brown ierr = PetscMalloc1(partial_sum,&all_factors[0]);CHKERRQ(ierr); 284674ae819SStefano Zampini for (i=0;i<pcis->n-1;i++) { 285674ae819SStefano Zampini j = mat_graph->count[i]; 286674ae819SStefano Zampini all_factors[i+1]=all_factors[i]+j; 287674ae819SStefano Zampini } 288674ae819SStefano Zampini /* scatter B scaling to N vec */ 289674ae819SStefano Zampini ierr = VecScatterBegin(pcis->N_to_B,pcis->D,pcis->vec1_N,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 290674ae819SStefano Zampini ierr = VecScatterEnd(pcis->N_to_B,pcis->D,pcis->vec1_N,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 291674ae819SStefano Zampini /* communications */ 2922b095fd8SStefano Zampini ierr = VecGetArrayRead(pcis->vec1_N,(const PetscScalar**)&array);CHKERRQ(ierr); 293674ae819SStefano Zampini for (i=1;i<pcis->n_neigh;i++) { 294674ae819SStefano Zampini for (j=0;j<pcis->n_shared[i];j++) { 295674ae819SStefano Zampini send_buffer[ptrs_buffer[i-1]+j]=array[pcis->shared[i][j]]; 296674ae819SStefano Zampini } 297674ae819SStefano Zampini ierr = PetscMPIIntCast(ptrs_buffer[i]-ptrs_buffer[i-1],&buf_size);CHKERRQ(ierr); 298674ae819SStefano Zampini ierr = PetscMPIIntCast(pcis->neigh[i],&neigh);CHKERRQ(ierr); 299674ae819SStefano Zampini ierr = MPI_Isend(&send_buffer[ptrs_buffer[i-1]],buf_size,MPIU_SCALAR,neigh,0,comm,&send_reqs[i-1]);CHKERRQ(ierr); 300674ae819SStefano Zampini ierr = MPI_Irecv(&recv_buffer[ptrs_buffer[i-1]],buf_size,MPIU_SCALAR,neigh,0,comm,&recv_reqs[i-1]);CHKERRQ(ierr); 301674ae819SStefano Zampini } 3022b095fd8SStefano Zampini ierr = VecRestoreArrayRead(pcis->vec1_N,(const PetscScalar**)&array);CHKERRQ(ierr); 3034b2aedd3SStefano Zampini if (pcis->n_neigh > 0) { 3044b2aedd3SStefano Zampini ierr = MPI_Waitall(pcis->n_neigh-1,recv_reqs,MPI_STATUSES_IGNORE);CHKERRQ(ierr); 3054b2aedd3SStefano Zampini } 306674ae819SStefano Zampini /* put values in correct places */ 307674ae819SStefano Zampini for (i=1;i<pcis->n_neigh;i++) { 308674ae819SStefano Zampini for (j=0;j<pcis->n_shared[i];j++) { 309674ae819SStefano Zampini k = pcis->shared[i][j]; 310674ae819SStefano Zampini neigh_position = 0; 311674ae819SStefano Zampini while(mat_graph->neighbours_set[k][neigh_position] != pcis->neigh[i]) {neigh_position++;} 312674ae819SStefano Zampini all_factors[k][neigh_position]=recv_buffer[ptrs_buffer[i-1]+j]; 313674ae819SStefano Zampini } 314674ae819SStefano Zampini } 3154b2aedd3SStefano Zampini if (pcis->n_neigh > 0) { 3164b2aedd3SStefano Zampini ierr = MPI_Waitall(pcis->n_neigh-1,send_reqs,MPI_STATUSES_IGNORE);CHKERRQ(ierr); 3174b2aedd3SStefano Zampini } 318674ae819SStefano Zampini ierr = PetscFree(send_reqs);CHKERRQ(ierr); 319674ae819SStefano Zampini ierr = PetscFree(recv_reqs);CHKERRQ(ierr); 320674ae819SStefano Zampini ierr = PetscFree(send_buffer);CHKERRQ(ierr); 321674ae819SStefano Zampini ierr = PetscFree(recv_buffer);CHKERRQ(ierr); 322674ae819SStefano Zampini ierr = PetscFree(ptrs_buffer);CHKERRQ(ierr); 323674ae819SStefano Zampini 324674ae819SStefano Zampini /* Compute B and B_delta (local actions) */ 325785e854fSJed Brown ierr = PetscMalloc1(pcis->n_neigh,&aux_sums);CHKERRQ(ierr); 326785e854fSJed Brown ierr = PetscMalloc1(n_local_lambda,&l2g_indices);CHKERRQ(ierr); 327785e854fSJed Brown ierr = PetscMalloc1(n_local_lambda,&vals_B_delta);CHKERRQ(ierr); 328785e854fSJed Brown ierr = PetscMalloc1(n_local_lambda,&cols_B_delta);CHKERRQ(ierr); 329785e854fSJed Brown ierr = PetscMalloc1(n_local_lambda,&scaling_factors);CHKERRQ(ierr); 330dc456d91SStefano Zampini ierr = ISGetIndices(subset_n,&aux_global_numbering);CHKERRQ(ierr); 331674ae819SStefano Zampini partial_sum=0; 332dc456d91SStefano Zampini cum = 0; 333674ae819SStefano Zampini for (i=0;i<dual_size;i++) { 334dc456d91SStefano Zampini n_global_lambda = aux_global_numbering[cum]; 335674ae819SStefano Zampini j = mat_graph->count[aux_local_numbering_1[i]]; 336674ae819SStefano Zampini aux_sums[0]=0; 337674ae819SStefano Zampini for (s=1;s<j;s++) { 338674ae819SStefano Zampini aux_sums[s]=aux_sums[s-1]+j-s+1; 339674ae819SStefano Zampini } 340674ae819SStefano Zampini array = all_factors[aux_local_numbering_1[i]]; 341674ae819SStefano Zampini n_neg_values = 0; 3422a7da448SStefano Zampini while(n_neg_values < j && mat_graph->neighbours_set[aux_local_numbering_1[i]][n_neg_values] < rank) {n_neg_values++;} 343674ae819SStefano Zampini n_pos_values = j - n_neg_values; 344674ae819SStefano Zampini if (fully_redundant) { 345674ae819SStefano Zampini for (s=0;s<n_neg_values;s++) { 346674ae819SStefano Zampini l2g_indices [partial_sum+s]=aux_sums[s]+n_neg_values-s-1+n_global_lambda; 347674ae819SStefano Zampini cols_B_delta [partial_sum+s]=dual_dofs_boundary_indices[i]; 348674ae819SStefano Zampini vals_B_delta [partial_sum+s]=-1.0; 349674ae819SStefano Zampini scaling_factors[partial_sum+s]=array[s]; 350674ae819SStefano Zampini } 351674ae819SStefano Zampini for (s=0;s<n_pos_values;s++) { 352674ae819SStefano Zampini l2g_indices [partial_sum+s+n_neg_values]=aux_sums[n_neg_values]+s+n_global_lambda; 353674ae819SStefano Zampini cols_B_delta [partial_sum+s+n_neg_values]=dual_dofs_boundary_indices[i]; 354674ae819SStefano Zampini vals_B_delta [partial_sum+s+n_neg_values]=1.0; 355674ae819SStefano Zampini scaling_factors[partial_sum+s+n_neg_values]=array[s+n_neg_values]; 356674ae819SStefano Zampini } 357674ae819SStefano Zampini partial_sum += j; 358674ae819SStefano Zampini } else { 359674ae819SStefano Zampini /* l2g_indices and default cols and vals of B_delta */ 360674ae819SStefano Zampini for (s=0;s<j;s++) { 361674ae819SStefano Zampini l2g_indices [partial_sum+s]=n_global_lambda+s; 362674ae819SStefano Zampini cols_B_delta [partial_sum+s]=dual_dofs_boundary_indices[i]; 363674ae819SStefano Zampini vals_B_delta [partial_sum+s]=0.0; 364674ae819SStefano Zampini } 365674ae819SStefano Zampini /* B_delta */ 366674ae819SStefano Zampini if ( n_neg_values > 0 ) { /* there's a rank next to me to the left */ 367674ae819SStefano Zampini vals_B_delta [partial_sum+n_neg_values-1]=-1.0; 368674ae819SStefano Zampini } 369674ae819SStefano Zampini if ( n_neg_values < j ) { /* there's a rank next to me to the right */ 370674ae819SStefano Zampini vals_B_delta [partial_sum+n_neg_values]=1.0; 371674ae819SStefano Zampini } 372674ae819SStefano Zampini /* scaling as in Klawonn-Widlund 1999*/ 373674ae819SStefano Zampini for (s=0;s<n_neg_values;s++) { 374674ae819SStefano Zampini scalar_value = 0.0; 375674ae819SStefano Zampini for (k=0;k<s+1;k++) { 376674ae819SStefano Zampini scalar_value += array[k]; 377674ae819SStefano Zampini } 378674ae819SStefano Zampini scaling_factors[partial_sum+s] = -scalar_value; 379674ae819SStefano Zampini } 380674ae819SStefano Zampini for (s=0;s<n_pos_values;s++) { 381674ae819SStefano Zampini scalar_value = 0.0; 382674ae819SStefano Zampini for (k=s+n_neg_values;k<j;k++) { 383674ae819SStefano Zampini scalar_value += array[k]; 384674ae819SStefano Zampini } 385674ae819SStefano Zampini scaling_factors[partial_sum+s+n_neg_values] = scalar_value; 386674ae819SStefano Zampini } 387674ae819SStefano Zampini partial_sum += j; 388674ae819SStefano Zampini } 389dc456d91SStefano Zampini cum += aux_local_numbering_2[i]; 390674ae819SStefano Zampini } 391dc456d91SStefano Zampini ierr = ISRestoreIndices(subset_n,&aux_global_numbering);CHKERRQ(ierr); 392dc456d91SStefano Zampini ierr = ISDestroy(&subset_mult);CHKERRQ(ierr); 393dc456d91SStefano Zampini ierr = ISDestroy(&subset_n);CHKERRQ(ierr); 394674ae819SStefano Zampini ierr = PetscFree(aux_sums);CHKERRQ(ierr); 395674ae819SStefano Zampini ierr = PetscFree(aux_local_numbering_1);CHKERRQ(ierr); 396674ae819SStefano Zampini ierr = PetscFree(dual_dofs_boundary_indices);CHKERRQ(ierr); 397674ae819SStefano Zampini ierr = PetscFree(all_factors[0]);CHKERRQ(ierr); 398674ae819SStefano Zampini ierr = PetscFree(all_factors);CHKERRQ(ierr); 399674ae819SStefano Zampini 400674ae819SStefano Zampini /* Create local part of B_delta */ 401302440fdSBarry Smith ierr = MatCreate(PETSC_COMM_SELF,&fetidpmat_ctx->B_delta);CHKERRQ(ierr); 402674ae819SStefano Zampini ierr = MatSetSizes(fetidpmat_ctx->B_delta,n_local_lambda,pcis->n_B,n_local_lambda,pcis->n_B);CHKERRQ(ierr); 403674ae819SStefano Zampini ierr = MatSetType(fetidpmat_ctx->B_delta,MATSEQAIJ);CHKERRQ(ierr); 404674ae819SStefano Zampini ierr = MatSeqAIJSetPreallocation(fetidpmat_ctx->B_delta,1,NULL);CHKERRQ(ierr); 405674ae819SStefano Zampini ierr = MatSetOption(fetidpmat_ctx->B_delta,MAT_IGNORE_ZERO_ENTRIES,PETSC_TRUE);CHKERRQ(ierr); 406674ae819SStefano Zampini for (i=0;i<n_local_lambda;i++) { 407674ae819SStefano Zampini ierr = MatSetValue(fetidpmat_ctx->B_delta,i,cols_B_delta[i],vals_B_delta[i],INSERT_VALUES);CHKERRQ(ierr); 408674ae819SStefano Zampini } 409674ae819SStefano Zampini ierr = PetscFree(vals_B_delta);CHKERRQ(ierr); 410674ae819SStefano Zampini ierr = MatAssemblyBegin(fetidpmat_ctx->B_delta,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 411674ae819SStefano Zampini ierr = MatAssemblyEnd (fetidpmat_ctx->B_delta,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 412674ae819SStefano Zampini 413674ae819SStefano Zampini if (fully_redundant) { 414302440fdSBarry Smith ierr = MatCreate(PETSC_COMM_SELF,&ScalingMat);CHKERRQ(ierr); 415674ae819SStefano Zampini ierr = MatSetSizes(ScalingMat,n_local_lambda,n_local_lambda,n_local_lambda,n_local_lambda);CHKERRQ(ierr); 416674ae819SStefano Zampini ierr = MatSetType(ScalingMat,MATSEQAIJ);CHKERRQ(ierr); 417674ae819SStefano Zampini ierr = MatSeqAIJSetPreallocation(ScalingMat,1,NULL);CHKERRQ(ierr); 418674ae819SStefano Zampini for (i=0;i<n_local_lambda;i++) { 419674ae819SStefano Zampini ierr = MatSetValue(ScalingMat,i,i,scaling_factors[i],INSERT_VALUES);CHKERRQ(ierr); 420674ae819SStefano Zampini } 421674ae819SStefano Zampini ierr = MatAssemblyBegin(ScalingMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 422674ae819SStefano Zampini ierr = MatAssemblyEnd (ScalingMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 423674ae819SStefano Zampini ierr = MatMatMult(ScalingMat,fetidpmat_ctx->B_delta,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&fetidpmat_ctx->B_Ddelta);CHKERRQ(ierr); 424674ae819SStefano Zampini ierr = MatDestroy(&ScalingMat);CHKERRQ(ierr); 425674ae819SStefano Zampini } else { 426302440fdSBarry Smith ierr = MatCreate(PETSC_COMM_SELF,&fetidpmat_ctx->B_Ddelta);CHKERRQ(ierr); 427674ae819SStefano Zampini ierr = MatSetSizes(fetidpmat_ctx->B_Ddelta,n_local_lambda,pcis->n_B,n_local_lambda,pcis->n_B);CHKERRQ(ierr); 428674ae819SStefano Zampini ierr = MatSetType(fetidpmat_ctx->B_Ddelta,MATSEQAIJ);CHKERRQ(ierr); 429674ae819SStefano Zampini ierr = MatSeqAIJSetPreallocation(fetidpmat_ctx->B_Ddelta,1,NULL);CHKERRQ(ierr); 430674ae819SStefano Zampini for (i=0;i<n_local_lambda;i++) { 431674ae819SStefano Zampini ierr = MatSetValue(fetidpmat_ctx->B_Ddelta,i,cols_B_delta[i],scaling_factors[i],INSERT_VALUES);CHKERRQ(ierr); 432674ae819SStefano Zampini } 433674ae819SStefano Zampini ierr = MatAssemblyBegin(fetidpmat_ctx->B_Ddelta,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 434674ae819SStefano Zampini ierr = MatAssemblyEnd (fetidpmat_ctx->B_Ddelta,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 435674ae819SStefano Zampini } 436674ae819SStefano Zampini ierr = PetscFree(scaling_factors);CHKERRQ(ierr); 437674ae819SStefano Zampini ierr = PetscFree(cols_B_delta);CHKERRQ(ierr); 438674ae819SStefano Zampini 439*457d4a33Sstefano_zampini /* Layout of multipliers */ 440*457d4a33Sstefano_zampini ierr = PetscLayoutCreate(comm,&llay);CHKERRQ(ierr); 441*457d4a33Sstefano_zampini ierr = PetscLayoutSetBlockSize(llay,1);CHKERRQ(ierr); 442*457d4a33Sstefano_zampini ierr = PetscLayoutSetSize(llay,fetidpmat_ctx->n_lambda);CHKERRQ(ierr); 443*457d4a33Sstefano_zampini ierr = PetscLayoutSetUp(llay);CHKERRQ(ierr); 444*457d4a33Sstefano_zampini ierr = PetscLayoutGetLocalSize(llay,&fetidpmat_ctx->n);CHKERRQ(ierr); 445*457d4a33Sstefano_zampini 446*457d4a33Sstefano_zampini /* fetidpmat sizes */ 447*457d4a33Sstefano_zampini fetidpmat_ctx->n += nPgl; 448*457d4a33Sstefano_zampini fetidpmat_ctx->N = fetidpmat_ctx->n_lambda+nPg; 449*457d4a33Sstefano_zampini 450*457d4a33Sstefano_zampini /* Local work vector of multipliers */ 451*457d4a33Sstefano_zampini ierr = VecCreate(PETSC_COMM_SELF,&fetidpmat_ctx->lambda_local);CHKERRQ(ierr); 452*457d4a33Sstefano_zampini ierr = VecSetSizes(fetidpmat_ctx->lambda_local,n_local_lambda,n_local_lambda);CHKERRQ(ierr); 453*457d4a33Sstefano_zampini ierr = VecSetType(fetidpmat_ctx->lambda_local,VECSEQ);CHKERRQ(ierr); 454*457d4a33Sstefano_zampini 455*457d4a33Sstefano_zampini /* Global vector for FETI-DP linear system */ 456*457d4a33Sstefano_zampini ierr = VecCreate(comm,&fetidp_global);CHKERRQ(ierr); 457*457d4a33Sstefano_zampini ierr = VecSetSizes(fetidp_global,fetidpmat_ctx->n,fetidpmat_ctx->N);CHKERRQ(ierr); 458*457d4a33Sstefano_zampini ierr = VecSetType(fetidp_global,VECMPI);CHKERRQ(ierr); 459*457d4a33Sstefano_zampini ierr = VecSetUp(fetidp_global);CHKERRQ(ierr); 460*457d4a33Sstefano_zampini 461*457d4a33Sstefano_zampini /* Decide layout fetidp dofs if it is a saddle point problem 462*457d4a33Sstefano_zampini pressures ordered first in the local part of the global vector 463*457d4a33Sstefano_zampini of the FETI-DP linear system */ 464*457d4a33Sstefano_zampini if (nPg) { 465*457d4a33Sstefano_zampini Mat nA; 466*457d4a33Sstefano_zampini Vec ov; 467*457d4a33Sstefano_zampini IS IS_l2g_p,ais,AmgP; 468*457d4a33Sstefano_zampini PetscLayout alay; 469*457d4a33Sstefano_zampini const PetscInt *idxs,*pranges,*aranges,*lranges; 470*457d4a33Sstefano_zampini PetscInt *l2g_indices_p,rst,ren; 471*457d4a33Sstefano_zampini 472*457d4a33Sstefano_zampini ierr = PetscMalloc1(nPl,&l2g_indices_p);CHKERRQ(ierr); 473*457d4a33Sstefano_zampini ierr = VecGetLayout(fetidp_global,&alay);CHKERRQ(ierr); 474*457d4a33Sstefano_zampini ierr = PetscLayoutGetRanges(alay,&aranges);CHKERRQ(ierr); 475*457d4a33Sstefano_zampini ierr = PetscLayoutGetRanges(play,&pranges);CHKERRQ(ierr); 476*457d4a33Sstefano_zampini ierr = PetscLayoutGetRanges(llay,&lranges);CHKERRQ(ierr); 477*457d4a33Sstefano_zampini ierr = ISLocalToGlobalMappingGetIndices(l2gmap_p,&idxs);CHKERRQ(ierr); 478*457d4a33Sstefano_zampini for (i=0;i<nPl;i++) { 479*457d4a33Sstefano_zampini PetscInt owner; 480*457d4a33Sstefano_zampini 481*457d4a33Sstefano_zampini ierr = PetscLayoutFindOwner(play,idxs[i],&owner);CHKERRQ(ierr); 482*457d4a33Sstefano_zampini l2g_indices_p[i] = idxs[i]-pranges[owner]+aranges[owner]; 483*457d4a33Sstefano_zampini } 484*457d4a33Sstefano_zampini ierr = ISLocalToGlobalMappingRestoreIndices(l2gmap_p,&idxs);CHKERRQ(ierr); 485*457d4a33Sstefano_zampini ierr = ISCreateGeneral(comm,nPl,l2g_indices_p,PETSC_OWN_POINTER,&IS_l2g_p);CHKERRQ(ierr); 486*457d4a33Sstefano_zampini /* local to global scatter for interface pressure */ 487*457d4a33Sstefano_zampini ierr = VecScatterCreate(fetidpmat_ctx->vP,NULL,fetidp_global,IS_l2g_p,&fetidpmat_ctx->l2g_p);CHKERRQ(ierr); 488*457d4a33Sstefano_zampini ierr = ISDestroy(&IS_l2g_p);CHKERRQ(ierr); 489*457d4a33Sstefano_zampini 490*457d4a33Sstefano_zampini /* shift local to global scatter for multipliers */ 491*457d4a33Sstefano_zampini for (i=0;i<n_local_lambda;i++) { 492*457d4a33Sstefano_zampini PetscInt owner,ps; 493*457d4a33Sstefano_zampini 494*457d4a33Sstefano_zampini ierr = PetscLayoutFindOwner(llay,l2g_indices[i],&owner);CHKERRQ(ierr); 495*457d4a33Sstefano_zampini ps = pranges[owner+1]-pranges[owner]; 496*457d4a33Sstefano_zampini l2g_indices[i] = l2g_indices[i]-lranges[owner]+aranges[owner]+ps; 497*457d4a33Sstefano_zampini } 498*457d4a33Sstefano_zampini 499*457d4a33Sstefano_zampini /* compute additional scatters needed */ 500*457d4a33Sstefano_zampini ierr = PCGetOperators(fetidpmat_ctx->pc,&nA,NULL);CHKERRQ(ierr); 501*457d4a33Sstefano_zampini ierr = MatCreateVecs(nA,&fetidpmat_ctx->work,NULL);CHKERRQ(ierr); 502*457d4a33Sstefano_zampini ierr = MatCreateVecs(oA,&ov,NULL);CHKERRQ(ierr); 503*457d4a33Sstefano_zampini /* scatter from alldofs to alldofs \setminus interfacepressures */ 504*457d4a33Sstefano_zampini ierr = MatGetOwnershipRange(oA,&rst,&ren);CHKERRQ(ierr); 505*457d4a33Sstefano_zampini ierr = ISCreateStride(comm,ren-rst,rst,1,&ais);CHKERRQ(ierr); 506*457d4a33Sstefano_zampini ierr = ISDifference(ais,pP,&AmgP);CHKERRQ(ierr); 507*457d4a33Sstefano_zampini ierr = ISDestroy(&ais);CHKERRQ(ierr); 508*457d4a33Sstefano_zampini ierr = VecScatterCreate(ov,AmgP,fetidpmat_ctx->work,NULL,&fetidpmat_ctx->g2g);CHKERRQ(ierr); 509*457d4a33Sstefano_zampini ierr = ISDestroy(&AmgP);CHKERRQ(ierr); 510*457d4a33Sstefano_zampini /* scatter from alldofs to local interface pressures */ 511*457d4a33Sstefano_zampini ierr = VecScatterCreate(ov,gP,fetidpmat_ctx->vP,NULL,&fetidpmat_ctx->g2l_p);CHKERRQ(ierr); 512*457d4a33Sstefano_zampini ierr = VecDestroy(&ov);CHKERRQ(ierr); 513*457d4a33Sstefano_zampini } 514*457d4a33Sstefano_zampini ierr = PetscLayoutDestroy(&llay);CHKERRQ(ierr); 515*457d4a33Sstefano_zampini ierr = PetscLayoutDestroy(&play);CHKERRQ(ierr); 516*457d4a33Sstefano_zampini ierr = ISCreateGeneral(comm,n_local_lambda,l2g_indices,PETSC_OWN_POINTER,&IS_l2g_lambda);CHKERRQ(ierr); 517*457d4a33Sstefano_zampini /* scatter for multipliers */ 518*457d4a33Sstefano_zampini ierr = VecScatterCreate(fetidpmat_ctx->lambda_local,NULL,fetidp_global,IS_l2g_lambda,&fetidpmat_ctx->l2g_lambda);CHKERRQ(ierr); 519*457d4a33Sstefano_zampini ierr = ISDestroy(&IS_l2g_lambda);CHKERRQ(ierr); 520*457d4a33Sstefano_zampini ierr = ISLocalToGlobalMappingDestroy(&l2gmap_p);CHKERRQ(ierr); 521*457d4a33Sstefano_zampini 522674ae819SStefano Zampini /* Create some vectors needed by fetidp */ 523674ae819SStefano Zampini ierr = VecDuplicate(pcis->vec1_B,&fetidpmat_ctx->temp_solution_B);CHKERRQ(ierr); 524674ae819SStefano Zampini ierr = VecDuplicate(pcis->vec1_D,&fetidpmat_ctx->temp_solution_D);CHKERRQ(ierr); 525674ae819SStefano Zampini 526674ae819SStefano Zampini test_fetidp = PETSC_FALSE; 527c5929fdfSBarry Smith ierr = PetscOptionsGetBool(NULL,NULL,"-fetidp_check",&test_fetidp,NULL);CHKERRQ(ierr); 528674ae819SStefano Zampini 529674ae819SStefano Zampini if (test_fetidp && !pcbddc->use_deluxe_scaling) { 530674ae819SStefano Zampini 5315b08dc53SStefano Zampini PetscReal real_value; 5325b08dc53SStefano Zampini 533674ae819SStefano Zampini ierr = PetscViewerASCIIGetStdout(comm,&viewer);CHKERRQ(ierr); 5341575c14dSBarry Smith ierr = PetscViewerASCIIPushSynchronized(viewer);CHKERRQ(ierr); 535*457d4a33Sstefano_zampini ierr = PetscViewerASCIIPrintf(viewer,"----------FETI_DP SIZES--------------\n");CHKERRQ(ierr); 536*457d4a33Sstefano_zampini ierr = PetscViewerASCIISynchronizedPrintf(viewer,"[%04d]: Sizes %D %D %D %D %D\n",rank,fetidpmat_ctx->n,fetidpmat_ctx->N,nPgl,nPg,nPl);CHKERRQ(ierr); 537*457d4a33Sstefano_zampini ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 538674ae819SStefano Zampini ierr = PetscViewerASCIIPrintf(viewer,"----------FETI_DP TESTS--------------\n");CHKERRQ(ierr); 539674ae819SStefano Zampini ierr = PetscViewerASCIIPrintf(viewer,"All tests should return zero!\n");CHKERRQ(ierr); 540674ae819SStefano Zampini ierr = PetscViewerASCIIPrintf(viewer,"FETIDP MAT context in the ");CHKERRQ(ierr); 541674ae819SStefano Zampini if (fully_redundant) { 542674ae819SStefano Zampini ierr = PetscViewerASCIIPrintf(viewer,"fully redundant case for lagrange multipliers.\n");CHKERRQ(ierr); 543674ae819SStefano Zampini } else { 544674ae819SStefano Zampini ierr = PetscViewerASCIIPrintf(viewer,"Non-fully redundant case for lagrange multiplier.\n");CHKERRQ(ierr); 545674ae819SStefano Zampini } 546674ae819SStefano Zampini ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 547674ae819SStefano Zampini 548674ae819SStefano Zampini /******************************************************************/ 549674ae819SStefano Zampini /* TEST A/B: Test numbering of global lambda dofs */ 550674ae819SStefano Zampini /******************************************************************/ 551674ae819SStefano Zampini 552674ae819SStefano Zampini ierr = VecDuplicate(fetidpmat_ctx->lambda_local,&test_vec);CHKERRQ(ierr); 553*457d4a33Sstefano_zampini ierr = VecSet(fetidp_global,1.0);CHKERRQ(ierr); 554674ae819SStefano Zampini ierr = VecSet(test_vec,1.0);CHKERRQ(ierr); 555*457d4a33Sstefano_zampini ierr = VecScatterBegin(fetidpmat_ctx->l2g_lambda,fetidp_global,fetidpmat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 556*457d4a33Sstefano_zampini ierr = VecScatterEnd (fetidpmat_ctx->l2g_lambda,fetidp_global,fetidpmat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 557674ae819SStefano Zampini scalar_value = -1.0; 558674ae819SStefano Zampini ierr = VecAXPY(test_vec,scalar_value,fetidpmat_ctx->lambda_local);CHKERRQ(ierr); 5595b08dc53SStefano Zampini ierr = VecNorm(test_vec,NORM_INFINITY,&real_value);CHKERRQ(ierr); 560674ae819SStefano Zampini ierr = VecDestroy(&test_vec);CHKERRQ(ierr); 5615b08dc53SStefano Zampini ierr = PetscViewerASCIISynchronizedPrintf(viewer,"A[%04d]: CHECK glob to loc: % 1.14e\n",rank,real_value);CHKERRQ(ierr); 562674ae819SStefano Zampini ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 563674ae819SStefano Zampini if (fully_redundant) { 564*457d4a33Sstefano_zampini ierr = VecSet(fetidp_global,0.0);CHKERRQ(ierr); 565674ae819SStefano Zampini ierr = VecSet(fetidpmat_ctx->lambda_local,0.5);CHKERRQ(ierr); 566*457d4a33Sstefano_zampini ierr = VecScatterBegin(fetidpmat_ctx->l2g_lambda,fetidpmat_ctx->lambda_local,fetidp_global,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 567*457d4a33Sstefano_zampini ierr = VecScatterEnd (fetidpmat_ctx->l2g_lambda,fetidpmat_ctx->lambda_local,fetidp_global,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 568*457d4a33Sstefano_zampini ierr = VecSum(fetidp_global,&scalar_value);CHKERRQ(ierr); 5695b08dc53SStefano Zampini ierr = PetscViewerASCIISynchronizedPrintf(viewer,"B[%04d]: CHECK loc to glob: % 1.14e\n",rank,PetscRealPart(scalar_value)-fetidpmat_ctx->n_lambda);CHKERRQ(ierr); 570674ae819SStefano Zampini ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 571674ae819SStefano Zampini } 572674ae819SStefano Zampini 573674ae819SStefano Zampini /******************************************************************/ 574674ae819SStefano Zampini /* TEST C: It should holds B_delta*w=0, w\in\widehat{W} */ 575674ae819SStefano Zampini /* This is the meaning of the B matrix */ 576674ae819SStefano Zampini /******************************************************************/ 577674ae819SStefano Zampini 578674ae819SStefano Zampini ierr = VecSetRandom(pcis->vec1_N,NULL);CHKERRQ(ierr); 579674ae819SStefano Zampini ierr = VecSet(pcis->vec1_global,0.0);CHKERRQ(ierr); 580e176bc59SStefano Zampini ierr = VecScatterBegin(matis->rctx,pcis->vec1_N,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 581e176bc59SStefano Zampini ierr = VecScatterEnd(matis->rctx,pcis->vec1_N,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 582e176bc59SStefano Zampini ierr = VecScatterBegin(matis->rctx,pcis->vec1_global,pcis->vec1_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 583e176bc59SStefano Zampini ierr = VecScatterEnd(matis->rctx,pcis->vec1_global,pcis->vec1_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 584674ae819SStefano Zampini ierr = VecScatterBegin(pcis->N_to_B,pcis->vec1_N,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 585674ae819SStefano Zampini ierr = VecScatterEnd (pcis->N_to_B,pcis->vec1_N,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 586674ae819SStefano Zampini /* Action of B_delta */ 587674ae819SStefano Zampini ierr = MatMult(fetidpmat_ctx->B_delta,pcis->vec1_B,fetidpmat_ctx->lambda_local);CHKERRQ(ierr); 588*457d4a33Sstefano_zampini ierr = VecSet(fetidp_global,0.0);CHKERRQ(ierr); 589*457d4a33Sstefano_zampini ierr = VecScatterBegin(fetidpmat_ctx->l2g_lambda,fetidpmat_ctx->lambda_local,fetidp_global,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 590*457d4a33Sstefano_zampini ierr = VecScatterEnd (fetidpmat_ctx->l2g_lambda,fetidpmat_ctx->lambda_local,fetidp_global,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 591*457d4a33Sstefano_zampini ierr = VecNorm(fetidp_global,NORM_INFINITY,&real_value);CHKERRQ(ierr); 5925b08dc53SStefano Zampini ierr = PetscViewerASCIIPrintf(viewer,"C[coll]: CHECK infty norm of B_delta*w (w continuous): % 1.14e\n",real_value);CHKERRQ(ierr); 593674ae819SStefano Zampini ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 594674ae819SStefano Zampini 595674ae819SStefano Zampini /******************************************************************/ 596674ae819SStefano Zampini /* TEST D: It should holds E_Dw = w - P_Dw w\in\widetilde{W} */ 597674ae819SStefano Zampini /* E_D = R_D^TR */ 598674ae819SStefano Zampini /* P_D = B_{D,delta}^T B_{delta} */ 599674ae819SStefano Zampini /* eq.44 Mandel Tezaur and Dohrmann 2005 */ 600674ae819SStefano Zampini /******************************************************************/ 601674ae819SStefano Zampini 602674ae819SStefano Zampini /* compute a random vector in \widetilde{W} */ 603674ae819SStefano Zampini ierr = VecSetRandom(pcis->vec1_N,NULL);CHKERRQ(ierr); 604674ae819SStefano Zampini scalar_value = 0.0; /* set zero at vertices */ 605674ae819SStefano Zampini ierr = VecGetArray(pcis->vec1_N,&array);CHKERRQ(ierr); 606674ae819SStefano Zampini for (i=0;i<n_vertices;i++) { array[vertex_indices[i]]=scalar_value; } 607674ae819SStefano Zampini ierr = VecRestoreArray(pcis->vec1_N,&array);CHKERRQ(ierr); 608674ae819SStefano Zampini /* store w for final comparison */ 609674ae819SStefano Zampini ierr = VecDuplicate(pcis->vec1_B,&test_vec);CHKERRQ(ierr); 610674ae819SStefano Zampini ierr = VecScatterBegin(pcis->N_to_B,pcis->vec1_N,test_vec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 611674ae819SStefano Zampini ierr = VecScatterEnd (pcis->N_to_B,pcis->vec1_N,test_vec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 612674ae819SStefano Zampini 613674ae819SStefano Zampini /* Jump operator P_D : results stored in pcis->vec1_B */ 614674ae819SStefano Zampini 615674ae819SStefano Zampini ierr = VecScatterBegin(pcis->N_to_B,pcis->vec1_N,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 616674ae819SStefano Zampini ierr = VecScatterEnd (pcis->N_to_B,pcis->vec1_N,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 617674ae819SStefano Zampini /* Action of B_delta */ 618674ae819SStefano Zampini ierr = MatMult(fetidpmat_ctx->B_delta,pcis->vec1_B,fetidpmat_ctx->lambda_local);CHKERRQ(ierr); 619*457d4a33Sstefano_zampini ierr = VecSet(fetidp_global,0.0);CHKERRQ(ierr); 620*457d4a33Sstefano_zampini ierr = VecScatterBegin(fetidpmat_ctx->l2g_lambda,fetidpmat_ctx->lambda_local,fetidp_global,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 621*457d4a33Sstefano_zampini ierr = VecScatterEnd (fetidpmat_ctx->l2g_lambda,fetidpmat_ctx->lambda_local,fetidp_global,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 622674ae819SStefano Zampini /* Action of B_Ddelta^T */ 623*457d4a33Sstefano_zampini ierr = VecScatterBegin(fetidpmat_ctx->l2g_lambda,fetidp_global,fetidpmat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 624*457d4a33Sstefano_zampini ierr = VecScatterEnd (fetidpmat_ctx->l2g_lambda,fetidp_global,fetidpmat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 625674ae819SStefano Zampini ierr = MatMultTranspose(fetidpmat_ctx->B_Ddelta,fetidpmat_ctx->lambda_local,pcis->vec1_B);CHKERRQ(ierr); 626674ae819SStefano Zampini 627674ae819SStefano Zampini /* Average operator E_D : results stored in pcis->vec2_B */ 628674ae819SStefano Zampini ierr = VecScatterBegin(pcis->N_to_B,pcis->vec1_N,pcis->vec2_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 629674ae819SStefano Zampini ierr = VecScatterEnd (pcis->N_to_B,pcis->vec1_N,pcis->vec2_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 630674ae819SStefano Zampini ierr = PCBDDCScalingExtension(fetidpmat_ctx->pc,pcis->vec2_B,pcis->vec1_global);CHKERRQ(ierr); 631674ae819SStefano Zampini ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_global,pcis->vec2_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 632674ae819SStefano Zampini ierr = VecScatterEnd (pcis->global_to_B,pcis->vec1_global,pcis->vec2_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 633674ae819SStefano Zampini 634674ae819SStefano Zampini /* test E_D=I-P_D */ 635674ae819SStefano Zampini scalar_value = 1.0; 636674ae819SStefano Zampini ierr = VecAXPY(pcis->vec1_B,scalar_value,pcis->vec2_B);CHKERRQ(ierr); 637674ae819SStefano Zampini scalar_value = -1.0; 638674ae819SStefano Zampini ierr = VecAXPY(pcis->vec1_B,scalar_value,test_vec);CHKERRQ(ierr); 6395b08dc53SStefano Zampini ierr = VecNorm(pcis->vec1_B,NORM_INFINITY,&real_value);CHKERRQ(ierr); 640674ae819SStefano Zampini ierr = VecDestroy(&test_vec);CHKERRQ(ierr); 6415b08dc53SStefano Zampini ierr = PetscViewerASCIISynchronizedPrintf(viewer,"D[%04d] CHECK infty norm of E_D + P_D - I: % 1.14e\n",rank,real_value);CHKERRQ(ierr); 642674ae819SStefano Zampini ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 643674ae819SStefano Zampini 644674ae819SStefano Zampini /******************************************************************/ 645674ae819SStefano Zampini /* TEST E: It should holds R_D^TP_Dw=0 w\in\widetilde{W} */ 646674ae819SStefano Zampini /* eq.48 Mandel Tezaur and Dohrmann 2005 */ 647674ae819SStefano Zampini /******************************************************************/ 648674ae819SStefano Zampini 649674ae819SStefano Zampini ierr = VecSetRandom(pcis->vec1_N,NULL);CHKERRQ(ierr); 650674ae819SStefano Zampini ierr = VecGetArray(pcis->vec1_N,&array);CHKERRQ(ierr); 651674ae819SStefano Zampini scalar_value = 0.0; /* set zero at vertices */ 652674ae819SStefano Zampini for (i=0;i<n_vertices;i++) { array[vertex_indices[i]]=scalar_value; } 653674ae819SStefano Zampini ierr = VecRestoreArray(pcis->vec1_N,&array);CHKERRQ(ierr); 654674ae819SStefano Zampini 655674ae819SStefano Zampini /* Jump operator P_D : results stored in pcis->vec1_B */ 656674ae819SStefano Zampini 657674ae819SStefano Zampini ierr = VecScatterBegin(pcis->N_to_B,pcis->vec1_N,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 658674ae819SStefano Zampini ierr = VecScatterEnd (pcis->N_to_B,pcis->vec1_N,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 659674ae819SStefano Zampini /* Action of B_delta */ 660674ae819SStefano Zampini ierr = MatMult(fetidpmat_ctx->B_delta,pcis->vec1_B,fetidpmat_ctx->lambda_local);CHKERRQ(ierr); 661*457d4a33Sstefano_zampini ierr = VecSet(fetidp_global,0.0);CHKERRQ(ierr); 662*457d4a33Sstefano_zampini ierr = VecScatterBegin(fetidpmat_ctx->l2g_lambda,fetidpmat_ctx->lambda_local,fetidp_global,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 663*457d4a33Sstefano_zampini ierr = VecScatterEnd (fetidpmat_ctx->l2g_lambda,fetidpmat_ctx->lambda_local,fetidp_global,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 664674ae819SStefano Zampini /* Action of B_Ddelta^T */ 665*457d4a33Sstefano_zampini ierr = VecScatterBegin(fetidpmat_ctx->l2g_lambda,fetidp_global,fetidpmat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 666*457d4a33Sstefano_zampini ierr = VecScatterEnd (fetidpmat_ctx->l2g_lambda,fetidp_global,fetidpmat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 667674ae819SStefano Zampini ierr = MatMultTranspose(fetidpmat_ctx->B_Ddelta,fetidpmat_ctx->lambda_local,pcis->vec1_B);CHKERRQ(ierr); 668674ae819SStefano Zampini /* scaling */ 669674ae819SStefano Zampini ierr = PCBDDCScalingExtension(fetidpmat_ctx->pc,pcis->vec1_B,pcis->vec1_global);CHKERRQ(ierr); 6705b08dc53SStefano Zampini ierr = VecNorm(pcis->vec1_global,NORM_INFINITY,&real_value);CHKERRQ(ierr); 6715b08dc53SStefano Zampini ierr = PetscViewerASCIIPrintf(viewer,"E[coll]: CHECK infty norm of R^T_D P_D: % 1.14e\n",real_value);CHKERRQ(ierr); 672674ae819SStefano Zampini ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 673674ae819SStefano Zampini 674674ae819SStefano Zampini if (!fully_redundant) { 675674ae819SStefano Zampini /******************************************************************/ 676674ae819SStefano Zampini /* TEST F: It should holds B_{delta}B^T_{D,delta}=I */ 677674ae819SStefano Zampini /* Corollary thm 14 Mandel Tezaur and Dohrmann 2005 */ 678674ae819SStefano Zampini /******************************************************************/ 679*457d4a33Sstefano_zampini ierr = VecDuplicate(fetidp_global,&test_vec);CHKERRQ(ierr); 680*457d4a33Sstefano_zampini ierr = VecSetRandom(fetidp_global,NULL);CHKERRQ(ierr); 681*457d4a33Sstefano_zampini if (fetidpmat_ctx->l2g_p) { 682*457d4a33Sstefano_zampini ierr = VecSet(fetidpmat_ctx->vP,0.);CHKERRQ(ierr); 683*457d4a33Sstefano_zampini ierr = VecScatterBegin(fetidpmat_ctx->l2g_p,fetidpmat_ctx->vP,fetidp_global,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 684*457d4a33Sstefano_zampini ierr = VecScatterEnd(fetidpmat_ctx->l2g_p,fetidpmat_ctx->vP,fetidp_global,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 685*457d4a33Sstefano_zampini } 686674ae819SStefano Zampini /* Action of B_Ddelta^T */ 687*457d4a33Sstefano_zampini ierr = VecScatterBegin(fetidpmat_ctx->l2g_lambda,fetidp_global,fetidpmat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 688*457d4a33Sstefano_zampini ierr = VecScatterEnd (fetidpmat_ctx->l2g_lambda,fetidp_global,fetidpmat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 689674ae819SStefano Zampini ierr = MatMultTranspose(fetidpmat_ctx->B_Ddelta,fetidpmat_ctx->lambda_local,pcis->vec1_B);CHKERRQ(ierr); 690674ae819SStefano Zampini /* Action of B_delta */ 691674ae819SStefano Zampini ierr = MatMult(fetidpmat_ctx->B_delta,pcis->vec1_B,fetidpmat_ctx->lambda_local);CHKERRQ(ierr); 692674ae819SStefano Zampini ierr = VecSet(test_vec,0.0);CHKERRQ(ierr); 693674ae819SStefano Zampini ierr = VecScatterBegin(fetidpmat_ctx->l2g_lambda,fetidpmat_ctx->lambda_local,test_vec,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 694674ae819SStefano Zampini ierr = VecScatterEnd (fetidpmat_ctx->l2g_lambda,fetidpmat_ctx->lambda_local,test_vec,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 695674ae819SStefano Zampini scalar_value = -1.0; 696*457d4a33Sstefano_zampini ierr = VecAXPY(fetidp_global,scalar_value,test_vec);CHKERRQ(ierr); 697*457d4a33Sstefano_zampini ierr = VecNorm(fetidp_global,NORM_INFINITY,&real_value);CHKERRQ(ierr); 6985b08dc53SStefano Zampini ierr = PetscViewerASCIIPrintf(viewer,"E[coll]: CHECK infty norm of P^T_D - I: % 1.14e\n",real_value);CHKERRQ(ierr); 699674ae819SStefano Zampini ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 700674ae819SStefano Zampini ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 701674ae819SStefano Zampini ierr = VecDestroy(&test_vec);CHKERRQ(ierr); 702674ae819SStefano Zampini } 703674ae819SStefano Zampini } 704674ae819SStefano Zampini /* final cleanup */ 705*457d4a33Sstefano_zampini ierr = VecDestroy(&fetidp_global);CHKERRQ(ierr); 706674ae819SStefano Zampini PetscFunctionReturn(0); 707674ae819SStefano Zampini } 708674ae819SStefano Zampini 709674ae819SStefano Zampini #undef __FUNCT__ 710674ae819SStefano Zampini #define __FUNCT__ "PCBDDCSetupFETIDPPCContext" 711674ae819SStefano Zampini PetscErrorCode PCBDDCSetupFETIDPPCContext(Mat fetimat, FETIDPPC_ctx fetidppc_ctx) 712674ae819SStefano Zampini { 713674ae819SStefano Zampini FETIDPMat_ctx mat_ctx; 714ed6c3d69SStefano Zampini PC_IS *pcis; 715f28b6018SStefano Zampini PetscBool lumped = PETSC_FALSE; 716674ae819SStefano Zampini PetscErrorCode ierr; 717674ae819SStefano Zampini 718674ae819SStefano Zampini PetscFunctionBegin; 719674ae819SStefano Zampini ierr = MatShellGetContext(fetimat,(void**)&mat_ctx);CHKERRQ(ierr); 720674ae819SStefano Zampini /* get references from objects created when setting up feti mat context */ 721674ae819SStefano Zampini ierr = PetscObjectReference((PetscObject)mat_ctx->lambda_local);CHKERRQ(ierr); 722674ae819SStefano Zampini fetidppc_ctx->lambda_local = mat_ctx->lambda_local; 723674ae819SStefano Zampini ierr = PetscObjectReference((PetscObject)mat_ctx->B_Ddelta);CHKERRQ(ierr); 724674ae819SStefano Zampini fetidppc_ctx->B_Ddelta = mat_ctx->B_Ddelta; 725674ae819SStefano Zampini ierr = PetscObjectReference((PetscObject)mat_ctx->l2g_lambda);CHKERRQ(ierr); 726674ae819SStefano Zampini fetidppc_ctx->l2g_lambda = mat_ctx->l2g_lambda; 727f28b6018SStefano Zampini /* create preconditioner */ 728ed6c3d69SStefano Zampini pcis = (PC_IS*)fetidppc_ctx->pc->data; 729f28b6018SStefano Zampini /* Dirichlet preconditioner */ 730f28b6018SStefano Zampini ierr = PetscOptionsGetBool(NULL,NULL,"-fetidp_lumped",&lumped,NULL);CHKERRQ(ierr); 731f28b6018SStefano Zampini if (!lumped) { 732ed6c3d69SStefano Zampini ierr = MatCreateSchurComplement(pcis->A_II,pcis->A_II,pcis->A_IB,pcis->A_BI,pcis->A_BB,&fetidppc_ctx->S_j);CHKERRQ(ierr); 733ed6c3d69SStefano Zampini ierr = MatSchurComplementSetKSP(fetidppc_ctx->S_j,pcis->ksp_D);CHKERRQ(ierr); 734f28b6018SStefano Zampini } else { 735f28b6018SStefano Zampini ierr = PetscObjectReference((PetscObject)pcis->A_BB);CHKERRQ(ierr); 736f28b6018SStefano Zampini fetidppc_ctx->S_j = pcis->A_BB; 737f28b6018SStefano Zampini } 738674ae819SStefano Zampini PetscFunctionReturn(0); 739674ae819SStefano Zampini } 740674ae819SStefano Zampini 741674ae819SStefano Zampini #undef __FUNCT__ 742674ae819SStefano Zampini #define __FUNCT__ "FETIDPMatMult" 743674ae819SStefano Zampini PetscErrorCode FETIDPMatMult(Mat fetimat, Vec x, Vec y) 744674ae819SStefano Zampini { 745674ae819SStefano Zampini FETIDPMat_ctx mat_ctx; 746617d11aeSStefano Zampini PC_BDDC *pcbddc; 747674ae819SStefano Zampini PC_IS *pcis; 748674ae819SStefano Zampini PetscErrorCode ierr; 749674ae819SStefano Zampini 750674ae819SStefano Zampini PetscFunctionBegin; 751674ae819SStefano Zampini ierr = MatShellGetContext(fetimat,(void**)&mat_ctx);CHKERRQ(ierr); 752674ae819SStefano Zampini pcis = (PC_IS*)mat_ctx->pc->data; 753617d11aeSStefano Zampini pcbddc = (PC_BDDC*)mat_ctx->pc->data; 754674ae819SStefano Zampini /* Application of B_delta^T */ 755674ae819SStefano Zampini ierr = VecScatterBegin(mat_ctx->l2g_lambda,x,mat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 756674ae819SStefano Zampini ierr = VecScatterEnd(mat_ctx->l2g_lambda,x,mat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 757674ae819SStefano Zampini ierr = MatMultTranspose(mat_ctx->B_delta,mat_ctx->lambda_local,pcis->vec1_B);CHKERRQ(ierr); 758674ae819SStefano Zampini /* Application of \widetilde{S}^-1 */ 759674ae819SStefano Zampini ierr = VecSet(pcis->vec1_D,0.0);CHKERRQ(ierr); 760617d11aeSStefano Zampini ierr = PetscMemzero(pcbddc->benign_p0,pcbddc->benign_n*sizeof(PetscScalar));CHKERRQ(ierr); 761dc359a40SStefano Zampini ierr = PCBDDCApplyInterfacePreconditioner(mat_ctx->pc,PETSC_FALSE);CHKERRQ(ierr); 762c7ffc8ceSStefano Zampini ierr = PetscMemzero(pcbddc->benign_p0,pcbddc->benign_n*sizeof(PetscScalar));CHKERRQ(ierr); 763674ae819SStefano Zampini /* Application of B_delta */ 764674ae819SStefano Zampini ierr = MatMult(mat_ctx->B_delta,pcis->vec1_B,mat_ctx->lambda_local);CHKERRQ(ierr); 765674ae819SStefano Zampini ierr = VecSet(y,0.0);CHKERRQ(ierr); 766674ae819SStefano Zampini ierr = VecScatterBegin(mat_ctx->l2g_lambda,mat_ctx->lambda_local,y,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 767674ae819SStefano Zampini ierr = VecScatterEnd(mat_ctx->l2g_lambda,mat_ctx->lambda_local,y,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 768674ae819SStefano Zampini PetscFunctionReturn(0); 769674ae819SStefano Zampini } 770674ae819SStefano Zampini 771674ae819SStefano Zampini #undef __FUNCT__ 772edf7251bSStefano Zampini #define __FUNCT__ "FETIDPMatMultTranspose" 773edf7251bSStefano Zampini PetscErrorCode FETIDPMatMultTranspose(Mat fetimat, Vec x, Vec y) 774edf7251bSStefano Zampini { 775edf7251bSStefano Zampini FETIDPMat_ctx mat_ctx; 776edf7251bSStefano Zampini PC_IS *pcis; 777edf7251bSStefano Zampini PetscErrorCode ierr; 778edf7251bSStefano Zampini 779edf7251bSStefano Zampini PetscFunctionBegin; 780edf7251bSStefano Zampini ierr = MatShellGetContext(fetimat,(void**)&mat_ctx);CHKERRQ(ierr); 781edf7251bSStefano Zampini pcis = (PC_IS*)mat_ctx->pc->data; 782edf7251bSStefano Zampini /* Application of B_delta^T */ 783edf7251bSStefano Zampini ierr = VecScatterBegin(mat_ctx->l2g_lambda,x,mat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 784edf7251bSStefano Zampini ierr = VecScatterEnd(mat_ctx->l2g_lambda,x,mat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 785edf7251bSStefano Zampini ierr = MatMultTranspose(mat_ctx->B_delta,mat_ctx->lambda_local,pcis->vec1_B);CHKERRQ(ierr); 786edf7251bSStefano Zampini /* Application of \widetilde{S}^-1 */ 787edf7251bSStefano Zampini ierr = VecSet(pcis->vec1_D,0.0);CHKERRQ(ierr); 788edf7251bSStefano Zampini ierr = PCBDDCApplyInterfacePreconditioner(mat_ctx->pc,PETSC_TRUE);CHKERRQ(ierr); 789edf7251bSStefano Zampini /* Application of B_delta */ 790edf7251bSStefano Zampini ierr = MatMult(mat_ctx->B_delta,pcis->vec1_B,mat_ctx->lambda_local);CHKERRQ(ierr); 791edf7251bSStefano Zampini ierr = VecSet(y,0.0);CHKERRQ(ierr); 792edf7251bSStefano Zampini ierr = VecScatterBegin(mat_ctx->l2g_lambda,mat_ctx->lambda_local,y,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 793edf7251bSStefano Zampini ierr = VecScatterEnd(mat_ctx->l2g_lambda,mat_ctx->lambda_local,y,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 794edf7251bSStefano Zampini PetscFunctionReturn(0); 795edf7251bSStefano Zampini } 796edf7251bSStefano Zampini 797edf7251bSStefano Zampini #undef __FUNCT__ 798674ae819SStefano Zampini #define __FUNCT__ "FETIDPPCApply" 799674ae819SStefano Zampini PetscErrorCode FETIDPPCApply(PC fetipc, Vec x, Vec y) 800674ae819SStefano Zampini { 801674ae819SStefano Zampini FETIDPPC_ctx pc_ctx; 802674ae819SStefano Zampini PC_IS *pcis; 803674ae819SStefano Zampini PetscErrorCode ierr; 804674ae819SStefano Zampini 805674ae819SStefano Zampini PetscFunctionBegin; 806302440fdSBarry Smith ierr = PCShellGetContext(fetipc,(void**)&pc_ctx);CHKERRQ(ierr); 807674ae819SStefano Zampini pcis = (PC_IS*)pc_ctx->pc->data; 808674ae819SStefano Zampini /* Application of B_Ddelta^T */ 809674ae819SStefano Zampini ierr = VecScatterBegin(pc_ctx->l2g_lambda,x,pc_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 810674ae819SStefano Zampini ierr = VecScatterEnd(pc_ctx->l2g_lambda,x,pc_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 811674ae819SStefano Zampini ierr = VecSet(pcis->vec2_B,0.0);CHKERRQ(ierr); 812674ae819SStefano Zampini ierr = MatMultTranspose(pc_ctx->B_Ddelta,pc_ctx->lambda_local,pcis->vec2_B);CHKERRQ(ierr); 813ed6c3d69SStefano Zampini /* Application of local Schur complement */ 814ed6c3d69SStefano Zampini ierr = MatMult(pc_ctx->S_j,pcis->vec2_B,pcis->vec1_B);CHKERRQ(ierr); 815edf7251bSStefano Zampini /* Application of B_Ddelta */ 816edf7251bSStefano Zampini ierr = MatMult(pc_ctx->B_Ddelta,pcis->vec1_B,pc_ctx->lambda_local);CHKERRQ(ierr); 817edf7251bSStefano Zampini ierr = VecSet(y,0.0);CHKERRQ(ierr); 818edf7251bSStefano Zampini ierr = VecScatterBegin(pc_ctx->l2g_lambda,pc_ctx->lambda_local,y,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 819edf7251bSStefano Zampini ierr = VecScatterEnd(pc_ctx->l2g_lambda,pc_ctx->lambda_local,y,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 820edf7251bSStefano Zampini PetscFunctionReturn(0); 821edf7251bSStefano Zampini } 822edf7251bSStefano Zampini 823edf7251bSStefano Zampini #undef __FUNCT__ 824edf7251bSStefano Zampini #define __FUNCT__ "FETIDPPCApplyTranspose" 825edf7251bSStefano Zampini PetscErrorCode FETIDPPCApplyTranspose(PC fetipc, Vec x, Vec y) 826edf7251bSStefano Zampini { 827edf7251bSStefano Zampini FETIDPPC_ctx pc_ctx; 828edf7251bSStefano Zampini PC_IS *pcis; 829edf7251bSStefano Zampini PetscErrorCode ierr; 830edf7251bSStefano Zampini 831edf7251bSStefano Zampini PetscFunctionBegin; 832302440fdSBarry Smith ierr = PCShellGetContext(fetipc,(void**)&pc_ctx);CHKERRQ(ierr); 833edf7251bSStefano Zampini pcis = (PC_IS*)pc_ctx->pc->data; 834edf7251bSStefano Zampini /* Application of B_Ddelta^T */ 835edf7251bSStefano Zampini ierr = VecScatterBegin(pc_ctx->l2g_lambda,x,pc_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 836edf7251bSStefano Zampini ierr = VecScatterEnd(pc_ctx->l2g_lambda,x,pc_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 837edf7251bSStefano Zampini ierr = VecSet(pcis->vec2_B,0.0);CHKERRQ(ierr); 838edf7251bSStefano Zampini ierr = MatMultTranspose(pc_ctx->B_Ddelta,pc_ctx->lambda_local,pcis->vec2_B);CHKERRQ(ierr); 839ed6c3d69SStefano Zampini /* Application of local Schur complement */ 840ed6c3d69SStefano Zampini ierr = MatMultTranspose(pc_ctx->S_j,pcis->vec2_B,pcis->vec1_B);CHKERRQ(ierr); 841674ae819SStefano Zampini /* Application of B_Ddelta */ 842674ae819SStefano Zampini ierr = MatMult(pc_ctx->B_Ddelta,pcis->vec1_B,pc_ctx->lambda_local);CHKERRQ(ierr); 843674ae819SStefano Zampini ierr = VecSet(y,0.0);CHKERRQ(ierr); 844674ae819SStefano Zampini ierr = VecScatterBegin(pc_ctx->l2g_lambda,pc_ctx->lambda_local,y,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 845674ae819SStefano Zampini ierr = VecScatterEnd(pc_ctx->l2g_lambda,pc_ctx->lambda_local,y,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 846674ae819SStefano Zampini PetscFunctionReturn(0); 847674ae819SStefano Zampini } 848