1*5e5bbd0aSStefano Zampini #include <petsc/private/pcbddcimpl.h> 2*5e5bbd0aSStefano Zampini #include <petsc/private/pcbddcprivateimpl.h> 325d06dbeSstefano_zampini #include <petscblaslapack.h> 425d06dbeSstefano_zampini 525d06dbeSstefano_zampini static PetscErrorCode MatMult_BDdelta_deluxe_nonred(Mat A, Vec x, Vec y) 625d06dbeSstefano_zampini { 725d06dbeSstefano_zampini BDdelta_DN ctx; 825d06dbeSstefano_zampini 925d06dbeSstefano_zampini PetscFunctionBegin; 109566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(A,&ctx)); 119566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(ctx->BD,x,ctx->work)); 129566063dSJacob Faibussowitsch PetscCall(KSPSolveTranspose(ctx->kBD,ctx->work,y)); 13c0decd05SBarry Smith /* No PC so cannot propagate up failure in KSPSolveTranspose() */ 1425d06dbeSstefano_zampini PetscFunctionReturn(0); 1525d06dbeSstefano_zampini } 1625d06dbeSstefano_zampini 1725d06dbeSstefano_zampini static PetscErrorCode MatMultTranspose_BDdelta_deluxe_nonred(Mat A, Vec x, Vec y) 1825d06dbeSstefano_zampini { 1925d06dbeSstefano_zampini BDdelta_DN ctx; 2025d06dbeSstefano_zampini 2125d06dbeSstefano_zampini PetscFunctionBegin; 229566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(A,&ctx)); 239566063dSJacob Faibussowitsch PetscCall(KSPSolve(ctx->kBD,x,ctx->work)); 24c0decd05SBarry Smith /* No PC so cannot propagate up failure in KSPSolve() */ 259566063dSJacob Faibussowitsch PetscCall(MatMult(ctx->BD,ctx->work,y)); 2625d06dbeSstefano_zampini PetscFunctionReturn(0); 2725d06dbeSstefano_zampini } 2825d06dbeSstefano_zampini 2925d06dbeSstefano_zampini static PetscErrorCode MatDestroy_BDdelta_deluxe_nonred(Mat A) 3025d06dbeSstefano_zampini { 3125d06dbeSstefano_zampini BDdelta_DN ctx; 3225d06dbeSstefano_zampini 3325d06dbeSstefano_zampini PetscFunctionBegin; 349566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(A,&ctx)); 359566063dSJacob Faibussowitsch PetscCall(MatDestroy(&ctx->BD)); 369566063dSJacob Faibussowitsch PetscCall(KSPDestroy(&ctx->kBD)); 379566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ctx->work)); 389566063dSJacob Faibussowitsch PetscCall(PetscFree(ctx)); 3925d06dbeSstefano_zampini PetscFunctionReturn(0); 4025d06dbeSstefano_zampini } 4125d06dbeSstefano_zampini 42674ae819SStefano Zampini PetscErrorCode PCBDDCCreateFETIDPMatContext(PC pc, FETIDPMat_ctx *fetidpmat_ctx) 43674ae819SStefano Zampini { 44674ae819SStefano Zampini FETIDPMat_ctx newctx; 45674ae819SStefano Zampini 46674ae819SStefano Zampini PetscFunctionBegin; 479566063dSJacob Faibussowitsch PetscCall(PetscNew(&newctx)); 48674ae819SStefano Zampini /* increase the reference count for BDDC preconditioner */ 499566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)pc)); 50674ae819SStefano Zampini newctx->pc = pc; 51674ae819SStefano Zampini *fetidpmat_ctx = newctx; 52674ae819SStefano Zampini PetscFunctionReturn(0); 53674ae819SStefano Zampini } 54674ae819SStefano Zampini 55674ae819SStefano Zampini PetscErrorCode PCBDDCCreateFETIDPPCContext(PC pc, FETIDPPC_ctx *fetidppc_ctx) 56674ae819SStefano Zampini { 57674ae819SStefano Zampini FETIDPPC_ctx newctx; 58674ae819SStefano Zampini 59674ae819SStefano Zampini PetscFunctionBegin; 609566063dSJacob Faibussowitsch PetscCall(PetscNew(&newctx)); 61674ae819SStefano Zampini /* increase the reference count for BDDC preconditioner */ 629566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)pc)); 63674ae819SStefano Zampini newctx->pc = pc; 64674ae819SStefano Zampini *fetidppc_ctx = newctx; 65674ae819SStefano Zampini PetscFunctionReturn(0); 66674ae819SStefano Zampini } 67674ae819SStefano Zampini 68674ae819SStefano Zampini PetscErrorCode PCBDDCDestroyFETIDPMat(Mat A) 69674ae819SStefano Zampini { 70674ae819SStefano Zampini FETIDPMat_ctx mat_ctx; 71674ae819SStefano Zampini 72674ae819SStefano Zampini PetscFunctionBegin; 739566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(A,&mat_ctx)); 749566063dSJacob Faibussowitsch PetscCall(VecDestroy(&mat_ctx->lambda_local)); 759566063dSJacob Faibussowitsch PetscCall(VecDestroy(&mat_ctx->temp_solution_D)); 769566063dSJacob Faibussowitsch PetscCall(VecDestroy(&mat_ctx->temp_solution_B)); 779566063dSJacob Faibussowitsch PetscCall(MatDestroy(&mat_ctx->B_delta)); 789566063dSJacob Faibussowitsch PetscCall(MatDestroy(&mat_ctx->B_Ddelta)); 799566063dSJacob Faibussowitsch PetscCall(MatDestroy(&mat_ctx->B_BB)); 809566063dSJacob Faibussowitsch PetscCall(MatDestroy(&mat_ctx->B_BI)); 819566063dSJacob Faibussowitsch PetscCall(MatDestroy(&mat_ctx->Bt_BB)); 829566063dSJacob Faibussowitsch PetscCall(MatDestroy(&mat_ctx->Bt_BI)); 839566063dSJacob Faibussowitsch PetscCall(MatDestroy(&mat_ctx->C)); 849566063dSJacob Faibussowitsch PetscCall(VecDestroy(&mat_ctx->rhs_flip)); 859566063dSJacob Faibussowitsch PetscCall(VecDestroy(&mat_ctx->vP)); 869566063dSJacob Faibussowitsch PetscCall(VecDestroy(&mat_ctx->xPg)); 879566063dSJacob Faibussowitsch PetscCall(VecDestroy(&mat_ctx->yPg)); 889566063dSJacob Faibussowitsch PetscCall(VecScatterDestroy(&mat_ctx->l2g_lambda)); 899566063dSJacob Faibussowitsch PetscCall(VecScatterDestroy(&mat_ctx->l2g_lambda_only)); 909566063dSJacob Faibussowitsch PetscCall(VecScatterDestroy(&mat_ctx->l2g_p)); 919566063dSJacob Faibussowitsch PetscCall(VecScatterDestroy(&mat_ctx->g2g_p)); 929566063dSJacob Faibussowitsch PetscCall(PCDestroy(&mat_ctx->pc)); /* decrease PCBDDC reference count */ 939566063dSJacob Faibussowitsch PetscCall(ISDestroy(&mat_ctx->pressure)); 949566063dSJacob Faibussowitsch PetscCall(ISDestroy(&mat_ctx->lagrange)); 959566063dSJacob Faibussowitsch PetscCall(PetscFree(mat_ctx)); 96674ae819SStefano Zampini PetscFunctionReturn(0); 97674ae819SStefano Zampini } 98674ae819SStefano Zampini 99674ae819SStefano Zampini PetscErrorCode PCBDDCDestroyFETIDPPC(PC pc) 100674ae819SStefano Zampini { 101674ae819SStefano Zampini FETIDPPC_ctx pc_ctx; 102674ae819SStefano Zampini 103674ae819SStefano Zampini PetscFunctionBegin; 1049566063dSJacob Faibussowitsch PetscCall(PCShellGetContext(pc,&pc_ctx)); 1059566063dSJacob Faibussowitsch PetscCall(VecDestroy(&pc_ctx->lambda_local)); 1069566063dSJacob Faibussowitsch PetscCall(MatDestroy(&pc_ctx->B_Ddelta)); 1079566063dSJacob Faibussowitsch PetscCall(VecScatterDestroy(&pc_ctx->l2g_lambda)); 1089566063dSJacob Faibussowitsch PetscCall(MatDestroy(&pc_ctx->S_j)); 1099566063dSJacob Faibussowitsch PetscCall(PCDestroy(&pc_ctx->pc)); /* decrease PCBDDC reference count */ 1109566063dSJacob Faibussowitsch PetscCall(VecDestroy(&pc_ctx->xPg)); 1119566063dSJacob Faibussowitsch PetscCall(VecDestroy(&pc_ctx->yPg)); 1129566063dSJacob Faibussowitsch PetscCall(PetscFree(pc_ctx)); 113674ae819SStefano Zampini PetscFunctionReturn(0); 114674ae819SStefano Zampini } 115674ae819SStefano Zampini 116674ae819SStefano Zampini PetscErrorCode PCBDDCSetupFETIDPMatContext(FETIDPMat_ctx fetidpmat_ctx) 117674ae819SStefano Zampini { 118674ae819SStefano Zampini PC_IS *pcis=(PC_IS*)fetidpmat_ctx->pc->data; 119674ae819SStefano Zampini PC_BDDC *pcbddc=(PC_BDDC*)fetidpmat_ctx->pc->data; 120674ae819SStefano Zampini PCBDDCGraph mat_graph=pcbddc->mat_graph; 121674ae819SStefano Zampini Mat_IS *matis = (Mat_IS*)fetidpmat_ctx->pc->pmat->data; 122674ae819SStefano Zampini MPI_Comm comm; 12325d06dbeSstefano_zampini Mat ScalingMat,BD1,BD2; 124457d4a33Sstefano_zampini Vec fetidp_global; 125674ae819SStefano Zampini IS IS_l2g_lambda; 126a1c0d0daSstefano_zampini IS subset,subset_mult,subset_n,isvert; 127674ae819SStefano Zampini PetscBool skip_node,fully_redundant; 128674ae819SStefano Zampini PetscInt i,j,k,s,n_boundary_dofs,n_global_lambda,n_vertices,partial_sum; 129dc456d91SStefano Zampini PetscInt cum,n_local_lambda,n_lambda_for_dof,dual_size,n_neg_values,n_pos_values; 13076ec1555SBarry Smith PetscMPIInt rank,size,buf_size,neigh; 131674ae819SStefano Zampini PetscScalar scalar_value; 132a1c0d0daSstefano_zampini const PetscInt *vertex_indices; 133dc456d91SStefano Zampini PetscInt *dual_dofs_boundary_indices,*aux_local_numbering_1; 134dc456d91SStefano Zampini const PetscInt *aux_global_numbering; 135674ae819SStefano Zampini PetscInt *aux_sums,*cols_B_delta,*l2g_indices; 136674ae819SStefano Zampini PetscScalar *array,*scaling_factors,*vals_B_delta; 13725d06dbeSstefano_zampini PetscScalar **all_factors; 138674ae819SStefano Zampini PetscInt *aux_local_numbering_2; 139457d4a33Sstefano_zampini PetscLayout llay; 140a1c0d0daSstefano_zampini 141457d4a33Sstefano_zampini /* saddlepoint */ 142457d4a33Sstefano_zampini ISLocalToGlobalMapping l2gmap_p; 143457d4a33Sstefano_zampini PetscLayout play; 144457d4a33Sstefano_zampini IS gP,pP; 145457d4a33Sstefano_zampini PetscInt nPl,nPg,nPgl; 146674ae819SStefano Zampini 147674ae819SStefano Zampini PetscFunctionBegin; 1489566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)(fetidpmat_ctx->pc),&comm)); 1499566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm,&rank)); 1509566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(comm,&size)); 151674ae819SStefano Zampini 152457d4a33Sstefano_zampini /* saddlepoint */ 153457d4a33Sstefano_zampini nPl = 0; 154457d4a33Sstefano_zampini nPg = 0; 155457d4a33Sstefano_zampini nPgl = 0; 156457d4a33Sstefano_zampini gP = NULL; 157457d4a33Sstefano_zampini pP = NULL; 158457d4a33Sstefano_zampini l2gmap_p = NULL; 159457d4a33Sstefano_zampini play = NULL; 1609566063dSJacob Faibussowitsch PetscCall(PetscObjectQuery((PetscObject)fetidpmat_ctx->pc,"__KSPFETIDP_pP",(PetscObject*)&pP)); 161022d8d2bSstefano_zampini if (pP) { /* saddle point */ 162457d4a33Sstefano_zampini /* subdomain pressures in global numbering */ 1639566063dSJacob Faibussowitsch PetscCall(PetscObjectQuery((PetscObject)fetidpmat_ctx->pc,"__KSPFETIDP_gP",(PetscObject*)&gP)); 16428b400f6SJacob Faibussowitsch PetscCheck(gP,PETSC_COMM_SELF,PETSC_ERR_PLIB,"gP not present"); 1659566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(gP,&nPl)); 1669566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_SELF,&fetidpmat_ctx->vP)); 1679566063dSJacob Faibussowitsch PetscCall(VecSetSizes(fetidpmat_ctx->vP,nPl,nPl)); 1689566063dSJacob Faibussowitsch PetscCall(VecSetType(fetidpmat_ctx->vP,VECSTANDARD)); 1699566063dSJacob Faibussowitsch PetscCall(VecSetUp(fetidpmat_ctx->vP)); 170457d4a33Sstefano_zampini 171e1214c54Sstefano_zampini /* pressure matrix */ 1729566063dSJacob Faibussowitsch PetscCall(PetscObjectQuery((PetscObject)fetidpmat_ctx->pc,"__KSPFETIDP_C",(PetscObject*)&fetidpmat_ctx->C)); 173e1214c54Sstefano_zampini if (!fetidpmat_ctx->C) { /* null pressure block, compute layout and global numbering for pressures */ 174457d4a33Sstefano_zampini IS Pg; 175457d4a33Sstefano_zampini 1769566063dSJacob Faibussowitsch PetscCall(ISRenumber(gP,NULL,&nPg,&Pg)); 1779566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingCreateIS(Pg,&l2gmap_p)); 1789566063dSJacob Faibussowitsch PetscCall(ISDestroy(&Pg)); 1799566063dSJacob Faibussowitsch PetscCall(PetscLayoutCreate(comm,&play)); 1809566063dSJacob Faibussowitsch PetscCall(PetscLayoutSetBlockSize(play,1)); 1819566063dSJacob Faibussowitsch PetscCall(PetscLayoutSetSize(play,nPg)); 1829566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(pP,&nPgl)); 1839566063dSJacob Faibussowitsch PetscCall(PetscLayoutSetLocalSize(play,nPgl)); 1849566063dSJacob Faibussowitsch PetscCall(PetscLayoutSetUp(play)); 185457d4a33Sstefano_zampini } else { 1869566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)fetidpmat_ctx->C)); 1879566063dSJacob Faibussowitsch PetscCall(MatISGetLocalToGlobalMapping(fetidpmat_ctx->C,&l2gmap_p,NULL)); 1889566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)l2gmap_p)); 1899566063dSJacob Faibussowitsch PetscCall(MatGetSize(fetidpmat_ctx->C,&nPg,NULL)); 1909566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(fetidpmat_ctx->C,NULL,&nPgl)); 1919566063dSJacob Faibussowitsch PetscCall(MatGetLayouts(fetidpmat_ctx->C,NULL,&llay)); 1929566063dSJacob Faibussowitsch PetscCall(PetscLayoutReference(llay,&play)); 193457d4a33Sstefano_zampini } 1949566063dSJacob Faibussowitsch PetscCall(VecCreateMPIWithArray(comm,1,nPgl,nPg,NULL,&fetidpmat_ctx->xPg)); 1959566063dSJacob Faibussowitsch PetscCall(VecCreateMPIWithArray(comm,1,nPgl,nPg,NULL,&fetidpmat_ctx->yPg)); 196457d4a33Sstefano_zampini 197e1214c54Sstefano_zampini /* import matrices for pressures coupling */ 1989566063dSJacob Faibussowitsch PetscCall(PetscObjectQuery((PetscObject)fetidpmat_ctx->pc,"__KSPFETIDP_B_BI",(PetscObject*)&fetidpmat_ctx->B_BI)); 19928b400f6SJacob Faibussowitsch PetscCheck(fetidpmat_ctx->B_BI,PETSC_COMM_SELF,PETSC_ERR_PLIB,"B_BI not present"); 2009566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)fetidpmat_ctx->B_BI)); 201457d4a33Sstefano_zampini 2029566063dSJacob Faibussowitsch PetscCall(PetscObjectQuery((PetscObject)fetidpmat_ctx->pc,"__KSPFETIDP_B_BB",(PetscObject*)&fetidpmat_ctx->B_BB)); 20328b400f6SJacob Faibussowitsch PetscCheck(fetidpmat_ctx->B_BB,PETSC_COMM_SELF,PETSC_ERR_PLIB,"B_BB not present"); 2049566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)fetidpmat_ctx->B_BB)); 205457d4a33Sstefano_zampini 2069566063dSJacob Faibussowitsch PetscCall(PetscObjectQuery((PetscObject)fetidpmat_ctx->pc,"__KSPFETIDP_Bt_BI",(PetscObject*)&fetidpmat_ctx->Bt_BI)); 20728b400f6SJacob Faibussowitsch PetscCheck(fetidpmat_ctx->Bt_BI,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Bt_BI not present"); 2089566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)fetidpmat_ctx->Bt_BI)); 209457d4a33Sstefano_zampini 2109566063dSJacob Faibussowitsch PetscCall(PetscObjectQuery((PetscObject)fetidpmat_ctx->pc,"__KSPFETIDP_Bt_BB",(PetscObject*)&fetidpmat_ctx->Bt_BB)); 21128b400f6SJacob Faibussowitsch PetscCheck(fetidpmat_ctx->Bt_BB,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Bt_BB not present"); 2129566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)fetidpmat_ctx->Bt_BB)); 2136cc1294bSstefano_zampini 2149566063dSJacob Faibussowitsch PetscCall(PetscObjectQuery((PetscObject)fetidpmat_ctx->pc,"__KSPFETIDP_flip" ,(PetscObject*)&fetidpmat_ctx->rhs_flip)); 2156cc1294bSstefano_zampini if (fetidpmat_ctx->rhs_flip) { 2169566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)fetidpmat_ctx->rhs_flip)); 2176cc1294bSstefano_zampini } 218457d4a33Sstefano_zampini } 219457d4a33Sstefano_zampini 220674ae819SStefano Zampini /* Default type of lagrange multipliers is non-redundant */ 221329cd39dSStefano Zampini fully_redundant = fetidpmat_ctx->fully_redundant; 222674ae819SStefano Zampini 223674ae819SStefano Zampini /* Evaluate local and global number of lagrange multipliers */ 2249566063dSJacob Faibussowitsch PetscCall(VecSet(pcis->vec1_N,0.0)); 225674ae819SStefano Zampini n_local_lambda = 0; 226674ae819SStefano Zampini partial_sum = 0; 227674ae819SStefano Zampini n_boundary_dofs = 0; 228674ae819SStefano Zampini s = 0; 229a1c0d0daSstefano_zampini 230674ae819SStefano Zampini /* Get Vertices used to define the BDDC */ 2319566063dSJacob Faibussowitsch PetscCall(PCBDDCGraphGetCandidatesIS(pcbddc->mat_graph,NULL,NULL,NULL,NULL,&isvert)); 2329566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(isvert,&n_vertices)); 2339566063dSJacob Faibussowitsch PetscCall(ISGetIndices(isvert,&vertex_indices)); 234a1c0d0daSstefano_zampini 235674ae819SStefano Zampini dual_size = pcis->n_B-n_vertices; 2369566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(dual_size,&dual_dofs_boundary_indices)); 2379566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(dual_size,&aux_local_numbering_1)); 2389566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(dual_size,&aux_local_numbering_2)); 239674ae819SStefano Zampini 2409566063dSJacob Faibussowitsch PetscCall(VecGetArray(pcis->vec1_N,&array)); 241674ae819SStefano Zampini for (i=0;i<pcis->n;i++) { 242674ae819SStefano Zampini j = mat_graph->count[i]; /* RECALL: mat_graph->count[i] does not count myself */ 2432d456bbaSstefano_zampini if (j > 0) n_boundary_dofs++; 244674ae819SStefano Zampini skip_node = PETSC_FALSE; 245674ae819SStefano Zampini if (s < n_vertices && vertex_indices[s] == i) { /* it works for a sorted set of vertices */ 246674ae819SStefano Zampini skip_node = PETSC_TRUE; 247674ae819SStefano Zampini s++; 248674ae819SStefano Zampini } 2492d456bbaSstefano_zampini if (j < 1) skip_node = PETSC_TRUE; 2502d456bbaSstefano_zampini if (mat_graph->special_dof[i] == PCBDDCGRAPH_DIRICHLET_MARK) skip_node = PETSC_TRUE; 251674ae819SStefano Zampini if (!skip_node) { 252674ae819SStefano Zampini if (fully_redundant) { 253674ae819SStefano Zampini /* fully redundant set of lagrange multipliers */ 254674ae819SStefano Zampini n_lambda_for_dof = (j*(j+1))/2; 255674ae819SStefano Zampini } else { 256674ae819SStefano Zampini n_lambda_for_dof = j; 257674ae819SStefano Zampini } 258674ae819SStefano Zampini n_local_lambda += j; 259674ae819SStefano Zampini /* needed to evaluate global number of lagrange multipliers */ 260674ae819SStefano Zampini array[i]=(1.0*n_lambda_for_dof)/(j+1.0); /* already scaled for the next global sum */ 261674ae819SStefano Zampini /* store some data needed */ 262674ae819SStefano Zampini dual_dofs_boundary_indices[partial_sum] = n_boundary_dofs-1; 263674ae819SStefano Zampini aux_local_numbering_1[partial_sum] = i; 264674ae819SStefano Zampini aux_local_numbering_2[partial_sum] = n_lambda_for_dof; 265674ae819SStefano Zampini partial_sum++; 266674ae819SStefano Zampini } 267674ae819SStefano Zampini } 2689566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(pcis->vec1_N,&array)); 2699566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(isvert,&vertex_indices)); 2709566063dSJacob Faibussowitsch PetscCall(PCBDDCGraphRestoreCandidatesIS(pcbddc->mat_graph,NULL,NULL,NULL,NULL,&isvert)); 2712d456bbaSstefano_zampini dual_size = partial_sum; 272674ae819SStefano Zampini 273674ae819SStefano Zampini /* compute global ordering of lagrange multipliers and associate l2g map */ 2749566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(comm,partial_sum,aux_local_numbering_1,PETSC_COPY_VALUES,&subset_n)); 2759566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingApplyIS(pcis->mapping,subset_n,&subset)); 2769566063dSJacob Faibussowitsch PetscCall(ISDestroy(&subset_n)); 2779566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(comm,partial_sum,aux_local_numbering_2,PETSC_OWN_POINTER,&subset_mult)); 2789566063dSJacob Faibussowitsch PetscCall(ISRenumber(subset,subset_mult,&fetidpmat_ctx->n_lambda,&subset_n)); 2799566063dSJacob Faibussowitsch PetscCall(ISDestroy(&subset)); 2803d996552SStefano Zampini 28176bd3646SJed Brown if (PetscDefined(USE_DEBUG)) { 2829566063dSJacob Faibussowitsch PetscCall(VecSet(pcis->vec1_global,0.0)); 2839566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(matis->rctx,pcis->vec1_N,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE)); 2849566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(matis->rctx,pcis->vec1_N,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE)); 2859566063dSJacob Faibussowitsch PetscCall(VecSum(pcis->vec1_global,&scalar_value)); 2864fe826edSStefano Zampini i = (PetscInt)PetscRealPart(scalar_value); 28763a3b9bcSJacob Faibussowitsch PetscCheck(i == fetidpmat_ctx->n_lambda,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"Global number of multipliers mismatch! (%" PetscInt_FMT " != %" PetscInt_FMT ")",fetidpmat_ctx->n_lambda,i); 28876bd3646SJed Brown } 289674ae819SStefano Zampini 290674ae819SStefano Zampini /* init data for scaling factors exchange */ 29125d06dbeSstefano_zampini if (!pcbddc->use_deluxe_scaling) { 29225d06dbeSstefano_zampini PetscInt *ptrs_buffer,neigh_position; 29325d06dbeSstefano_zampini PetscScalar *send_buffer,*recv_buffer; 29425d06dbeSstefano_zampini MPI_Request *send_reqs,*recv_reqs; 29525d06dbeSstefano_zampini 296674ae819SStefano Zampini partial_sum = 0; 2979566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(pcis->n_neigh,&ptrs_buffer)); 2989566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(PetscMax(pcis->n_neigh-1,0),&send_reqs)); 2999566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(PetscMax(pcis->n_neigh-1,0),&recv_reqs)); 3009566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(pcis->n+1,&all_factors)); 3014b2aedd3SStefano Zampini if (pcis->n_neigh > 0) ptrs_buffer[0]=0; 302674ae819SStefano Zampini for (i=1;i<pcis->n_neigh;i++) { 303674ae819SStefano Zampini partial_sum += pcis->n_shared[i]; 304674ae819SStefano Zampini ptrs_buffer[i] = ptrs_buffer[i-1]+pcis->n_shared[i]; 305674ae819SStefano Zampini } 3069566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(partial_sum,&send_buffer)); 3079566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(partial_sum,&recv_buffer)); 3089566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(partial_sum,&all_factors[0])); 309674ae819SStefano Zampini for (i=0;i<pcis->n-1;i++) { 310674ae819SStefano Zampini j = mat_graph->count[i]; 311674ae819SStefano Zampini all_factors[i+1]=all_factors[i]+j; 312674ae819SStefano Zampini } 31325d06dbeSstefano_zampini 314674ae819SStefano Zampini /* scatter B scaling to N vec */ 3159566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(pcis->N_to_B,pcis->D,pcis->vec1_N,INSERT_VALUES,SCATTER_REVERSE)); 3169566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(pcis->N_to_B,pcis->D,pcis->vec1_N,INSERT_VALUES,SCATTER_REVERSE)); 317674ae819SStefano Zampini /* communications */ 3189566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(pcis->vec1_N,(const PetscScalar**)&array)); 319674ae819SStefano Zampini for (i=1;i<pcis->n_neigh;i++) { 320674ae819SStefano Zampini for (j=0;j<pcis->n_shared[i];j++) { 321674ae819SStefano Zampini send_buffer[ptrs_buffer[i-1]+j]=array[pcis->shared[i][j]]; 322674ae819SStefano Zampini } 3239566063dSJacob Faibussowitsch PetscCall(PetscMPIIntCast(ptrs_buffer[i]-ptrs_buffer[i-1],&buf_size)); 3249566063dSJacob Faibussowitsch PetscCall(PetscMPIIntCast(pcis->neigh[i],&neigh)); 3259566063dSJacob Faibussowitsch PetscCallMPI(MPI_Isend(&send_buffer[ptrs_buffer[i-1]],buf_size,MPIU_SCALAR,neigh,0,comm,&send_reqs[i-1])); 3269566063dSJacob Faibussowitsch PetscCallMPI(MPI_Irecv(&recv_buffer[ptrs_buffer[i-1]],buf_size,MPIU_SCALAR,neigh,0,comm,&recv_reqs[i-1])); 327674ae819SStefano Zampini } 3289566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(pcis->vec1_N,(const PetscScalar**)&array)); 3294b2aedd3SStefano Zampini if (pcis->n_neigh > 0) { 3309566063dSJacob Faibussowitsch PetscCallMPI(MPI_Waitall(pcis->n_neigh-1,recv_reqs,MPI_STATUSES_IGNORE)); 3314b2aedd3SStefano Zampini } 332674ae819SStefano Zampini /* put values in correct places */ 333674ae819SStefano Zampini for (i=1;i<pcis->n_neigh;i++) { 334674ae819SStefano Zampini for (j=0;j<pcis->n_shared[i];j++) { 335674ae819SStefano Zampini k = pcis->shared[i][j]; 336674ae819SStefano Zampini neigh_position = 0; 337674ae819SStefano Zampini while (mat_graph->neighbours_set[k][neigh_position] != pcis->neigh[i]) {neigh_position++;} 338674ae819SStefano Zampini all_factors[k][neigh_position]=recv_buffer[ptrs_buffer[i-1]+j]; 339674ae819SStefano Zampini } 340674ae819SStefano Zampini } 3414b2aedd3SStefano Zampini if (pcis->n_neigh > 0) { 3429566063dSJacob Faibussowitsch PetscCallMPI(MPI_Waitall(pcis->n_neigh-1,send_reqs,MPI_STATUSES_IGNORE)); 3434b2aedd3SStefano Zampini } 3449566063dSJacob Faibussowitsch PetscCall(PetscFree(send_reqs)); 3459566063dSJacob Faibussowitsch PetscCall(PetscFree(recv_reqs)); 3469566063dSJacob Faibussowitsch PetscCall(PetscFree(send_buffer)); 3479566063dSJacob Faibussowitsch PetscCall(PetscFree(recv_buffer)); 3489566063dSJacob Faibussowitsch PetscCall(PetscFree(ptrs_buffer)); 34925d06dbeSstefano_zampini } 350674ae819SStefano Zampini 351674ae819SStefano Zampini /* Compute B and B_delta (local actions) */ 3529566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(pcis->n_neigh,&aux_sums)); 3539566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n_local_lambda,&l2g_indices)); 3549566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n_local_lambda,&vals_B_delta)); 3559566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n_local_lambda,&cols_B_delta)); 35625d06dbeSstefano_zampini if (!pcbddc->use_deluxe_scaling) { 3579566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n_local_lambda,&scaling_factors)); 35825d06dbeSstefano_zampini } else { 35925d06dbeSstefano_zampini scaling_factors = NULL; 36025d06dbeSstefano_zampini all_factors = NULL; 36125d06dbeSstefano_zampini } 3629566063dSJacob Faibussowitsch PetscCall(ISGetIndices(subset_n,&aux_global_numbering)); 363674ae819SStefano Zampini partial_sum=0; 364dc456d91SStefano Zampini cum = 0; 365674ae819SStefano Zampini for (i=0;i<dual_size;i++) { 366dc456d91SStefano Zampini n_global_lambda = aux_global_numbering[cum]; 367674ae819SStefano Zampini j = mat_graph->count[aux_local_numbering_1[i]]; 368674ae819SStefano Zampini aux_sums[0]=0; 369674ae819SStefano Zampini for (s=1;s<j;s++) { 370674ae819SStefano Zampini aux_sums[s]=aux_sums[s-1]+j-s+1; 371674ae819SStefano Zampini } 37225d06dbeSstefano_zampini if (all_factors) array = all_factors[aux_local_numbering_1[i]]; 373674ae819SStefano Zampini n_neg_values = 0; 3742a7da448SStefano Zampini while (n_neg_values < j && mat_graph->neighbours_set[aux_local_numbering_1[i]][n_neg_values] < rank) {n_neg_values++;} 375674ae819SStefano Zampini n_pos_values = j - n_neg_values; 376674ae819SStefano Zampini if (fully_redundant) { 377674ae819SStefano Zampini for (s=0;s<n_neg_values;s++) { 378674ae819SStefano Zampini l2g_indices [partial_sum+s]=aux_sums[s]+n_neg_values-s-1+n_global_lambda; 379674ae819SStefano Zampini cols_B_delta [partial_sum+s]=dual_dofs_boundary_indices[i]; 380674ae819SStefano Zampini vals_B_delta [partial_sum+s]=-1.0; 38125d06dbeSstefano_zampini if (!pcbddc->use_deluxe_scaling) scaling_factors[partial_sum+s]=array[s]; 382674ae819SStefano Zampini } 383674ae819SStefano Zampini for (s=0;s<n_pos_values;s++) { 384674ae819SStefano Zampini l2g_indices [partial_sum+s+n_neg_values]=aux_sums[n_neg_values]+s+n_global_lambda; 385674ae819SStefano Zampini cols_B_delta [partial_sum+s+n_neg_values]=dual_dofs_boundary_indices[i]; 386674ae819SStefano Zampini vals_B_delta [partial_sum+s+n_neg_values]=1.0; 38725d06dbeSstefano_zampini if (!pcbddc->use_deluxe_scaling) scaling_factors[partial_sum+s+n_neg_values]=array[s+n_neg_values]; 388674ae819SStefano Zampini } 389674ae819SStefano Zampini partial_sum += j; 390674ae819SStefano Zampini } else { 391674ae819SStefano Zampini /* l2g_indices and default cols and vals of B_delta */ 392674ae819SStefano Zampini for (s=0;s<j;s++) { 393674ae819SStefano Zampini l2g_indices [partial_sum+s]=n_global_lambda+s; 394674ae819SStefano Zampini cols_B_delta [partial_sum+s]=dual_dofs_boundary_indices[i]; 395674ae819SStefano Zampini vals_B_delta [partial_sum+s]=0.0; 396674ae819SStefano Zampini } 397674ae819SStefano Zampini /* B_delta */ 398674ae819SStefano Zampini if (n_neg_values > 0) { /* there's a rank next to me to the left */ 399674ae819SStefano Zampini vals_B_delta [partial_sum+n_neg_values-1]=-1.0; 400674ae819SStefano Zampini } 401674ae819SStefano Zampini if (n_neg_values < j) { /* there's a rank next to me to the right */ 402674ae819SStefano Zampini vals_B_delta [partial_sum+n_neg_values]=1.0; 403674ae819SStefano Zampini } 404674ae819SStefano Zampini /* scaling as in Klawonn-Widlund 1999 */ 40525d06dbeSstefano_zampini if (!pcbddc->use_deluxe_scaling) { 406674ae819SStefano Zampini for (s=0;s<n_neg_values;s++) { 407674ae819SStefano Zampini scalar_value = 0.0; 40825d06dbeSstefano_zampini for (k=0;k<s+1;k++) scalar_value += array[k]; 409674ae819SStefano Zampini scaling_factors[partial_sum+s] = -scalar_value; 410674ae819SStefano Zampini } 411674ae819SStefano Zampini for (s=0;s<n_pos_values;s++) { 412674ae819SStefano Zampini scalar_value = 0.0; 41325d06dbeSstefano_zampini for (k=s+n_neg_values;k<j;k++) scalar_value += array[k]; 414674ae819SStefano Zampini scaling_factors[partial_sum+s+n_neg_values] = scalar_value; 415674ae819SStefano Zampini } 41625d06dbeSstefano_zampini } 417674ae819SStefano Zampini partial_sum += j; 418674ae819SStefano Zampini } 419dc456d91SStefano Zampini cum += aux_local_numbering_2[i]; 420674ae819SStefano Zampini } 4219566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(subset_n,&aux_global_numbering)); 4229566063dSJacob Faibussowitsch PetscCall(ISDestroy(&subset_mult)); 4239566063dSJacob Faibussowitsch PetscCall(ISDestroy(&subset_n)); 4249566063dSJacob Faibussowitsch PetscCall(PetscFree(aux_sums)); 4259566063dSJacob Faibussowitsch PetscCall(PetscFree(aux_local_numbering_1)); 4269566063dSJacob Faibussowitsch PetscCall(PetscFree(dual_dofs_boundary_indices)); 42725d06dbeSstefano_zampini if (all_factors) { 4289566063dSJacob Faibussowitsch PetscCall(PetscFree(all_factors[0])); 4299566063dSJacob Faibussowitsch PetscCall(PetscFree(all_factors)); 43025d06dbeSstefano_zampini } 431674ae819SStefano Zampini 432674ae819SStefano Zampini /* Create local part of B_delta */ 4339566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_SELF,&fetidpmat_ctx->B_delta)); 4349566063dSJacob Faibussowitsch PetscCall(MatSetSizes(fetidpmat_ctx->B_delta,n_local_lambda,pcis->n_B,n_local_lambda,pcis->n_B)); 4359566063dSJacob Faibussowitsch PetscCall(MatSetType(fetidpmat_ctx->B_delta,MATSEQAIJ)); 4369566063dSJacob Faibussowitsch PetscCall(MatSeqAIJSetPreallocation(fetidpmat_ctx->B_delta,1,NULL)); 4379566063dSJacob Faibussowitsch PetscCall(MatSetOption(fetidpmat_ctx->B_delta,MAT_IGNORE_ZERO_ENTRIES,PETSC_TRUE)); 438674ae819SStefano Zampini for (i=0;i<n_local_lambda;i++) { 4399566063dSJacob Faibussowitsch PetscCall(MatSetValue(fetidpmat_ctx->B_delta,i,cols_B_delta[i],vals_B_delta[i],INSERT_VALUES)); 440674ae819SStefano Zampini } 4419566063dSJacob Faibussowitsch PetscCall(PetscFree(vals_B_delta)); 4429566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(fetidpmat_ctx->B_delta,MAT_FINAL_ASSEMBLY)); 4439566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(fetidpmat_ctx->B_delta,MAT_FINAL_ASSEMBLY)); 444674ae819SStefano Zampini 44525d06dbeSstefano_zampini BD1 = NULL; 44625d06dbeSstefano_zampini BD2 = NULL; 447674ae819SStefano Zampini if (fully_redundant) { 44828b400f6SJacob Faibussowitsch PetscCheck(!pcbddc->use_deluxe_scaling,comm,PETSC_ERR_SUP,"Deluxe FETIDP with fully-redundant multipliers to be implemented"); 4499566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_SELF,&ScalingMat)); 4509566063dSJacob Faibussowitsch PetscCall(MatSetSizes(ScalingMat,n_local_lambda,n_local_lambda,n_local_lambda,n_local_lambda)); 4519566063dSJacob Faibussowitsch PetscCall(MatSetType(ScalingMat,MATSEQAIJ)); 4529566063dSJacob Faibussowitsch PetscCall(MatSeqAIJSetPreallocation(ScalingMat,1,NULL)); 453674ae819SStefano Zampini for (i=0;i<n_local_lambda;i++) { 4549566063dSJacob Faibussowitsch PetscCall(MatSetValue(ScalingMat,i,i,scaling_factors[i],INSERT_VALUES)); 455674ae819SStefano Zampini } 4569566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(ScalingMat,MAT_FINAL_ASSEMBLY)); 4579566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(ScalingMat,MAT_FINAL_ASSEMBLY)); 4589566063dSJacob Faibussowitsch PetscCall(MatMatMult(ScalingMat,fetidpmat_ctx->B_delta,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&fetidpmat_ctx->B_Ddelta)); 4599566063dSJacob Faibussowitsch PetscCall(MatDestroy(&ScalingMat)); 460674ae819SStefano Zampini } else { 4619566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_SELF,&fetidpmat_ctx->B_Ddelta)); 4629566063dSJacob Faibussowitsch PetscCall(MatSetSizes(fetidpmat_ctx->B_Ddelta,n_local_lambda,pcis->n_B,n_local_lambda,pcis->n_B)); 463524221abSStefano Zampini if (!pcbddc->use_deluxe_scaling || !pcbddc->sub_schurs) { 4649566063dSJacob Faibussowitsch PetscCall(MatSetType(fetidpmat_ctx->B_Ddelta,MATSEQAIJ)); 4659566063dSJacob Faibussowitsch PetscCall(MatSeqAIJSetPreallocation(fetidpmat_ctx->B_Ddelta,1,NULL)); 466674ae819SStefano Zampini for (i=0;i<n_local_lambda;i++) { 4679566063dSJacob Faibussowitsch PetscCall(MatSetValue(fetidpmat_ctx->B_Ddelta,i,cols_B_delta[i],scaling_factors[i],INSERT_VALUES)); 468674ae819SStefano Zampini } 4699566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(fetidpmat_ctx->B_Ddelta,MAT_FINAL_ASSEMBLY)); 4709566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(fetidpmat_ctx->B_Ddelta,MAT_FINAL_ASSEMBLY)); 47125d06dbeSstefano_zampini } else { 47225d06dbeSstefano_zampini /* scaling as in Klawonn-Widlund 1999 */ 47325d06dbeSstefano_zampini PCBDDCDeluxeScaling deluxe_ctx = pcbddc->deluxe_ctx; 47425d06dbeSstefano_zampini PCBDDCSubSchurs sub_schurs = pcbddc->sub_schurs; 47525d06dbeSstefano_zampini Mat T; 47625d06dbeSstefano_zampini PetscScalar *W,lwork,*Bwork; 477451a39c7SStefano Zampini const PetscInt *idxs = NULL; 47825d06dbeSstefano_zampini PetscInt cum,mss,*nnz; 47925d06dbeSstefano_zampini PetscBLASInt *pivots,B_lwork,B_N,B_ierr; 48025d06dbeSstefano_zampini 48128b400f6SJacob Faibussowitsch PetscCheck(pcbddc->deluxe_singlemat,comm,PETSC_ERR_USER,"Cannot compute B_Ddelta! rerun with -pc_bddc_deluxe_singlemat"); 48225d06dbeSstefano_zampini mss = 0; 4839566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(pcis->n_B,&nnz)); 48425d06dbeSstefano_zampini if (sub_schurs->is_Ej_all) { 4859566063dSJacob Faibussowitsch PetscCall(ISGetIndices(sub_schurs->is_Ej_all,&idxs)); 48625d06dbeSstefano_zampini for (i=0,cum=0;i<sub_schurs->n_subs;i++) { 48725d06dbeSstefano_zampini PetscInt subset_size; 48825d06dbeSstefano_zampini 4899566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(sub_schurs->is_subs[i],&subset_size)); 49025d06dbeSstefano_zampini for (j=0;j<subset_size;j++) nnz[idxs[j+cum]] = subset_size; 49125d06dbeSstefano_zampini mss = PetscMax(mss,subset_size); 49225d06dbeSstefano_zampini cum += subset_size; 49325d06dbeSstefano_zampini } 49425d06dbeSstefano_zampini } 4959566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_SELF,&T)); 4969566063dSJacob Faibussowitsch PetscCall(MatSetSizes(T,pcis->n_B,pcis->n_B,pcis->n_B,pcis->n_B)); 4979566063dSJacob Faibussowitsch PetscCall(MatSetType(T,MATSEQAIJ)); 4989566063dSJacob Faibussowitsch PetscCall(MatSeqAIJSetPreallocation(T,0,nnz)); 4999566063dSJacob Faibussowitsch PetscCall(PetscFree(nnz)); 50025d06dbeSstefano_zampini 50125d06dbeSstefano_zampini /* workspace allocation */ 5027f00194aSstefano_zampini B_lwork = 0; 5037f00194aSstefano_zampini if (mss) { 504be8bcbd0Sstefano_zampini PetscScalar dummy = 1; 505be8bcbd0Sstefano_zampini 50625d06dbeSstefano_zampini B_lwork = -1; 5079566063dSJacob Faibussowitsch PetscCall(PetscBLASIntCast(mss,&B_N)); 5089566063dSJacob Faibussowitsch PetscCall(PetscFPTrapPush(PETSC_FP_TRAP_OFF)); 509be8bcbd0Sstefano_zampini PetscStackCallBLAS("LAPACKgetri",LAPACKgetri_(&B_N,&dummy,&B_N,&B_N,&lwork,&B_lwork,&B_ierr)); 5109566063dSJacob Faibussowitsch PetscCall(PetscFPTrapPop()); 51128b400f6SJacob Faibussowitsch PetscCheck(!B_ierr,PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in query to GETRI Lapack routine %d",(int)B_ierr); 5129566063dSJacob Faibussowitsch PetscCall(PetscBLASIntCast((PetscInt)PetscRealPart(lwork),&B_lwork)); 5137f00194aSstefano_zampini } 5149566063dSJacob Faibussowitsch PetscCall(PetscMalloc3(mss*mss,&W,mss,&pivots,B_lwork,&Bwork)); 51525d06dbeSstefano_zampini 51625d06dbeSstefano_zampini for (i=0,cum=0;i<sub_schurs->n_subs;i++) { 5171683a169SBarry Smith const PetscScalar *M; 51825d06dbeSstefano_zampini PetscInt subset_size; 51925d06dbeSstefano_zampini 5209566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(sub_schurs->is_subs[i],&subset_size)); 5219566063dSJacob Faibussowitsch PetscCall(PetscBLASIntCast(subset_size,&B_N)); 5229566063dSJacob Faibussowitsch PetscCall(MatDenseGetArrayRead(deluxe_ctx->seq_mat[i],&M)); 5239566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(W,M,subset_size*subset_size)); 5249566063dSJacob Faibussowitsch PetscCall(MatDenseRestoreArrayRead(deluxe_ctx->seq_mat[i],&M)); 5259566063dSJacob Faibussowitsch PetscCall(PetscFPTrapPush(PETSC_FP_TRAP_OFF)); 52625d06dbeSstefano_zampini PetscStackCallBLAS("LAPACKgetrf",LAPACKgetrf_(&B_N,&B_N,W,&B_N,pivots,&B_ierr)); 52728b400f6SJacob Faibussowitsch PetscCheck(!B_ierr,PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRF Lapack routine %d",(int)B_ierr); 52825d06dbeSstefano_zampini PetscStackCallBLAS("LAPACKgetri",LAPACKgetri_(&B_N,W,&B_N,pivots,Bwork,&B_lwork,&B_ierr)); 52928b400f6SJacob Faibussowitsch PetscCheck(!B_ierr,PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRI Lapack routine %d",(int)B_ierr); 5309566063dSJacob Faibussowitsch PetscCall(PetscFPTrapPop()); 531451a39c7SStefano Zampini /* silent static analyzer */ 53228b400f6SJacob Faibussowitsch PetscCheck(idxs,PETSC_COMM_SELF,PETSC_ERR_PLIB,"IDXS not present"); 5339566063dSJacob Faibussowitsch PetscCall(MatSetValues(T,subset_size,idxs+cum,subset_size,idxs+cum,W,INSERT_VALUES)); 53425d06dbeSstefano_zampini cum += subset_size; 53525d06dbeSstefano_zampini } 5369566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(T,MAT_FINAL_ASSEMBLY)); 5379566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(T,MAT_FINAL_ASSEMBLY)); 5389566063dSJacob Faibussowitsch PetscCall(MatMatTransposeMult(T,fetidpmat_ctx->B_delta,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&BD1)); 5399566063dSJacob Faibussowitsch PetscCall(MatMatMult(fetidpmat_ctx->B_delta,BD1,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&BD2)); 5409566063dSJacob Faibussowitsch PetscCall(MatDestroy(&T)); 5419566063dSJacob Faibussowitsch PetscCall(PetscFree3(W,pivots,Bwork)); 54225d06dbeSstefano_zampini if (sub_schurs->is_Ej_all) { 5439566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(sub_schurs->is_Ej_all,&idxs)); 54425d06dbeSstefano_zampini } 54525d06dbeSstefano_zampini } 546674ae819SStefano Zampini } 5479566063dSJacob Faibussowitsch PetscCall(PetscFree(scaling_factors)); 5489566063dSJacob Faibussowitsch PetscCall(PetscFree(cols_B_delta)); 549674ae819SStefano Zampini 550457d4a33Sstefano_zampini /* Layout of multipliers */ 5519566063dSJacob Faibussowitsch PetscCall(PetscLayoutCreate(comm,&llay)); 5529566063dSJacob Faibussowitsch PetscCall(PetscLayoutSetBlockSize(llay,1)); 5539566063dSJacob Faibussowitsch PetscCall(PetscLayoutSetSize(llay,fetidpmat_ctx->n_lambda)); 5549566063dSJacob Faibussowitsch PetscCall(PetscLayoutSetUp(llay)); 5559566063dSJacob Faibussowitsch PetscCall(PetscLayoutGetLocalSize(llay,&fetidpmat_ctx->n)); 556457d4a33Sstefano_zampini 557457d4a33Sstefano_zampini /* Local work vector of multipliers */ 5589566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_SELF,&fetidpmat_ctx->lambda_local)); 5599566063dSJacob Faibussowitsch PetscCall(VecSetSizes(fetidpmat_ctx->lambda_local,n_local_lambda,n_local_lambda)); 5609566063dSJacob Faibussowitsch PetscCall(VecSetType(fetidpmat_ctx->lambda_local,VECSEQ)); 561457d4a33Sstefano_zampini 56225d06dbeSstefano_zampini if (BD2) { 56325d06dbeSstefano_zampini ISLocalToGlobalMapping l2g; 56425d06dbeSstefano_zampini Mat T,TA,*pT; 56525d06dbeSstefano_zampini IS is; 56625d06dbeSstefano_zampini PetscInt nl,N; 56725d06dbeSstefano_zampini BDdelta_DN ctx; 56825d06dbeSstefano_zampini 5699566063dSJacob Faibussowitsch PetscCall(PetscLayoutGetLocalSize(llay,&nl)); 5709566063dSJacob Faibussowitsch PetscCall(PetscLayoutGetSize(llay,&N)); 5719566063dSJacob Faibussowitsch PetscCall(MatCreate(comm,&T)); 5729566063dSJacob Faibussowitsch PetscCall(MatSetSizes(T,nl,nl,N,N)); 5739566063dSJacob Faibussowitsch PetscCall(MatSetType(T,MATIS)); 5749566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingCreate(comm,1,n_local_lambda,l2g_indices,PETSC_COPY_VALUES,&l2g)); 5759566063dSJacob Faibussowitsch PetscCall(MatSetLocalToGlobalMapping(T,l2g,l2g)); 5769566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingDestroy(&l2g)); 5779566063dSJacob Faibussowitsch PetscCall(MatISSetLocalMat(T,BD2)); 5789566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(T,MAT_FINAL_ASSEMBLY)); 5799566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(T,MAT_FINAL_ASSEMBLY)); 5809566063dSJacob Faibussowitsch PetscCall(MatDestroy(&BD2)); 5819566063dSJacob Faibussowitsch PetscCall(MatConvert(T,MATAIJ,MAT_INITIAL_MATRIX,&TA)); 5829566063dSJacob Faibussowitsch PetscCall(MatDestroy(&T)); 5839566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(comm,n_local_lambda,l2g_indices,PETSC_USE_POINTER,&is)); 5849566063dSJacob Faibussowitsch PetscCall(MatCreateSubMatrices(TA,1,&is,&is,MAT_INITIAL_MATRIX,&pT)); 5859566063dSJacob Faibussowitsch PetscCall(MatDestroy(&TA)); 5869566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is)); 58725d06dbeSstefano_zampini BD2 = pT[0]; 5889566063dSJacob Faibussowitsch PetscCall(PetscFree(pT)); 58925d06dbeSstefano_zampini 59025d06dbeSstefano_zampini /* B_Ddelta for non-redundant multipliers with deluxe scaling */ 5919566063dSJacob Faibussowitsch PetscCall(PetscNew(&ctx)); 5929566063dSJacob Faibussowitsch PetscCall(MatSetType(fetidpmat_ctx->B_Ddelta,MATSHELL)); 5939566063dSJacob Faibussowitsch PetscCall(MatShellSetContext(fetidpmat_ctx->B_Ddelta,ctx)); 5949566063dSJacob Faibussowitsch PetscCall(MatShellSetOperation(fetidpmat_ctx->B_Ddelta,MATOP_MULT,(void (*)(void))MatMult_BDdelta_deluxe_nonred)); 5959566063dSJacob Faibussowitsch PetscCall(MatShellSetOperation(fetidpmat_ctx->B_Ddelta,MATOP_MULT_TRANSPOSE,(void (*)(void))MatMultTranspose_BDdelta_deluxe_nonred)); 5969566063dSJacob Faibussowitsch PetscCall(MatShellSetOperation(fetidpmat_ctx->B_Ddelta,MATOP_DESTROY,(void (*)(void))MatDestroy_BDdelta_deluxe_nonred)); 5979566063dSJacob Faibussowitsch PetscCall(MatSetUp(fetidpmat_ctx->B_Ddelta)); 59825d06dbeSstefano_zampini 5999566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)BD1)); 60025d06dbeSstefano_zampini ctx->BD = BD1; 6019566063dSJacob Faibussowitsch PetscCall(KSPCreate(PETSC_COMM_SELF,&ctx->kBD)); 6029566063dSJacob Faibussowitsch PetscCall(KSPSetOperators(ctx->kBD,BD2,BD2)); 6039566063dSJacob Faibussowitsch PetscCall(VecDuplicate(fetidpmat_ctx->lambda_local,&ctx->work)); 60425d06dbeSstefano_zampini fetidpmat_ctx->deluxe_nonred = PETSC_TRUE; 60525d06dbeSstefano_zampini } 6069566063dSJacob Faibussowitsch PetscCall(MatDestroy(&BD1)); 6079566063dSJacob Faibussowitsch PetscCall(MatDestroy(&BD2)); 60825d06dbeSstefano_zampini 60925d06dbeSstefano_zampini /* fetidpmat sizes */ 61025d06dbeSstefano_zampini fetidpmat_ctx->n += nPgl; 61125d06dbeSstefano_zampini fetidpmat_ctx->N = fetidpmat_ctx->n_lambda+nPg; 61225d06dbeSstefano_zampini 613457d4a33Sstefano_zampini /* Global vector for FETI-DP linear system */ 6149566063dSJacob Faibussowitsch PetscCall(VecCreate(comm,&fetidp_global)); 6159566063dSJacob Faibussowitsch PetscCall(VecSetSizes(fetidp_global,fetidpmat_ctx->n,fetidpmat_ctx->N)); 6169566063dSJacob Faibussowitsch PetscCall(VecSetType(fetidp_global,VECMPI)); 6179566063dSJacob Faibussowitsch PetscCall(VecSetUp(fetidp_global)); 618457d4a33Sstefano_zampini 6199eec4de8Sstefano_zampini /* Decide layout for fetidp dofs: if it is a saddle point problem 6209eec4de8Sstefano_zampini pressure is ordered first in the local part of the global vector 621457d4a33Sstefano_zampini of the FETI-DP linear system */ 622457d4a33Sstefano_zampini if (nPg) { 623e1214c54Sstefano_zampini Vec v; 624af140850Sstefano_zampini IS IS_l2g_p,ais; 625457d4a33Sstefano_zampini PetscLayout alay; 626457d4a33Sstefano_zampini const PetscInt *idxs,*pranges,*aranges,*lranges; 627af140850Sstefano_zampini PetscInt *l2g_indices_p,rst; 628e1214c54Sstefano_zampini PetscMPIInt rank; 629457d4a33Sstefano_zampini 6309566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nPl,&l2g_indices_p)); 6319566063dSJacob Faibussowitsch PetscCall(VecGetLayout(fetidp_global,&alay)); 6329566063dSJacob Faibussowitsch PetscCall(PetscLayoutGetRanges(alay,&aranges)); 6339566063dSJacob Faibussowitsch PetscCall(PetscLayoutGetRanges(play,&pranges)); 6349566063dSJacob Faibussowitsch PetscCall(PetscLayoutGetRanges(llay,&lranges)); 635e1214c54Sstefano_zampini 6369566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)fetidp_global),&rank)); 6379566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PetscObjectComm((PetscObject)fetidp_global),pranges[rank+1]-pranges[rank],aranges[rank],1,&fetidpmat_ctx->pressure)); 6389566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)fetidpmat_ctx->pressure,"F_P")); 6399566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PetscObjectComm((PetscObject)fetidp_global),lranges[rank+1]-lranges[rank],aranges[rank]+pranges[rank+1]-pranges[rank],1,&fetidpmat_ctx->lagrange)); 6409566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)fetidpmat_ctx->lagrange,"F_L")); 6419566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingGetIndices(l2gmap_p,&idxs)); 642af140850Sstefano_zampini /* shift local to global indices for pressure */ 643457d4a33Sstefano_zampini for (i=0;i<nPl;i++) { 644131c27b5Sprj- PetscMPIInt owner; 645457d4a33Sstefano_zampini 6469566063dSJacob Faibussowitsch PetscCall(PetscLayoutFindOwner(play,idxs[i],&owner)); 647457d4a33Sstefano_zampini l2g_indices_p[i] = idxs[i]-pranges[owner]+aranges[owner]; 648457d4a33Sstefano_zampini } 6499566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingRestoreIndices(l2gmap_p,&idxs)); 6509566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(comm,nPl,l2g_indices_p,PETSC_OWN_POINTER,&IS_l2g_p)); 651af140850Sstefano_zampini 652e1214c54Sstefano_zampini /* local to global scatter for pressure */ 6539566063dSJacob Faibussowitsch PetscCall(VecScatterCreate(fetidpmat_ctx->vP,NULL,fetidp_global,IS_l2g_p,&fetidpmat_ctx->l2g_p)); 6549566063dSJacob Faibussowitsch PetscCall(ISDestroy(&IS_l2g_p)); 655457d4a33Sstefano_zampini 656e1214c54Sstefano_zampini /* scatter for lagrange multipliers only */ 6579566063dSJacob Faibussowitsch PetscCall(VecCreate(comm,&v)); 6589566063dSJacob Faibussowitsch PetscCall(VecSetType(v,VECSTANDARD)); 6599566063dSJacob Faibussowitsch PetscCall(VecSetLayout(v,llay)); 6609566063dSJacob Faibussowitsch PetscCall(VecSetUp(v)); 6619566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(comm,n_local_lambda,l2g_indices,PETSC_COPY_VALUES,&ais)); 6629566063dSJacob Faibussowitsch PetscCall(VecScatterCreate(fetidpmat_ctx->lambda_local,NULL,v,ais,&fetidpmat_ctx->l2g_lambda_only)); 6639566063dSJacob Faibussowitsch PetscCall(ISDestroy(&ais)); 6649566063dSJacob Faibussowitsch PetscCall(VecDestroy(&v)); 665e1214c54Sstefano_zampini 666af140850Sstefano_zampini /* shift local to global indices for multipliers */ 667457d4a33Sstefano_zampini for (i=0;i<n_local_lambda;i++) { 668131c27b5Sprj- PetscInt ps; 669131c27b5Sprj- PetscMPIInt owner; 670457d4a33Sstefano_zampini 6719566063dSJacob Faibussowitsch PetscCall(PetscLayoutFindOwner(llay,l2g_indices[i],&owner)); 672457d4a33Sstefano_zampini ps = pranges[owner+1]-pranges[owner]; 673457d4a33Sstefano_zampini l2g_indices[i] = l2g_indices[i]-lranges[owner]+aranges[owner]+ps; 674457d4a33Sstefano_zampini } 675457d4a33Sstefano_zampini 676e1214c54Sstefano_zampini /* scatter from alldofs to pressures global fetidp vector */ 6779566063dSJacob Faibussowitsch PetscCall(PetscLayoutGetRange(alay,&rst,NULL)); 6789566063dSJacob Faibussowitsch PetscCall(ISCreateStride(comm,nPgl,rst,1,&ais)); 6799566063dSJacob Faibussowitsch PetscCall(VecScatterCreate(pcis->vec1_global,pP,fetidp_global,ais,&fetidpmat_ctx->g2g_p)); 6809566063dSJacob Faibussowitsch PetscCall(ISDestroy(&ais)); 681457d4a33Sstefano_zampini } 6829566063dSJacob Faibussowitsch PetscCall(PetscLayoutDestroy(&llay)); 6839566063dSJacob Faibussowitsch PetscCall(PetscLayoutDestroy(&play)); 6849566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(comm,n_local_lambda,l2g_indices,PETSC_OWN_POINTER,&IS_l2g_lambda)); 685a1c0d0daSstefano_zampini 6869cc7774eSstefano_zampini /* scatter from local to global multipliers */ 6879566063dSJacob Faibussowitsch PetscCall(VecScatterCreate(fetidpmat_ctx->lambda_local,NULL,fetidp_global,IS_l2g_lambda,&fetidpmat_ctx->l2g_lambda)); 6889566063dSJacob Faibussowitsch PetscCall(ISDestroy(&IS_l2g_lambda)); 6899566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingDestroy(&l2gmap_p)); 6909566063dSJacob Faibussowitsch PetscCall(VecDestroy(&fetidp_global)); 691457d4a33Sstefano_zampini 692a1c0d0daSstefano_zampini /* Create some work vectors needed by fetidp */ 6939566063dSJacob Faibussowitsch PetscCall(VecDuplicate(pcis->vec1_B,&fetidpmat_ctx->temp_solution_B)); 6949566063dSJacob Faibussowitsch PetscCall(VecDuplicate(pcis->vec1_D,&fetidpmat_ctx->temp_solution_D)); 695674ae819SStefano Zampini PetscFunctionReturn(0); 696674ae819SStefano Zampini } 697674ae819SStefano Zampini 698674ae819SStefano Zampini PetscErrorCode PCBDDCSetupFETIDPPCContext(Mat fetimat, FETIDPPC_ctx fetidppc_ctx) 699674ae819SStefano Zampini { 700674ae819SStefano Zampini FETIDPMat_ctx mat_ctx; 70125d06dbeSstefano_zampini PC_BDDC *pcbddc = (PC_BDDC*)fetidppc_ctx->pc->data; 70225d06dbeSstefano_zampini PC_IS *pcis = (PC_IS*)fetidppc_ctx->pc->data; 703f28b6018SStefano Zampini PetscBool lumped = PETSC_FALSE; 704674ae819SStefano Zampini 705674ae819SStefano Zampini PetscFunctionBegin; 7069566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(fetimat,&mat_ctx)); 707674ae819SStefano Zampini /* get references from objects created when setting up feti mat context */ 7089566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)mat_ctx->lambda_local)); 709674ae819SStefano Zampini fetidppc_ctx->lambda_local = mat_ctx->lambda_local; 7109566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)mat_ctx->B_Ddelta)); 711674ae819SStefano Zampini fetidppc_ctx->B_Ddelta = mat_ctx->B_Ddelta; 71225d06dbeSstefano_zampini if (mat_ctx->deluxe_nonred) { 71325d06dbeSstefano_zampini PC pc,mpc; 71425d06dbeSstefano_zampini BDdelta_DN ctx; 7153ca39a21SBarry Smith MatSolverType solver; 71625d06dbeSstefano_zampini const char *prefix; 71725d06dbeSstefano_zampini 7189566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(mat_ctx->B_Ddelta,&ctx)); 7199566063dSJacob Faibussowitsch PetscCall(KSPSetType(ctx->kBD,KSPPREONLY)); 7209566063dSJacob Faibussowitsch PetscCall(KSPGetPC(ctx->kBD,&mpc)); 7219566063dSJacob Faibussowitsch PetscCall(KSPGetPC(pcbddc->ksp_D,&pc)); 7229566063dSJacob Faibussowitsch PetscCall(PCSetType(mpc,PCLU)); 7239566063dSJacob Faibussowitsch PetscCall(PCFactorGetMatSolverType(pc,(MatSolverType*)&solver)); 72425d06dbeSstefano_zampini if (solver) { 7259566063dSJacob Faibussowitsch PetscCall(PCFactorSetMatSolverType(mpc,solver)); 72625d06dbeSstefano_zampini } 7279566063dSJacob Faibussowitsch PetscCall(MatGetOptionsPrefix(fetimat,&prefix)); 7289566063dSJacob Faibussowitsch PetscCall(KSPSetOptionsPrefix(ctx->kBD,prefix)); 7299566063dSJacob Faibussowitsch PetscCall(KSPAppendOptionsPrefix(ctx->kBD,"bddelta_")); 7309566063dSJacob Faibussowitsch PetscCall(KSPSetFromOptions(ctx->kBD)); 73125d06dbeSstefano_zampini } 73225d06dbeSstefano_zampini 733e1214c54Sstefano_zampini if (mat_ctx->l2g_lambda_only) { 7349566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)mat_ctx->l2g_lambda_only)); 735e1214c54Sstefano_zampini fetidppc_ctx->l2g_lambda = mat_ctx->l2g_lambda_only; 736e1214c54Sstefano_zampini } else { 7379566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)mat_ctx->l2g_lambda)); 738674ae819SStefano Zampini fetidppc_ctx->l2g_lambda = mat_ctx->l2g_lambda; 739e1214c54Sstefano_zampini } 740f28b6018SStefano Zampini /* Dirichlet preconditioner */ 7419566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL,((PetscObject)fetimat)->prefix,"-pc_lumped",&lumped,NULL)); 742f28b6018SStefano Zampini if (!lumped) { 7439feb908dSstefano_zampini IS iV; 7449c2d02cdSstefano_zampini PetscBool discrete_harmonic = PETSC_FALSE; 7459c2d02cdSstefano_zampini 7469566063dSJacob Faibussowitsch PetscCall(PetscObjectQuery((PetscObject)fetidppc_ctx->pc,"__KSPFETIDP_iV",(PetscObject*)&iV)); 7479feb908dSstefano_zampini if (iV) { 7489566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL,((PetscObject)fetimat)->prefix,"-pc_discrete_harmonic",&discrete_harmonic,NULL)); 7499c2d02cdSstefano_zampini } 7509c2d02cdSstefano_zampini if (discrete_harmonic) { 7519c2d02cdSstefano_zampini KSP sksp; 7529c2d02cdSstefano_zampini PC pc; 753111315fdSstefano_zampini PCBDDCSubSchurs sub_schurs = pcbddc->sub_schurs; 7549c2d02cdSstefano_zampini Mat A_II,A_IB,A_BI; 755111315fdSstefano_zampini IS iP = NULL; 756111315fdSstefano_zampini PetscBool isshell,reuse = PETSC_FALSE; 7579c2d02cdSstefano_zampini KSPType ksptype; 7589c2d02cdSstefano_zampini const char *prefix; 7599c2d02cdSstefano_zampini 7609c2d02cdSstefano_zampini /* 7619c2d02cdSstefano_zampini We constructs a Schur complement for 7629c2d02cdSstefano_zampini 7639c2d02cdSstefano_zampini | A_II A_ID | 7649c2d02cdSstefano_zampini | A_DI A_DD | 7659c2d02cdSstefano_zampini 7669c2d02cdSstefano_zampini instead of 7679c2d02cdSstefano_zampini 7689c2d02cdSstefano_zampini | A_II B^t_II A_ID | 7699c2d02cdSstefano_zampini | B_II -C_II B_ID | 7709c2d02cdSstefano_zampini | A_DI B^t_ID A_DD | 7719c2d02cdSstefano_zampini 7729c2d02cdSstefano_zampini */ 773111315fdSstefano_zampini if (sub_schurs && sub_schurs->reuse_solver) { 7749566063dSJacob Faibussowitsch PetscCall(PetscObjectQuery((PetscObject)sub_schurs->A,"__KSPFETIDP_iP",(PetscObject*)&iP)); 775111315fdSstefano_zampini if (iP) reuse = PETSC_TRUE; 776111315fdSstefano_zampini } 777111315fdSstefano_zampini if (!reuse) { 778111315fdSstefano_zampini IS aB; 779111315fdSstefano_zampini PetscInt nb; 7809566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(pcis->is_B_local,&nb)); 7819566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PetscObjectComm((PetscObject)pcis->A_II),nb,0,1,&aB)); 7829566063dSJacob Faibussowitsch PetscCall(MatCreateSubMatrix(pcis->A_II,iV,iV,MAT_INITIAL_MATRIX,&A_II)); 7839566063dSJacob Faibussowitsch PetscCall(MatCreateSubMatrix(pcis->A_IB,iV,aB,MAT_INITIAL_MATRIX,&A_IB)); 7849566063dSJacob Faibussowitsch PetscCall(MatCreateSubMatrix(pcis->A_BI,aB,iV,MAT_INITIAL_MATRIX,&A_BI)); 7859566063dSJacob Faibussowitsch PetscCall(ISDestroy(&aB)); 786111315fdSstefano_zampini } else { 7879566063dSJacob Faibussowitsch PetscCall(MatCreateSubMatrix(sub_schurs->A,pcis->is_I_local,pcis->is_B_local,MAT_INITIAL_MATRIX,&A_IB)); 7889566063dSJacob Faibussowitsch PetscCall(MatCreateSubMatrix(sub_schurs->A,pcis->is_B_local,pcis->is_I_local,MAT_INITIAL_MATRIX,&A_BI)); 7899566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)pcis->A_II)); 790111315fdSstefano_zampini A_II = pcis->A_II; 791111315fdSstefano_zampini } 7929566063dSJacob Faibussowitsch PetscCall(MatCreateSchurComplement(A_II,A_II,A_IB,A_BI,pcis->A_BB,&fetidppc_ctx->S_j)); 7939c2d02cdSstefano_zampini 7949c2d02cdSstefano_zampini /* propagate settings of solver */ 7959566063dSJacob Faibussowitsch PetscCall(MatSchurComplementGetKSP(fetidppc_ctx->S_j,&sksp)); 7969566063dSJacob Faibussowitsch PetscCall(KSPGetType(pcis->ksp_D,&ksptype)); 7979566063dSJacob Faibussowitsch PetscCall(KSPSetType(sksp,ksptype)); 7989566063dSJacob Faibussowitsch PetscCall(KSPGetPC(pcis->ksp_D,&pc)); 7999566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)pc,PCSHELL,&isshell)); 8005334bea6Sstefano_zampini if (!isshell) { 8013ca39a21SBarry Smith MatSolverType solver; 8025334bea6Sstefano_zampini PCType pctype; 8035334bea6Sstefano_zampini 8049566063dSJacob Faibussowitsch PetscCall(PCGetType(pc,&pctype)); 8059566063dSJacob Faibussowitsch PetscCall(PCFactorGetMatSolverType(pc,(MatSolverType*)&solver)); 8069566063dSJacob Faibussowitsch PetscCall(KSPGetPC(sksp,&pc)); 8079566063dSJacob Faibussowitsch PetscCall(PCSetType(pc,pctype)); 8089c2d02cdSstefano_zampini if (solver) { 8099566063dSJacob Faibussowitsch PetscCall(PCFactorSetMatSolverType(pc,solver)); 8109c2d02cdSstefano_zampini } 8115334bea6Sstefano_zampini } else { 8129566063dSJacob Faibussowitsch PetscCall(KSPGetPC(sksp,&pc)); 8139566063dSJacob Faibussowitsch PetscCall(PCSetType(pc,PCLU)); 8145334bea6Sstefano_zampini } 8159566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A_II)); 8169566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A_IB)); 8179566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A_BI)); 8189566063dSJacob Faibussowitsch PetscCall(MatGetOptionsPrefix(fetimat,&prefix)); 8199566063dSJacob Faibussowitsch PetscCall(KSPSetOptionsPrefix(sksp,prefix)); 8209566063dSJacob Faibussowitsch PetscCall(KSPAppendOptionsPrefix(sksp,"harmonic_")); 8219566063dSJacob Faibussowitsch PetscCall(KSPSetFromOptions(sksp)); 822111315fdSstefano_zampini if (reuse) { 8239566063dSJacob Faibussowitsch PetscCall(KSPSetPC(sksp,sub_schurs->reuse_solver->interior_solver)); 8249566063dSJacob Faibussowitsch PetscCall(PetscObjectIncrementTabLevel((PetscObject)sub_schurs->reuse_solver->interior_solver,(PetscObject)sksp,0)); 825111315fdSstefano_zampini } 8269c2d02cdSstefano_zampini } else { /* default Dirichlet preconditioner is pde-harmonic */ 8279566063dSJacob Faibussowitsch PetscCall(MatCreateSchurComplement(pcis->A_II,pcis->A_II,pcis->A_IB,pcis->A_BI,pcis->A_BB,&fetidppc_ctx->S_j)); 8289566063dSJacob Faibussowitsch PetscCall(MatSchurComplementSetKSP(fetidppc_ctx->S_j,pcis->ksp_D)); 8299c2d02cdSstefano_zampini } 830f28b6018SStefano Zampini } else { 8319566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)pcis->A_BB)); 832f28b6018SStefano Zampini fetidppc_ctx->S_j = pcis->A_BB; 833f28b6018SStefano Zampini } 834af140850Sstefano_zampini /* saddle-point */ 835af140850Sstefano_zampini if (mat_ctx->xPg) { 8369566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)mat_ctx->xPg)); 837af140850Sstefano_zampini fetidppc_ctx->xPg = mat_ctx->xPg; 8389566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)mat_ctx->yPg)); 839af140850Sstefano_zampini fetidppc_ctx->yPg = mat_ctx->yPg; 840af140850Sstefano_zampini } 841674ae819SStefano Zampini PetscFunctionReturn(0); 842674ae819SStefano Zampini } 843674ae819SStefano Zampini 844ef425aeaSstefano_zampini PetscErrorCode FETIDPMatMult_Kernel(Mat fetimat, Vec x, Vec y, PetscBool trans) 845674ae819SStefano Zampini { 846674ae819SStefano Zampini FETIDPMat_ctx mat_ctx; 847617d11aeSStefano Zampini PC_BDDC *pcbddc; 848674ae819SStefano Zampini PC_IS *pcis; 849674ae819SStefano Zampini 850674ae819SStefano Zampini PetscFunctionBegin; 8519566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(fetimat,&mat_ctx)); 852674ae819SStefano Zampini pcis = (PC_IS*)mat_ctx->pc->data; 853617d11aeSStefano Zampini pcbddc = (PC_BDDC*)mat_ctx->pc->data; 854674ae819SStefano Zampini /* Application of B_delta^T */ 8559566063dSJacob Faibussowitsch PetscCall(VecSet(pcis->vec1_B,0.)); 8569566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(mat_ctx->l2g_lambda,x,mat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE)); 8579566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(mat_ctx->l2g_lambda,x,mat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE)); 8589566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(mat_ctx->B_delta,mat_ctx->lambda_local,pcis->vec1_B)); 859af140850Sstefano_zampini 860af140850Sstefano_zampini /* Add contribution from saddle point */ 861af140850Sstefano_zampini if (mat_ctx->l2g_p) { 8629566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(mat_ctx->l2g_p,x,mat_ctx->vP,INSERT_VALUES,SCATTER_REVERSE)); 8639566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(mat_ctx->l2g_p,x,mat_ctx->vP,INSERT_VALUES,SCATTER_REVERSE)); 864af140850Sstefano_zampini if (pcbddc->switch_static) { 865ef425aeaSstefano_zampini if (trans) { 8669566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(mat_ctx->B_BI,mat_ctx->vP,pcis->vec1_D)); 867ef425aeaSstefano_zampini } else { 8689566063dSJacob Faibussowitsch PetscCall(MatMult(mat_ctx->Bt_BI,mat_ctx->vP,pcis->vec1_D)); 869af140850Sstefano_zampini } 870ef425aeaSstefano_zampini } 871ef425aeaSstefano_zampini if (trans) { 8729566063dSJacob Faibussowitsch PetscCall(MatMultTransposeAdd(mat_ctx->B_BB,mat_ctx->vP,pcis->vec1_B,pcis->vec1_B)); 873ef425aeaSstefano_zampini } else { 8749566063dSJacob Faibussowitsch PetscCall(MatMultAdd(mat_ctx->Bt_BB,mat_ctx->vP,pcis->vec1_B,pcis->vec1_B)); 875ef425aeaSstefano_zampini } 876af140850Sstefano_zampini } else { 877af140850Sstefano_zampini if (pcbddc->switch_static) { 8789566063dSJacob Faibussowitsch PetscCall(VecSet(pcis->vec1_D,0.0)); 879af140850Sstefano_zampini } 880af140850Sstefano_zampini } 881af140850Sstefano_zampini /* Application of \widetilde{S}^-1 */ 8829566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(pcbddc->benign_p0,pcbddc->benign_n)); 8839566063dSJacob Faibussowitsch PetscCall(PCBDDCApplyInterfacePreconditioner(mat_ctx->pc,trans)); 8849566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(pcbddc->benign_p0,pcbddc->benign_n)); 8859566063dSJacob Faibussowitsch PetscCall(VecSet(y,0.0)); 886674ae819SStefano Zampini /* Application of B_delta */ 8879566063dSJacob Faibussowitsch PetscCall(MatMult(mat_ctx->B_delta,pcis->vec1_B,mat_ctx->lambda_local)); 888af140850Sstefano_zampini /* Contribution from boundary pressures */ 889af140850Sstefano_zampini if (mat_ctx->C) { 890af140850Sstefano_zampini const PetscScalar *lx; 891af140850Sstefano_zampini PetscScalar *ly; 892af140850Sstefano_zampini 89315579a77SStefano Zampini /* pressure ordered first in the local part of x and y */ 8949566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(x,&lx)); 8959566063dSJacob Faibussowitsch PetscCall(VecGetArray(y,&ly)); 8969566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(mat_ctx->xPg,lx)); 8979566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(mat_ctx->yPg,ly)); 898ef425aeaSstefano_zampini if (trans) { 8999566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(mat_ctx->C,mat_ctx->xPg,mat_ctx->yPg)); 900ef425aeaSstefano_zampini } else { 9019566063dSJacob Faibussowitsch PetscCall(MatMult(mat_ctx->C,mat_ctx->xPg,mat_ctx->yPg)); 902ef425aeaSstefano_zampini } 9039566063dSJacob Faibussowitsch PetscCall(VecResetArray(mat_ctx->xPg)); 9049566063dSJacob Faibussowitsch PetscCall(VecResetArray(mat_ctx->yPg)); 9059566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(x,&lx)); 9069566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(y,&ly)); 907af140850Sstefano_zampini } 908af140850Sstefano_zampini /* Add contribution from saddle point */ 909af140850Sstefano_zampini if (mat_ctx->l2g_p) { 910ef425aeaSstefano_zampini if (trans) { 9119566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(mat_ctx->Bt_BB,pcis->vec1_B,mat_ctx->vP)); 912ef425aeaSstefano_zampini } else { 9139566063dSJacob Faibussowitsch PetscCall(MatMult(mat_ctx->B_BB,pcis->vec1_B,mat_ctx->vP)); 914ef425aeaSstefano_zampini } 915af140850Sstefano_zampini if (pcbddc->switch_static) { 916ef425aeaSstefano_zampini if (trans) { 9179566063dSJacob Faibussowitsch PetscCall(MatMultTransposeAdd(mat_ctx->Bt_BI,pcis->vec1_D,mat_ctx->vP,mat_ctx->vP)); 918ef425aeaSstefano_zampini } else { 9199566063dSJacob Faibussowitsch PetscCall(MatMultAdd(mat_ctx->B_BI,pcis->vec1_D,mat_ctx->vP,mat_ctx->vP)); 920af140850Sstefano_zampini } 921ef425aeaSstefano_zampini } 9229566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(mat_ctx->l2g_p,mat_ctx->vP,y,ADD_VALUES,SCATTER_FORWARD)); 9239566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(mat_ctx->l2g_p,mat_ctx->vP,y,ADD_VALUES,SCATTER_FORWARD)); 924af140850Sstefano_zampini } 9259566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(mat_ctx->l2g_lambda,mat_ctx->lambda_local,y,ADD_VALUES,SCATTER_FORWARD)); 9269566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(mat_ctx->l2g_lambda,mat_ctx->lambda_local,y,ADD_VALUES,SCATTER_FORWARD)); 927674ae819SStefano Zampini PetscFunctionReturn(0); 928674ae819SStefano Zampini } 929674ae819SStefano Zampini 930ef425aeaSstefano_zampini PetscErrorCode FETIDPMatMult(Mat fetimat, Vec x, Vec y) 931edf7251bSStefano Zampini { 932edf7251bSStefano Zampini PetscFunctionBegin; 9339566063dSJacob Faibussowitsch PetscCall(FETIDPMatMult_Kernel(fetimat,x,y,PETSC_FALSE)); 934ef425aeaSstefano_zampini PetscFunctionReturn(0); 935ef425aeaSstefano_zampini } 936ef425aeaSstefano_zampini 937ef425aeaSstefano_zampini PetscErrorCode FETIDPMatMultTranspose(Mat fetimat, Vec x, Vec y) 938ef425aeaSstefano_zampini { 939ef425aeaSstefano_zampini PetscFunctionBegin; 9409566063dSJacob Faibussowitsch PetscCall(FETIDPMatMult_Kernel(fetimat,x,y,PETSC_TRUE)); 941edf7251bSStefano Zampini PetscFunctionReturn(0); 942edf7251bSStefano Zampini } 943edf7251bSStefano Zampini 94416597391Sstefano_zampini PetscErrorCode FETIDPPCApply_Kernel(PC fetipc, Vec x, Vec y, PetscBool trans) 945674ae819SStefano Zampini { 946674ae819SStefano Zampini FETIDPPC_ctx pc_ctx; 947674ae819SStefano Zampini PC_IS *pcis; 948674ae819SStefano Zampini 949674ae819SStefano Zampini PetscFunctionBegin; 9509566063dSJacob Faibussowitsch PetscCall(PCShellGetContext(fetipc,&pc_ctx)); 951674ae819SStefano Zampini pcis = (PC_IS*)pc_ctx->pc->data; 952674ae819SStefano Zampini /* Application of B_Ddelta^T */ 9539566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(pc_ctx->l2g_lambda,x,pc_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE)); 9549566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(pc_ctx->l2g_lambda,x,pc_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE)); 9559566063dSJacob Faibussowitsch PetscCall(VecSet(pcis->vec2_B,0.0)); 9569566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(pc_ctx->B_Ddelta,pc_ctx->lambda_local,pcis->vec2_B)); 957ed6c3d69SStefano Zampini /* Application of local Schur complement */ 95816597391Sstefano_zampini if (trans) { 9599566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(pc_ctx->S_j,pcis->vec2_B,pcis->vec1_B)); 96016597391Sstefano_zampini } else { 9619566063dSJacob Faibussowitsch PetscCall(MatMult(pc_ctx->S_j,pcis->vec2_B,pcis->vec1_B)); 96216597391Sstefano_zampini } 963edf7251bSStefano Zampini /* Application of B_Ddelta */ 9649566063dSJacob Faibussowitsch PetscCall(MatMult(pc_ctx->B_Ddelta,pcis->vec1_B,pc_ctx->lambda_local)); 9659566063dSJacob Faibussowitsch PetscCall(VecSet(y,0.0)); 9669566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(pc_ctx->l2g_lambda,pc_ctx->lambda_local,y,ADD_VALUES,SCATTER_FORWARD)); 9679566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(pc_ctx->l2g_lambda,pc_ctx->lambda_local,y,ADD_VALUES,SCATTER_FORWARD)); 968edf7251bSStefano Zampini PetscFunctionReturn(0); 969edf7251bSStefano Zampini } 970edf7251bSStefano Zampini 97116597391Sstefano_zampini PetscErrorCode FETIDPPCApply(PC pc, Vec x, Vec y) 972edf7251bSStefano Zampini { 973edf7251bSStefano Zampini PetscFunctionBegin; 9749566063dSJacob Faibussowitsch PetscCall(FETIDPPCApply_Kernel(pc,x,y,PETSC_FALSE)); 97516597391Sstefano_zampini PetscFunctionReturn(0); 97616597391Sstefano_zampini } 97716597391Sstefano_zampini 97816597391Sstefano_zampini PetscErrorCode FETIDPPCApplyTranspose(PC pc, Vec x, Vec y) 97916597391Sstefano_zampini { 98016597391Sstefano_zampini PetscFunctionBegin; 9819566063dSJacob Faibussowitsch PetscCall(FETIDPPCApply_Kernel(pc,x,y,PETSC_TRUE)); 982674ae819SStefano Zampini PetscFunctionReturn(0); 983674ae819SStefano Zampini } 984c45b8d2dSstefano_zampini 985c45b8d2dSstefano_zampini PetscErrorCode FETIDPPCView(PC pc, PetscViewer viewer) 986c45b8d2dSstefano_zampini { 987c45b8d2dSstefano_zampini FETIDPPC_ctx pc_ctx; 988c45b8d2dSstefano_zampini PetscBool iascii; 989c45b8d2dSstefano_zampini PetscViewer sviewer; 990c45b8d2dSstefano_zampini 991c45b8d2dSstefano_zampini PetscFunctionBegin; 9929566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii)); 993c45b8d2dSstefano_zampini if (iascii) { 994c45b8d2dSstefano_zampini PetscMPIInt rank; 99525d06dbeSstefano_zampini PetscBool isschur,isshell; 996c45b8d2dSstefano_zampini 9979566063dSJacob Faibussowitsch PetscCall(PCShellGetContext(pc,&pc_ctx)); 9989566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)pc),&rank)); 9999566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)pc_ctx->S_j,MATSCHURCOMPLEMENT,&isschur)); 1000c45b8d2dSstefano_zampini if (isschur) { 10019566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer," Dirichlet preconditioner (just from rank 0)\n")); 1002c45b8d2dSstefano_zampini } else { 10039566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer," Lumped preconditioner (just from rank 0)\n")); 1004c45b8d2dSstefano_zampini } 10059566063dSJacob Faibussowitsch PetscCall(PetscViewerGetSubViewer(viewer,PetscObjectComm((PetscObject)pc_ctx->S_j),&sviewer)); 1006dd400576SPatrick Sanan if (rank == 0) { 10079566063dSJacob Faibussowitsch PetscCall(PetscViewerPushFormat(sviewer,PETSC_VIEWER_ASCII_INFO)); 10089566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushTab(sviewer)); 10099566063dSJacob Faibussowitsch PetscCall(MatView(pc_ctx->S_j,sviewer)); 10109566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopTab(sviewer)); 10119566063dSJacob Faibussowitsch PetscCall(PetscViewerPopFormat(sviewer)); 1012c45b8d2dSstefano_zampini } 10139566063dSJacob Faibussowitsch PetscCall(PetscViewerRestoreSubViewer(viewer,PetscObjectComm((PetscObject)pc_ctx->S_j),&sviewer)); 10149566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)pc_ctx->B_Ddelta,MATSHELL,&isshell)); 101525d06dbeSstefano_zampini if (isshell) { 101625d06dbeSstefano_zampini BDdelta_DN ctx; 10179566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer," FETI-DP BDdelta: DB^t * (B D^-1 B^t)^-1 for deluxe scaling (just from rank 0)\n")); 10189566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(pc_ctx->B_Ddelta,&ctx)); 10199566063dSJacob Faibussowitsch PetscCall(PetscViewerGetSubViewer(viewer,PetscObjectComm((PetscObject)pc_ctx->S_j),&sviewer)); 1020dd400576SPatrick Sanan if (rank == 0) { 102115579a77SStefano Zampini PetscInt tl; 102215579a77SStefano Zampini 10239566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIGetTab(sviewer,&tl)); 10249566063dSJacob Faibussowitsch PetscCall(PetscObjectSetTabLevel((PetscObject)ctx->kBD,tl)); 10259566063dSJacob Faibussowitsch PetscCall(KSPView(ctx->kBD,sviewer)); 10269566063dSJacob Faibussowitsch PetscCall(PetscViewerPushFormat(sviewer,PETSC_VIEWER_ASCII_INFO)); 10279566063dSJacob Faibussowitsch PetscCall(MatView(ctx->BD,sviewer)); 10289566063dSJacob Faibussowitsch PetscCall(PetscViewerPopFormat(sviewer)); 102925d06dbeSstefano_zampini } 10309566063dSJacob Faibussowitsch PetscCall(PetscViewerRestoreSubViewer(viewer,PetscObjectComm((PetscObject)pc_ctx->S_j),&sviewer)); 1031e5afcf28SBarry Smith } 10329566063dSJacob Faibussowitsch PetscCall(PetscViewerFlush(viewer)); 1033c45b8d2dSstefano_zampini } 1034c45b8d2dSstefano_zampini PetscFunctionReturn(0); 1035c45b8d2dSstefano_zampini } 1036