1ab5c6b0cSJed Brown #include <../src/ksp/pc/impls/bddc/bddc.h> 2ab5c6b0cSJed Brown #include <../src/ksp/pc/impls/bddc/bddcprivate.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 PetscErrorCode ierr; 925d06dbeSstefano_zampini 1025d06dbeSstefano_zampini PetscFunctionBegin; 1125d06dbeSstefano_zampini ierr = MatShellGetContext(A,(void**)&ctx);CHKERRQ(ierr); 1225d06dbeSstefano_zampini ierr = MatMultTranspose(ctx->BD,x,ctx->work);CHKERRQ(ierr); 1325d06dbeSstefano_zampini ierr = KSPSolveTranspose(ctx->kBD,ctx->work,y);CHKERRQ(ierr); 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 PetscErrorCode ierr; 2125d06dbeSstefano_zampini 2225d06dbeSstefano_zampini PetscFunctionBegin; 2325d06dbeSstefano_zampini ierr = MatShellGetContext(A,(void**)&ctx);CHKERRQ(ierr); 2425d06dbeSstefano_zampini ierr = KSPSolve(ctx->kBD,x,ctx->work);CHKERRQ(ierr); 2525d06dbeSstefano_zampini ierr = MatMult(ctx->BD,ctx->work,y);CHKERRQ(ierr); 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 PetscErrorCode ierr; 3325d06dbeSstefano_zampini 3425d06dbeSstefano_zampini PetscFunctionBegin; 3525d06dbeSstefano_zampini ierr = MatShellGetContext(A,(void**)&ctx);CHKERRQ(ierr); 3625d06dbeSstefano_zampini ierr = MatDestroy(&ctx->BD);CHKERRQ(ierr); 3725d06dbeSstefano_zampini ierr = KSPDestroy(&ctx->kBD);CHKERRQ(ierr); 3825d06dbeSstefano_zampini ierr = VecDestroy(&ctx->work);CHKERRQ(ierr); 3925d06dbeSstefano_zampini ierr = PetscFree(ctx);CHKERRQ(ierr); 4025d06dbeSstefano_zampini PetscFunctionReturn(0); 4125d06dbeSstefano_zampini } 4225d06dbeSstefano_zampini 43674ae819SStefano Zampini 44674ae819SStefano Zampini PetscErrorCode PCBDDCCreateFETIDPMatContext(PC pc, FETIDPMat_ctx *fetidpmat_ctx) 45674ae819SStefano Zampini { 46674ae819SStefano Zampini FETIDPMat_ctx newctx; 47674ae819SStefano Zampini PetscErrorCode ierr; 48674ae819SStefano Zampini 49674ae819SStefano Zampini PetscFunctionBegin; 50854ce69bSBarry Smith ierr = PetscNew(&newctx);CHKERRQ(ierr); 51674ae819SStefano Zampini /* increase the reference count for BDDC preconditioner */ 52674ae819SStefano Zampini ierr = PetscObjectReference((PetscObject)pc);CHKERRQ(ierr); 53674ae819SStefano Zampini newctx->pc = pc; 54674ae819SStefano Zampini *fetidpmat_ctx = newctx; 55674ae819SStefano Zampini PetscFunctionReturn(0); 56674ae819SStefano Zampini } 57674ae819SStefano Zampini 58674ae819SStefano Zampini PetscErrorCode PCBDDCCreateFETIDPPCContext(PC pc, FETIDPPC_ctx *fetidppc_ctx) 59674ae819SStefano Zampini { 60674ae819SStefano Zampini FETIDPPC_ctx newctx; 61674ae819SStefano Zampini PetscErrorCode ierr; 62674ae819SStefano Zampini 63674ae819SStefano Zampini PetscFunctionBegin; 64854ce69bSBarry Smith ierr = PetscNew(&newctx);CHKERRQ(ierr); 65674ae819SStefano Zampini /* increase the reference count for BDDC preconditioner */ 66674ae819SStefano Zampini ierr = PetscObjectReference((PetscObject)pc);CHKERRQ(ierr); 67674ae819SStefano Zampini newctx->pc = pc; 68674ae819SStefano Zampini *fetidppc_ctx = newctx; 69674ae819SStefano Zampini PetscFunctionReturn(0); 70674ae819SStefano Zampini } 71674ae819SStefano Zampini 72674ae819SStefano Zampini PetscErrorCode PCBDDCDestroyFETIDPMat(Mat A) 73674ae819SStefano Zampini { 74674ae819SStefano Zampini FETIDPMat_ctx mat_ctx; 75674ae819SStefano Zampini PetscErrorCode ierr; 76674ae819SStefano Zampini 77674ae819SStefano Zampini PetscFunctionBegin; 78674ae819SStefano Zampini ierr = MatShellGetContext(A,(void**)&mat_ctx);CHKERRQ(ierr); 79674ae819SStefano Zampini ierr = VecDestroy(&mat_ctx->lambda_local);CHKERRQ(ierr); 80674ae819SStefano Zampini ierr = VecDestroy(&mat_ctx->temp_solution_D);CHKERRQ(ierr); 81674ae819SStefano Zampini ierr = VecDestroy(&mat_ctx->temp_solution_B);CHKERRQ(ierr); 82674ae819SStefano Zampini ierr = MatDestroy(&mat_ctx->B_delta);CHKERRQ(ierr); 83674ae819SStefano Zampini ierr = MatDestroy(&mat_ctx->B_Ddelta);CHKERRQ(ierr); 84457d4a33Sstefano_zampini ierr = MatDestroy(&mat_ctx->B_BB);CHKERRQ(ierr); 85457d4a33Sstefano_zampini ierr = MatDestroy(&mat_ctx->B_BI);CHKERRQ(ierr); 86457d4a33Sstefano_zampini ierr = MatDestroy(&mat_ctx->Bt_BB);CHKERRQ(ierr); 87457d4a33Sstefano_zampini ierr = MatDestroy(&mat_ctx->Bt_BI);CHKERRQ(ierr); 88457d4a33Sstefano_zampini ierr = MatDestroy(&mat_ctx->C);CHKERRQ(ierr); 896cc1294bSstefano_zampini ierr = VecDestroy(&mat_ctx->rhs_flip);CHKERRQ(ierr); 90457d4a33Sstefano_zampini ierr = VecDestroy(&mat_ctx->vP);CHKERRQ(ierr); 91457d4a33Sstefano_zampini ierr = VecDestroy(&mat_ctx->xPg);CHKERRQ(ierr); 92457d4a33Sstefano_zampini ierr = VecDestroy(&mat_ctx->yPg);CHKERRQ(ierr); 93674ae819SStefano Zampini ierr = VecScatterDestroy(&mat_ctx->l2g_lambda);CHKERRQ(ierr); 94457d4a33Sstefano_zampini ierr = VecScatterDestroy(&mat_ctx->l2g_p);CHKERRQ(ierr); 959cc7774eSstefano_zampini ierr = VecScatterDestroy(&mat_ctx->g2g_p);CHKERRQ(ierr); 96674ae819SStefano Zampini ierr = PCDestroy(&mat_ctx->pc);CHKERRQ(ierr); /* decrease PCBDDC reference count */ 97674ae819SStefano Zampini ierr = PetscFree(mat_ctx);CHKERRQ(ierr); 98674ae819SStefano Zampini PetscFunctionReturn(0); 99674ae819SStefano Zampini } 100674ae819SStefano Zampini 101674ae819SStefano Zampini PetscErrorCode PCBDDCDestroyFETIDPPC(PC pc) 102674ae819SStefano Zampini { 103674ae819SStefano Zampini FETIDPPC_ctx pc_ctx; 104674ae819SStefano Zampini PetscErrorCode ierr; 105674ae819SStefano Zampini 106674ae819SStefano Zampini PetscFunctionBegin; 107674ae819SStefano Zampini ierr = PCShellGetContext(pc,(void**)&pc_ctx);CHKERRQ(ierr); 108674ae819SStefano Zampini ierr = VecDestroy(&pc_ctx->lambda_local);CHKERRQ(ierr); 109674ae819SStefano Zampini ierr = MatDestroy(&pc_ctx->B_Ddelta);CHKERRQ(ierr); 110674ae819SStefano Zampini ierr = VecScatterDestroy(&pc_ctx->l2g_lambda);CHKERRQ(ierr); 11168668abeSStefano Zampini ierr = MatDestroy(&pc_ctx->S_j);CHKERRQ(ierr); 112674ae819SStefano Zampini ierr = PCDestroy(&pc_ctx->pc);CHKERRQ(ierr); /* decrease PCBDDC reference count */ 113457d4a33Sstefano_zampini ierr = KSPDestroy(&pc_ctx->kP);CHKERRQ(ierr); 114457d4a33Sstefano_zampini ierr = VecDestroy(&pc_ctx->xPg);CHKERRQ(ierr); 115457d4a33Sstefano_zampini ierr = VecDestroy(&pc_ctx->yPg);CHKERRQ(ierr); 116674ae819SStefano Zampini ierr = PetscFree(pc_ctx);CHKERRQ(ierr); 117674ae819SStefano Zampini PetscFunctionReturn(0); 118674ae819SStefano Zampini } 119674ae819SStefano Zampini 120674ae819SStefano Zampini PetscErrorCode PCBDDCSetupFETIDPMatContext(FETIDPMat_ctx fetidpmat_ctx ) 121674ae819SStefano Zampini { 122674ae819SStefano Zampini PetscErrorCode ierr; 123674ae819SStefano Zampini PC_IS *pcis=(PC_IS*)fetidpmat_ctx->pc->data; 124674ae819SStefano Zampini PC_BDDC *pcbddc=(PC_BDDC*)fetidpmat_ctx->pc->data; 125674ae819SStefano Zampini PCBDDCGraph mat_graph=pcbddc->mat_graph; 126be8bcbd0Sstefano_zampini #if defined(PETSC_USE_DEBUG) 127674ae819SStefano Zampini Mat_IS *matis = (Mat_IS*)fetidpmat_ctx->pc->pmat->data; 128be8bcbd0Sstefano_zampini #endif 129674ae819SStefano Zampini MPI_Comm comm; 13025d06dbeSstefano_zampini Mat ScalingMat,BD1,BD2; 131457d4a33Sstefano_zampini Vec fetidp_global; 132674ae819SStefano Zampini IS IS_l2g_lambda; 133a1c0d0daSstefano_zampini IS subset,subset_mult,subset_n,isvert; 134674ae819SStefano Zampini PetscBool skip_node,fully_redundant; 135674ae819SStefano Zampini PetscInt i,j,k,s,n_boundary_dofs,n_global_lambda,n_vertices,partial_sum; 136dc456d91SStefano Zampini PetscInt cum,n_local_lambda,n_lambda_for_dof,dual_size,n_neg_values,n_pos_values; 13776ec1555SBarry Smith PetscMPIInt rank,size,buf_size,neigh; 138674ae819SStefano Zampini PetscScalar scalar_value; 139a1c0d0daSstefano_zampini const PetscInt *vertex_indices; 140dc456d91SStefano Zampini PetscInt *dual_dofs_boundary_indices,*aux_local_numbering_1; 141dc456d91SStefano Zampini const PetscInt *aux_global_numbering; 142674ae819SStefano Zampini PetscInt *aux_sums,*cols_B_delta,*l2g_indices; 143674ae819SStefano Zampini PetscScalar *array,*scaling_factors,*vals_B_delta; 14425d06dbeSstefano_zampini PetscScalar **all_factors; 145674ae819SStefano Zampini PetscInt *aux_local_numbering_2; 146457d4a33Sstefano_zampini PetscLayout llay; 147a1c0d0daSstefano_zampini 148457d4a33Sstefano_zampini /* saddlepoint */ 149457d4a33Sstefano_zampini ISLocalToGlobalMapping l2gmap_p; 150457d4a33Sstefano_zampini PetscLayout play; 151457d4a33Sstefano_zampini IS gP,pP; 152457d4a33Sstefano_zampini PetscInt nPl,nPg,nPgl; 153674ae819SStefano Zampini 154674ae819SStefano Zampini PetscFunctionBegin; 155674ae819SStefano Zampini ierr = PetscObjectGetComm((PetscObject)(fetidpmat_ctx->pc),&comm);CHKERRQ(ierr); 156674ae819SStefano Zampini ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 15776ec1555SBarry Smith ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 158674ae819SStefano Zampini 159457d4a33Sstefano_zampini /* saddlepoint */ 160457d4a33Sstefano_zampini nPl = 0; 161457d4a33Sstefano_zampini nPg = 0; 162457d4a33Sstefano_zampini nPgl = 0; 163457d4a33Sstefano_zampini gP = NULL; 164457d4a33Sstefano_zampini pP = NULL; 165457d4a33Sstefano_zampini l2gmap_p = NULL; 166457d4a33Sstefano_zampini play = NULL; 167457d4a33Sstefano_zampini ierr = PetscObjectQuery((PetscObject)fetidpmat_ctx->pc,"__KSPFETIDP_pP",(PetscObject*)&pP);CHKERRQ(ierr); 168022d8d2bSstefano_zampini if (pP) { /* saddle point */ 169457d4a33Sstefano_zampini /* subdomain pressures in global numbering */ 170457d4a33Sstefano_zampini ierr = PetscObjectQuery((PetscObject)fetidpmat_ctx->pc,"__KSPFETIDP_gP",(PetscObject*)&gP);CHKERRQ(ierr); 171457d4a33Sstefano_zampini if (!gP) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"gP not present"); 172457d4a33Sstefano_zampini ierr = ISGetLocalSize(gP,&nPl);CHKERRQ(ierr); 173457d4a33Sstefano_zampini ierr = VecCreate(PETSC_COMM_SELF,&fetidpmat_ctx->vP);CHKERRQ(ierr); 174457d4a33Sstefano_zampini ierr = VecSetSizes(fetidpmat_ctx->vP,nPl,nPl);CHKERRQ(ierr); 175457d4a33Sstefano_zampini ierr = VecSetType(fetidpmat_ctx->vP,VECSTANDARD);CHKERRQ(ierr); 176457d4a33Sstefano_zampini ierr = VecSetUp(fetidpmat_ctx->vP);CHKERRQ(ierr); 177457d4a33Sstefano_zampini 178457d4a33Sstefano_zampini /* interface pressure matrix */ 179457d4a33Sstefano_zampini ierr = PetscObjectQuery((PetscObject)fetidpmat_ctx->pc,"__KSPFETIDP_C",(PetscObject*)&fetidpmat_ctx->C);CHKERRQ(ierr); 180457d4a33Sstefano_zampini if (!fetidpmat_ctx->C) { /* null pressure block, compute layout and global numbering for interface pressures */ 181457d4a33Sstefano_zampini IS Pg; 182457d4a33Sstefano_zampini 183457d4a33Sstefano_zampini ierr = ISRenumber(gP,NULL,&nPg,&Pg);CHKERRQ(ierr); 184457d4a33Sstefano_zampini ierr = ISLocalToGlobalMappingCreateIS(Pg,&l2gmap_p);CHKERRQ(ierr); 185457d4a33Sstefano_zampini ierr = ISDestroy(&Pg);CHKERRQ(ierr); 186457d4a33Sstefano_zampini ierr = PetscLayoutCreate(comm,&play);CHKERRQ(ierr); 187457d4a33Sstefano_zampini ierr = PetscLayoutSetBlockSize(play,1);CHKERRQ(ierr); 188457d4a33Sstefano_zampini ierr = PetscLayoutSetSize(play,nPg);CHKERRQ(ierr); 189457d4a33Sstefano_zampini ierr = ISGetLocalSize(pP,&nPgl);CHKERRQ(ierr); 190457d4a33Sstefano_zampini ierr = PetscLayoutSetLocalSize(play,nPgl);CHKERRQ(ierr); 191457d4a33Sstefano_zampini ierr = PetscLayoutSetUp(play);CHKERRQ(ierr); 192457d4a33Sstefano_zampini } else { 193457d4a33Sstefano_zampini ierr = PetscObjectReference((PetscObject)fetidpmat_ctx->C);CHKERRQ(ierr); 194457d4a33Sstefano_zampini ierr = MatGetLocalToGlobalMapping(fetidpmat_ctx->C,&l2gmap_p,NULL);CHKERRQ(ierr); 195457d4a33Sstefano_zampini ierr = PetscObjectReference((PetscObject)l2gmap_p);CHKERRQ(ierr); 196457d4a33Sstefano_zampini ierr = MatGetSize(fetidpmat_ctx->C,&nPg,NULL);CHKERRQ(ierr); 197457d4a33Sstefano_zampini ierr = MatGetLocalSize(fetidpmat_ctx->C,NULL,&nPgl);CHKERRQ(ierr); 198457d4a33Sstefano_zampini ierr = MatGetLayouts(fetidpmat_ctx->C,NULL,&llay);CHKERRQ(ierr); 199457d4a33Sstefano_zampini ierr = PetscLayoutReference(llay,&play);CHKERRQ(ierr); 200457d4a33Sstefano_zampini } 201457d4a33Sstefano_zampini ierr = VecCreateMPIWithArray(comm,1,nPgl,nPg,NULL,&fetidpmat_ctx->xPg);CHKERRQ(ierr); 202457d4a33Sstefano_zampini ierr = VecCreateMPIWithArray(comm,1,nPgl,nPg,NULL,&fetidpmat_ctx->yPg);CHKERRQ(ierr); 203457d4a33Sstefano_zampini 204457d4a33Sstefano_zampini /* import matrices for interface pressures coupling */ 205457d4a33Sstefano_zampini ierr = PetscObjectQuery((PetscObject)fetidpmat_ctx->pc,"__KSPFETIDP_B_BI",(PetscObject*)&fetidpmat_ctx->B_BI);CHKERRQ(ierr); 206457d4a33Sstefano_zampini if (!fetidpmat_ctx->B_BI) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"B_BI not present"); 207457d4a33Sstefano_zampini ierr = PetscObjectReference((PetscObject)fetidpmat_ctx->B_BI);CHKERRQ(ierr); 208457d4a33Sstefano_zampini 209457d4a33Sstefano_zampini ierr = PetscObjectQuery((PetscObject)fetidpmat_ctx->pc,"__KSPFETIDP_B_BB",(PetscObject*)&fetidpmat_ctx->B_BB);CHKERRQ(ierr); 210457d4a33Sstefano_zampini if (!fetidpmat_ctx->B_BB) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"B_BB not present"); 211457d4a33Sstefano_zampini ierr = PetscObjectReference((PetscObject)fetidpmat_ctx->B_BB);CHKERRQ(ierr); 212457d4a33Sstefano_zampini 213457d4a33Sstefano_zampini ierr = PetscObjectQuery((PetscObject)fetidpmat_ctx->pc,"__KSPFETIDP_Bt_BI",(PetscObject*)&fetidpmat_ctx->Bt_BI);CHKERRQ(ierr); 214457d4a33Sstefano_zampini if (!fetidpmat_ctx->Bt_BI) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Bt_BI not present"); 215457d4a33Sstefano_zampini ierr = PetscObjectReference((PetscObject)fetidpmat_ctx->Bt_BI);CHKERRQ(ierr); 216457d4a33Sstefano_zampini 217457d4a33Sstefano_zampini ierr = PetscObjectQuery((PetscObject)fetidpmat_ctx->pc,"__KSPFETIDP_Bt_BB",(PetscObject*)&fetidpmat_ctx->Bt_BB);CHKERRQ(ierr); 218457d4a33Sstefano_zampini if (!fetidpmat_ctx->Bt_BB) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Bt_BB not present"); 219457d4a33Sstefano_zampini ierr = PetscObjectReference((PetscObject)fetidpmat_ctx->Bt_BB);CHKERRQ(ierr); 2206cc1294bSstefano_zampini 2216cc1294bSstefano_zampini ierr = PetscObjectQuery((PetscObject)fetidpmat_ctx->pc,"__KSPFETIDP_flip" ,(PetscObject*)&fetidpmat_ctx->rhs_flip);CHKERRQ(ierr); 2226cc1294bSstefano_zampini if (fetidpmat_ctx->rhs_flip) { 2236cc1294bSstefano_zampini ierr = PetscObjectReference((PetscObject)fetidpmat_ctx->rhs_flip);CHKERRQ(ierr); 2246cc1294bSstefano_zampini } 225457d4a33Sstefano_zampini } 226457d4a33Sstefano_zampini 227674ae819SStefano Zampini /* Default type of lagrange multipliers is non-redundant */ 228329cd39dSStefano Zampini fully_redundant = fetidpmat_ctx->fully_redundant; 229674ae819SStefano Zampini 230674ae819SStefano Zampini /* Evaluate local and global number of lagrange multipliers */ 231674ae819SStefano Zampini ierr = VecSet(pcis->vec1_N,0.0);CHKERRQ(ierr); 232674ae819SStefano Zampini n_local_lambda = 0; 233674ae819SStefano Zampini partial_sum = 0; 234674ae819SStefano Zampini n_boundary_dofs = 0; 235674ae819SStefano Zampini s = 0; 236a1c0d0daSstefano_zampini 237674ae819SStefano Zampini /* Get Vertices used to define the BDDC */ 238a1c0d0daSstefano_zampini ierr = PCBDDCGraphGetCandidatesIS(pcbddc->mat_graph,NULL,NULL,NULL,NULL,&isvert);CHKERRQ(ierr); 239a1c0d0daSstefano_zampini ierr = ISGetLocalSize(isvert,&n_vertices);CHKERRQ(ierr); 240a1c0d0daSstefano_zampini ierr = ISGetIndices(isvert,&vertex_indices);CHKERRQ(ierr); 241a1c0d0daSstefano_zampini 242674ae819SStefano Zampini dual_size = pcis->n_B-n_vertices; 243785e854fSJed Brown ierr = PetscMalloc1(dual_size,&dual_dofs_boundary_indices);CHKERRQ(ierr); 244785e854fSJed Brown ierr = PetscMalloc1(dual_size,&aux_local_numbering_1);CHKERRQ(ierr); 245785e854fSJed Brown ierr = PetscMalloc1(dual_size,&aux_local_numbering_2);CHKERRQ(ierr); 246674ae819SStefano Zampini 247674ae819SStefano Zampini ierr = VecGetArray(pcis->vec1_N,&array);CHKERRQ(ierr); 248674ae819SStefano Zampini for (i=0;i<pcis->n;i++){ 249674ae819SStefano Zampini j = mat_graph->count[i]; /* RECALL: mat_graph->count[i] does not count myself */ 2502d456bbaSstefano_zampini if (j > 0) n_boundary_dofs++; 251674ae819SStefano Zampini skip_node = PETSC_FALSE; 252674ae819SStefano Zampini if (s < n_vertices && vertex_indices[s] == i) { /* it works for a sorted set of vertices */ 253674ae819SStefano Zampini skip_node = PETSC_TRUE; 254674ae819SStefano Zampini s++; 255674ae819SStefano Zampini } 2562d456bbaSstefano_zampini if (j < 1) skip_node = PETSC_TRUE; 2572d456bbaSstefano_zampini if (mat_graph->special_dof[i] == PCBDDCGRAPH_DIRICHLET_MARK) skip_node = PETSC_TRUE; 258674ae819SStefano Zampini if (!skip_node) { 259674ae819SStefano Zampini if (fully_redundant) { 260674ae819SStefano Zampini /* fully redundant set of lagrange multipliers */ 261674ae819SStefano Zampini n_lambda_for_dof = (j*(j+1))/2; 262674ae819SStefano Zampini } else { 263674ae819SStefano Zampini n_lambda_for_dof = j; 264674ae819SStefano Zampini } 265674ae819SStefano Zampini n_local_lambda += j; 266674ae819SStefano Zampini /* needed to evaluate global number of lagrange multipliers */ 267674ae819SStefano Zampini array[i]=(1.0*n_lambda_for_dof)/(j+1.0); /* already scaled for the next global sum */ 268674ae819SStefano Zampini /* store some data needed */ 269674ae819SStefano Zampini dual_dofs_boundary_indices[partial_sum] = n_boundary_dofs-1; 270674ae819SStefano Zampini aux_local_numbering_1[partial_sum] = i; 271674ae819SStefano Zampini aux_local_numbering_2[partial_sum] = n_lambda_for_dof; 272674ae819SStefano Zampini partial_sum++; 273674ae819SStefano Zampini } 274674ae819SStefano Zampini } 275674ae819SStefano Zampini ierr = VecRestoreArray(pcis->vec1_N,&array);CHKERRQ(ierr); 276a1c0d0daSstefano_zampini ierr = ISRestoreIndices(isvert,&vertex_indices);CHKERRQ(ierr); 277a1c0d0daSstefano_zampini ierr = PCBDDCGraphRestoreCandidatesIS(pcbddc->mat_graph,NULL,NULL,NULL,NULL,&isvert);CHKERRQ(ierr); 2782d456bbaSstefano_zampini dual_size = partial_sum; 279674ae819SStefano Zampini 280674ae819SStefano Zampini /* compute global ordering of lagrange multipliers and associate l2g map */ 281dc456d91SStefano Zampini ierr = ISCreateGeneral(comm,partial_sum,aux_local_numbering_1,PETSC_COPY_VALUES,&subset_n);CHKERRQ(ierr); 2823bbff08aSStefano Zampini ierr = ISLocalToGlobalMappingApplyIS(pcis->mapping,subset_n,&subset);CHKERRQ(ierr); 283dc456d91SStefano Zampini ierr = ISDestroy(&subset_n);CHKERRQ(ierr); 284dc456d91SStefano Zampini ierr = ISCreateGeneral(comm,partial_sum,aux_local_numbering_2,PETSC_OWN_POINTER,&subset_mult);CHKERRQ(ierr); 2853d996552SStefano Zampini ierr = ISRenumber(subset,subset_mult,&fetidpmat_ctx->n_lambda,&subset_n);CHKERRQ(ierr); 286dc456d91SStefano Zampini ierr = ISDestroy(&subset);CHKERRQ(ierr); 2873d996552SStefano Zampini 2884fe826edSStefano Zampini #if defined(PETSC_USE_DEBUG) 2894fe826edSStefano Zampini ierr = VecSet(pcis->vec1_global,0.0);CHKERRQ(ierr); 2904fe826edSStefano Zampini ierr = VecScatterBegin(matis->rctx,pcis->vec1_N,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 2914fe826edSStefano Zampini ierr = VecScatterEnd(matis->rctx,pcis->vec1_N,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 2924fe826edSStefano Zampini ierr = VecSum(pcis->vec1_global,&scalar_value);CHKERRQ(ierr); 2934fe826edSStefano Zampini i = (PetscInt)PetscRealPart(scalar_value); 2946c4ed002SBarry 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); 2954fe826edSStefano Zampini #endif 296674ae819SStefano Zampini 297674ae819SStefano Zampini /* init data for scaling factors exchange */ 29825d06dbeSstefano_zampini if (!pcbddc->use_deluxe_scaling) { 29925d06dbeSstefano_zampini PetscInt *ptrs_buffer,neigh_position; 30025d06dbeSstefano_zampini PetscScalar *send_buffer,*recv_buffer; 30125d06dbeSstefano_zampini MPI_Request *send_reqs,*recv_reqs; 30225d06dbeSstefano_zampini 303674ae819SStefano Zampini partial_sum = 0; 304785e854fSJed Brown ierr = PetscMalloc1(pcis->n_neigh,&ptrs_buffer);CHKERRQ(ierr); 3054b2aedd3SStefano Zampini ierr = PetscMalloc1(PetscMax(pcis->n_neigh-1,0),&send_reqs);CHKERRQ(ierr); 3064b2aedd3SStefano Zampini ierr = PetscMalloc1(PetscMax(pcis->n_neigh-1,0),&recv_reqs);CHKERRQ(ierr); 30719c16490Sstefano_zampini ierr = PetscMalloc1(pcis->n+1,&all_factors);CHKERRQ(ierr); 3084b2aedd3SStefano Zampini if (pcis->n_neigh > 0) ptrs_buffer[0]=0; 309674ae819SStefano Zampini for (i=1;i<pcis->n_neigh;i++) { 310674ae819SStefano Zampini partial_sum += pcis->n_shared[i]; 311674ae819SStefano Zampini ptrs_buffer[i] = ptrs_buffer[i-1]+pcis->n_shared[i]; 312674ae819SStefano Zampini } 313785e854fSJed Brown ierr = PetscMalloc1(partial_sum,&send_buffer);CHKERRQ(ierr); 314785e854fSJed Brown ierr = PetscMalloc1(partial_sum,&recv_buffer);CHKERRQ(ierr); 315785e854fSJed Brown ierr = PetscMalloc1(partial_sum,&all_factors[0]);CHKERRQ(ierr); 316674ae819SStefano Zampini for (i=0;i<pcis->n-1;i++) { 317674ae819SStefano Zampini j = mat_graph->count[i]; 318674ae819SStefano Zampini all_factors[i+1]=all_factors[i]+j; 319674ae819SStefano Zampini } 32025d06dbeSstefano_zampini 321674ae819SStefano Zampini /* scatter B scaling to N vec */ 322674ae819SStefano Zampini ierr = VecScatterBegin(pcis->N_to_B,pcis->D,pcis->vec1_N,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 323674ae819SStefano Zampini ierr = VecScatterEnd(pcis->N_to_B,pcis->D,pcis->vec1_N,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 324674ae819SStefano Zampini /* communications */ 3252b095fd8SStefano Zampini ierr = VecGetArrayRead(pcis->vec1_N,(const PetscScalar**)&array);CHKERRQ(ierr); 326674ae819SStefano Zampini for (i=1;i<pcis->n_neigh;i++) { 327674ae819SStefano Zampini for (j=0;j<pcis->n_shared[i];j++) { 328674ae819SStefano Zampini send_buffer[ptrs_buffer[i-1]+j]=array[pcis->shared[i][j]]; 329674ae819SStefano Zampini } 330674ae819SStefano Zampini ierr = PetscMPIIntCast(ptrs_buffer[i]-ptrs_buffer[i-1],&buf_size);CHKERRQ(ierr); 331674ae819SStefano Zampini ierr = PetscMPIIntCast(pcis->neigh[i],&neigh);CHKERRQ(ierr); 332674ae819SStefano Zampini ierr = MPI_Isend(&send_buffer[ptrs_buffer[i-1]],buf_size,MPIU_SCALAR,neigh,0,comm,&send_reqs[i-1]);CHKERRQ(ierr); 333674ae819SStefano Zampini ierr = MPI_Irecv(&recv_buffer[ptrs_buffer[i-1]],buf_size,MPIU_SCALAR,neigh,0,comm,&recv_reqs[i-1]);CHKERRQ(ierr); 334674ae819SStefano Zampini } 3352b095fd8SStefano Zampini ierr = VecRestoreArrayRead(pcis->vec1_N,(const PetscScalar**)&array);CHKERRQ(ierr); 3364b2aedd3SStefano Zampini if (pcis->n_neigh > 0) { 3374b2aedd3SStefano Zampini ierr = MPI_Waitall(pcis->n_neigh-1,recv_reqs,MPI_STATUSES_IGNORE);CHKERRQ(ierr); 3384b2aedd3SStefano Zampini } 339674ae819SStefano Zampini /* put values in correct places */ 340674ae819SStefano Zampini for (i=1;i<pcis->n_neigh;i++) { 341674ae819SStefano Zampini for (j=0;j<pcis->n_shared[i];j++) { 342674ae819SStefano Zampini k = pcis->shared[i][j]; 343674ae819SStefano Zampini neigh_position = 0; 344674ae819SStefano Zampini while(mat_graph->neighbours_set[k][neigh_position] != pcis->neigh[i]) {neigh_position++;} 345674ae819SStefano Zampini all_factors[k][neigh_position]=recv_buffer[ptrs_buffer[i-1]+j]; 346674ae819SStefano Zampini } 347674ae819SStefano Zampini } 3484b2aedd3SStefano Zampini if (pcis->n_neigh > 0) { 3494b2aedd3SStefano Zampini ierr = MPI_Waitall(pcis->n_neigh-1,send_reqs,MPI_STATUSES_IGNORE);CHKERRQ(ierr); 3504b2aedd3SStefano Zampini } 351674ae819SStefano Zampini ierr = PetscFree(send_reqs);CHKERRQ(ierr); 352674ae819SStefano Zampini ierr = PetscFree(recv_reqs);CHKERRQ(ierr); 353674ae819SStefano Zampini ierr = PetscFree(send_buffer);CHKERRQ(ierr); 354674ae819SStefano Zampini ierr = PetscFree(recv_buffer);CHKERRQ(ierr); 355674ae819SStefano Zampini ierr = PetscFree(ptrs_buffer);CHKERRQ(ierr); 35625d06dbeSstefano_zampini } 357674ae819SStefano Zampini 358674ae819SStefano Zampini /* Compute B and B_delta (local actions) */ 359785e854fSJed Brown ierr = PetscMalloc1(pcis->n_neigh,&aux_sums);CHKERRQ(ierr); 360785e854fSJed Brown ierr = PetscMalloc1(n_local_lambda,&l2g_indices);CHKERRQ(ierr); 361785e854fSJed Brown ierr = PetscMalloc1(n_local_lambda,&vals_B_delta);CHKERRQ(ierr); 362785e854fSJed Brown ierr = PetscMalloc1(n_local_lambda,&cols_B_delta);CHKERRQ(ierr); 36325d06dbeSstefano_zampini if (!pcbddc->use_deluxe_scaling) { 364785e854fSJed Brown ierr = PetscMalloc1(n_local_lambda,&scaling_factors);CHKERRQ(ierr); 36525d06dbeSstefano_zampini } else { 36625d06dbeSstefano_zampini scaling_factors = NULL; 36725d06dbeSstefano_zampini all_factors = NULL; 36825d06dbeSstefano_zampini } 369dc456d91SStefano Zampini ierr = ISGetIndices(subset_n,&aux_global_numbering);CHKERRQ(ierr); 370674ae819SStefano Zampini partial_sum=0; 371dc456d91SStefano Zampini cum = 0; 372674ae819SStefano Zampini for (i=0;i<dual_size;i++) { 373dc456d91SStefano Zampini n_global_lambda = aux_global_numbering[cum]; 374674ae819SStefano Zampini j = mat_graph->count[aux_local_numbering_1[i]]; 375674ae819SStefano Zampini aux_sums[0]=0; 376674ae819SStefano Zampini for (s=1;s<j;s++) { 377674ae819SStefano Zampini aux_sums[s]=aux_sums[s-1]+j-s+1; 378674ae819SStefano Zampini } 37925d06dbeSstefano_zampini if (all_factors) array = all_factors[aux_local_numbering_1[i]]; 380674ae819SStefano Zampini n_neg_values = 0; 3812a7da448SStefano Zampini while(n_neg_values < j && mat_graph->neighbours_set[aux_local_numbering_1[i]][n_neg_values] < rank) {n_neg_values++;} 382674ae819SStefano Zampini n_pos_values = j - n_neg_values; 383674ae819SStefano Zampini if (fully_redundant) { 384674ae819SStefano Zampini for (s=0;s<n_neg_values;s++) { 385674ae819SStefano Zampini l2g_indices [partial_sum+s]=aux_sums[s]+n_neg_values-s-1+n_global_lambda; 386674ae819SStefano Zampini cols_B_delta [partial_sum+s]=dual_dofs_boundary_indices[i]; 387674ae819SStefano Zampini vals_B_delta [partial_sum+s]=-1.0; 38825d06dbeSstefano_zampini if (!pcbddc->use_deluxe_scaling) scaling_factors[partial_sum+s]=array[s]; 389674ae819SStefano Zampini } 390674ae819SStefano Zampini for (s=0;s<n_pos_values;s++) { 391674ae819SStefano Zampini l2g_indices [partial_sum+s+n_neg_values]=aux_sums[n_neg_values]+s+n_global_lambda; 392674ae819SStefano Zampini cols_B_delta [partial_sum+s+n_neg_values]=dual_dofs_boundary_indices[i]; 393674ae819SStefano Zampini vals_B_delta [partial_sum+s+n_neg_values]=1.0; 39425d06dbeSstefano_zampini if (!pcbddc->use_deluxe_scaling) scaling_factors[partial_sum+s+n_neg_values]=array[s+n_neg_values]; 395674ae819SStefano Zampini } 396674ae819SStefano Zampini partial_sum += j; 397674ae819SStefano Zampini } else { 398674ae819SStefano Zampini /* l2g_indices and default cols and vals of B_delta */ 399674ae819SStefano Zampini for (s=0;s<j;s++) { 400674ae819SStefano Zampini l2g_indices [partial_sum+s]=n_global_lambda+s; 401674ae819SStefano Zampini cols_B_delta [partial_sum+s]=dual_dofs_boundary_indices[i]; 402674ae819SStefano Zampini vals_B_delta [partial_sum+s]=0.0; 403674ae819SStefano Zampini } 404674ae819SStefano Zampini /* B_delta */ 405674ae819SStefano Zampini if ( n_neg_values > 0 ) { /* there's a rank next to me to the left */ 406674ae819SStefano Zampini vals_B_delta [partial_sum+n_neg_values-1]=-1.0; 407674ae819SStefano Zampini } 408674ae819SStefano Zampini if ( n_neg_values < j ) { /* there's a rank next to me to the right */ 409674ae819SStefano Zampini vals_B_delta [partial_sum+n_neg_values]=1.0; 410674ae819SStefano Zampini } 411674ae819SStefano Zampini /* scaling as in Klawonn-Widlund 1999 */ 41225d06dbeSstefano_zampini if (!pcbddc->use_deluxe_scaling) { 413674ae819SStefano Zampini for (s=0;s<n_neg_values;s++) { 414674ae819SStefano Zampini scalar_value = 0.0; 41525d06dbeSstefano_zampini for (k=0;k<s+1;k++) scalar_value += array[k]; 416674ae819SStefano Zampini scaling_factors[partial_sum+s] = -scalar_value; 417674ae819SStefano Zampini } 418674ae819SStefano Zampini for (s=0;s<n_pos_values;s++) { 419674ae819SStefano Zampini scalar_value = 0.0; 42025d06dbeSstefano_zampini for (k=s+n_neg_values;k<j;k++) scalar_value += array[k]; 421674ae819SStefano Zampini scaling_factors[partial_sum+s+n_neg_values] = scalar_value; 422674ae819SStefano Zampini } 42325d06dbeSstefano_zampini } 424674ae819SStefano Zampini partial_sum += j; 425674ae819SStefano Zampini } 426dc456d91SStefano Zampini cum += aux_local_numbering_2[i]; 427674ae819SStefano Zampini } 428dc456d91SStefano Zampini ierr = ISRestoreIndices(subset_n,&aux_global_numbering);CHKERRQ(ierr); 429dc456d91SStefano Zampini ierr = ISDestroy(&subset_mult);CHKERRQ(ierr); 430dc456d91SStefano Zampini ierr = ISDestroy(&subset_n);CHKERRQ(ierr); 431674ae819SStefano Zampini ierr = PetscFree(aux_sums);CHKERRQ(ierr); 432674ae819SStefano Zampini ierr = PetscFree(aux_local_numbering_1);CHKERRQ(ierr); 433674ae819SStefano Zampini ierr = PetscFree(dual_dofs_boundary_indices);CHKERRQ(ierr); 43425d06dbeSstefano_zampini if (all_factors) { 435674ae819SStefano Zampini ierr = PetscFree(all_factors[0]);CHKERRQ(ierr); 436674ae819SStefano Zampini ierr = PetscFree(all_factors);CHKERRQ(ierr); 43725d06dbeSstefano_zampini } 438674ae819SStefano Zampini 439674ae819SStefano Zampini /* Create local part of B_delta */ 440302440fdSBarry Smith ierr = MatCreate(PETSC_COMM_SELF,&fetidpmat_ctx->B_delta);CHKERRQ(ierr); 441674ae819SStefano Zampini ierr = MatSetSizes(fetidpmat_ctx->B_delta,n_local_lambda,pcis->n_B,n_local_lambda,pcis->n_B);CHKERRQ(ierr); 442674ae819SStefano Zampini ierr = MatSetType(fetidpmat_ctx->B_delta,MATSEQAIJ);CHKERRQ(ierr); 443674ae819SStefano Zampini ierr = MatSeqAIJSetPreallocation(fetidpmat_ctx->B_delta,1,NULL);CHKERRQ(ierr); 444674ae819SStefano Zampini ierr = MatSetOption(fetidpmat_ctx->B_delta,MAT_IGNORE_ZERO_ENTRIES,PETSC_TRUE);CHKERRQ(ierr); 445674ae819SStefano Zampini for (i=0;i<n_local_lambda;i++) { 446674ae819SStefano Zampini ierr = MatSetValue(fetidpmat_ctx->B_delta,i,cols_B_delta[i],vals_B_delta[i],INSERT_VALUES);CHKERRQ(ierr); 447674ae819SStefano Zampini } 448674ae819SStefano Zampini ierr = PetscFree(vals_B_delta);CHKERRQ(ierr); 449674ae819SStefano Zampini ierr = MatAssemblyBegin(fetidpmat_ctx->B_delta,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 450674ae819SStefano Zampini ierr = MatAssemblyEnd(fetidpmat_ctx->B_delta,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 451674ae819SStefano Zampini 45225d06dbeSstefano_zampini BD1 = NULL; 45325d06dbeSstefano_zampini BD2 = NULL; 454674ae819SStefano Zampini if (fully_redundant) { 45525d06dbeSstefano_zampini if (pcbddc->use_deluxe_scaling) SETERRQ(comm,PETSC_ERR_SUP,"Deluxe FETIDP with fully-redundant multipliers to be implemented"); 456302440fdSBarry Smith ierr = MatCreate(PETSC_COMM_SELF,&ScalingMat);CHKERRQ(ierr); 457674ae819SStefano Zampini ierr = MatSetSizes(ScalingMat,n_local_lambda,n_local_lambda,n_local_lambda,n_local_lambda);CHKERRQ(ierr); 458674ae819SStefano Zampini ierr = MatSetType(ScalingMat,MATSEQAIJ);CHKERRQ(ierr); 459674ae819SStefano Zampini ierr = MatSeqAIJSetPreallocation(ScalingMat,1,NULL);CHKERRQ(ierr); 460674ae819SStefano Zampini for (i=0;i<n_local_lambda;i++) { 461674ae819SStefano Zampini ierr = MatSetValue(ScalingMat,i,i,scaling_factors[i],INSERT_VALUES);CHKERRQ(ierr); 462674ae819SStefano Zampini } 463674ae819SStefano Zampini ierr = MatAssemblyBegin(ScalingMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 464674ae819SStefano Zampini ierr = MatAssemblyEnd(ScalingMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 465674ae819SStefano Zampini ierr = MatMatMult(ScalingMat,fetidpmat_ctx->B_delta,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&fetidpmat_ctx->B_Ddelta);CHKERRQ(ierr); 466674ae819SStefano Zampini ierr = MatDestroy(&ScalingMat);CHKERRQ(ierr); 467674ae819SStefano Zampini } else { 468302440fdSBarry Smith ierr = MatCreate(PETSC_COMM_SELF,&fetidpmat_ctx->B_Ddelta);CHKERRQ(ierr); 469674ae819SStefano Zampini ierr = MatSetSizes(fetidpmat_ctx->B_Ddelta,n_local_lambda,pcis->n_B,n_local_lambda,pcis->n_B);CHKERRQ(ierr); 47025d06dbeSstefano_zampini if (!pcbddc->use_deluxe_scaling) { 471674ae819SStefano Zampini ierr = MatSetType(fetidpmat_ctx->B_Ddelta,MATSEQAIJ);CHKERRQ(ierr); 472674ae819SStefano Zampini ierr = MatSeqAIJSetPreallocation(fetidpmat_ctx->B_Ddelta,1,NULL);CHKERRQ(ierr); 473674ae819SStefano Zampini for (i=0;i<n_local_lambda;i++) { 474674ae819SStefano Zampini ierr = MatSetValue(fetidpmat_ctx->B_Ddelta,i,cols_B_delta[i],scaling_factors[i],INSERT_VALUES);CHKERRQ(ierr); 475674ae819SStefano Zampini } 476674ae819SStefano Zampini ierr = MatAssemblyBegin(fetidpmat_ctx->B_Ddelta,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 477674ae819SStefano Zampini ierr = MatAssemblyEnd(fetidpmat_ctx->B_Ddelta,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 47825d06dbeSstefano_zampini } else { 47925d06dbeSstefano_zampini /* scaling as in Klawonn-Widlund 1999 */ 48025d06dbeSstefano_zampini PCBDDCDeluxeScaling deluxe_ctx = pcbddc->deluxe_ctx; 48125d06dbeSstefano_zampini PCBDDCSubSchurs sub_schurs = pcbddc->sub_schurs; 48225d06dbeSstefano_zampini Mat T; 48325d06dbeSstefano_zampini PetscScalar *W,lwork,*Bwork; 48425d06dbeSstefano_zampini const PetscInt *idxs; 48525d06dbeSstefano_zampini PetscInt cum,mss,*nnz; 48625d06dbeSstefano_zampini PetscBLASInt *pivots,B_lwork,B_N,B_ierr; 48725d06dbeSstefano_zampini 48825d06dbeSstefano_zampini if (!pcbddc->deluxe_singlemat) SETERRQ(comm,PETSC_ERR_USER,"Cannot compute B_Ddelta! rerun with -pc_bddc_deluxe_singlemat"); 48925d06dbeSstefano_zampini if (deluxe_ctx->change) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Cannot compute B_Ddelta with deluxe scaling with active change context"); 49025d06dbeSstefano_zampini 49125d06dbeSstefano_zampini mss = 0; 49225d06dbeSstefano_zampini ierr = PetscCalloc1(pcis->n_B,&nnz);CHKERRQ(ierr); 49325d06dbeSstefano_zampini if (sub_schurs->is_Ej_all) { 49425d06dbeSstefano_zampini ierr = ISGetIndices(sub_schurs->is_Ej_all,&idxs);CHKERRQ(ierr); 49525d06dbeSstefano_zampini for (i=0,cum=0;i<sub_schurs->n_subs;i++) { 49625d06dbeSstefano_zampini PetscInt subset_size; 49725d06dbeSstefano_zampini 49825d06dbeSstefano_zampini ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr); 49925d06dbeSstefano_zampini for (j=0;j<subset_size;j++) nnz[idxs[j+cum]] = subset_size; 50025d06dbeSstefano_zampini mss = PetscMax(mss,subset_size); 50125d06dbeSstefano_zampini cum += subset_size; 50225d06dbeSstefano_zampini } 50325d06dbeSstefano_zampini } 50425d06dbeSstefano_zampini ierr = MatCreate(PETSC_COMM_SELF,&T);CHKERRQ(ierr); 50525d06dbeSstefano_zampini ierr = MatSetSizes(T,pcis->n_B,pcis->n_B,pcis->n_B,pcis->n_B);CHKERRQ(ierr); 50625d06dbeSstefano_zampini ierr = MatSetType(T,MATSEQAIJ);CHKERRQ(ierr); 50725d06dbeSstefano_zampini ierr = MatSeqAIJSetPreallocation(T,0,nnz);CHKERRQ(ierr); 50825d06dbeSstefano_zampini ierr = PetscFree(nnz);CHKERRQ(ierr); 50925d06dbeSstefano_zampini 51025d06dbeSstefano_zampini /* workspace allocation */ 5117f00194aSstefano_zampini B_lwork = 0; 5127f00194aSstefano_zampini if (mss) { 513be8bcbd0Sstefano_zampini PetscScalar dummy = 1; 514be8bcbd0Sstefano_zampini 51525d06dbeSstefano_zampini B_lwork = -1; 51625d06dbeSstefano_zampini ierr = PetscBLASIntCast(mss,&B_N);CHKERRQ(ierr); 51725d06dbeSstefano_zampini ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr); 518be8bcbd0Sstefano_zampini PetscStackCallBLAS("LAPACKgetri",LAPACKgetri_(&B_N,&dummy,&B_N,&B_N,&lwork,&B_lwork,&B_ierr)); 51925d06dbeSstefano_zampini ierr = PetscFPTrapPop();CHKERRQ(ierr); 52025d06dbeSstefano_zampini if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in query to GETRI Lapack routine %d",(int)B_ierr); 52125d06dbeSstefano_zampini ierr = PetscBLASIntCast((PetscInt)PetscRealPart(lwork),&B_lwork);CHKERRQ(ierr); 5227f00194aSstefano_zampini } 52325d06dbeSstefano_zampini ierr = PetscMalloc3(mss*mss,&W,mss,&pivots,B_lwork,&Bwork);CHKERRQ(ierr); 52425d06dbeSstefano_zampini 52525d06dbeSstefano_zampini for (i=0,cum=0;i<sub_schurs->n_subs;i++) { 52625d06dbeSstefano_zampini PetscScalar *M; 52725d06dbeSstefano_zampini PetscInt subset_size; 52825d06dbeSstefano_zampini 52925d06dbeSstefano_zampini ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr); 53025d06dbeSstefano_zampini ierr = PetscBLASIntCast(subset_size,&B_N);CHKERRQ(ierr); 53125d06dbeSstefano_zampini ierr = MatDenseGetArray(deluxe_ctx->seq_mat[i],&M);CHKERRQ(ierr); 53225d06dbeSstefano_zampini ierr = PetscMemcpy(W,M,subset_size*subset_size*sizeof(PetscScalar));CHKERRQ(ierr); 53325d06dbeSstefano_zampini ierr = MatDenseRestoreArray(deluxe_ctx->seq_mat[i],&M);CHKERRQ(ierr); 53425d06dbeSstefano_zampini ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr); 53525d06dbeSstefano_zampini PetscStackCallBLAS("LAPACKgetrf",LAPACKgetrf_(&B_N,&B_N,W,&B_N,pivots,&B_ierr)); 53625d06dbeSstefano_zampini if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRF Lapack routine %d",(int)B_ierr); 53725d06dbeSstefano_zampini PetscStackCallBLAS("LAPACKgetri",LAPACKgetri_(&B_N,W,&B_N,pivots,Bwork,&B_lwork,&B_ierr)); 53825d06dbeSstefano_zampini if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRI Lapack routine %d",(int)B_ierr); 53925d06dbeSstefano_zampini ierr = PetscFPTrapPop();CHKERRQ(ierr); 54025d06dbeSstefano_zampini ierr = MatSetValues(T,subset_size,idxs+cum,subset_size,idxs+cum,W,INSERT_VALUES);CHKERRQ(ierr); 54125d06dbeSstefano_zampini cum += subset_size; 54225d06dbeSstefano_zampini } 54325d06dbeSstefano_zampini ierr = MatAssemblyBegin(T,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 54425d06dbeSstefano_zampini ierr = MatAssemblyEnd(T,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 54525d06dbeSstefano_zampini ierr = MatMatTransposeMult(T,fetidpmat_ctx->B_delta,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&BD1);CHKERRQ(ierr); 54625d06dbeSstefano_zampini ierr = MatMatMult(fetidpmat_ctx->B_delta,BD1,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&BD2);CHKERRQ(ierr); 54725d06dbeSstefano_zampini ierr = MatDestroy(&T);CHKERRQ(ierr); 54825d06dbeSstefano_zampini ierr = PetscFree3(W,pivots,Bwork);CHKERRQ(ierr); 54925d06dbeSstefano_zampini if (sub_schurs->is_Ej_all) { 55025d06dbeSstefano_zampini ierr = ISRestoreIndices(sub_schurs->is_Ej_all,&idxs);CHKERRQ(ierr); 55125d06dbeSstefano_zampini } 55225d06dbeSstefano_zampini } 553674ae819SStefano Zampini } 554674ae819SStefano Zampini ierr = PetscFree(scaling_factors);CHKERRQ(ierr); 555674ae819SStefano Zampini ierr = PetscFree(cols_B_delta);CHKERRQ(ierr); 556674ae819SStefano Zampini 557457d4a33Sstefano_zampini /* Layout of multipliers */ 558457d4a33Sstefano_zampini ierr = PetscLayoutCreate(comm,&llay);CHKERRQ(ierr); 559457d4a33Sstefano_zampini ierr = PetscLayoutSetBlockSize(llay,1);CHKERRQ(ierr); 560457d4a33Sstefano_zampini ierr = PetscLayoutSetSize(llay,fetidpmat_ctx->n_lambda);CHKERRQ(ierr); 561457d4a33Sstefano_zampini ierr = PetscLayoutSetUp(llay);CHKERRQ(ierr); 562457d4a33Sstefano_zampini ierr = PetscLayoutGetLocalSize(llay,&fetidpmat_ctx->n);CHKERRQ(ierr); 563457d4a33Sstefano_zampini 564457d4a33Sstefano_zampini /* Local work vector of multipliers */ 565457d4a33Sstefano_zampini ierr = VecCreate(PETSC_COMM_SELF,&fetidpmat_ctx->lambda_local);CHKERRQ(ierr); 566457d4a33Sstefano_zampini ierr = VecSetSizes(fetidpmat_ctx->lambda_local,n_local_lambda,n_local_lambda);CHKERRQ(ierr); 567457d4a33Sstefano_zampini ierr = VecSetType(fetidpmat_ctx->lambda_local,VECSEQ);CHKERRQ(ierr); 568457d4a33Sstefano_zampini 56925d06dbeSstefano_zampini if (BD2) { 57025d06dbeSstefano_zampini ISLocalToGlobalMapping l2g; 57125d06dbeSstefano_zampini Mat T,TA,*pT; 57225d06dbeSstefano_zampini IS is; 57325d06dbeSstefano_zampini PetscInt nl,N; 57425d06dbeSstefano_zampini BDdelta_DN ctx; 57525d06dbeSstefano_zampini 57625d06dbeSstefano_zampini ierr = PetscLayoutGetLocalSize(llay,&nl);CHKERRQ(ierr); 57725d06dbeSstefano_zampini ierr = PetscLayoutGetSize(llay,&N);CHKERRQ(ierr); 57825d06dbeSstefano_zampini ierr = MatCreate(comm,&T);CHKERRQ(ierr); 57925d06dbeSstefano_zampini ierr = MatSetSizes(T,nl,nl,N,N);CHKERRQ(ierr); 58025d06dbeSstefano_zampini ierr = MatSetType(T,MATIS);CHKERRQ(ierr); 58125d06dbeSstefano_zampini ierr = ISLocalToGlobalMappingCreate(comm,1,n_local_lambda,l2g_indices,PETSC_COPY_VALUES,&l2g);CHKERRQ(ierr); 58225d06dbeSstefano_zampini ierr = MatSetLocalToGlobalMapping(T,l2g,l2g);CHKERRQ(ierr); 58325d06dbeSstefano_zampini ierr = ISLocalToGlobalMappingDestroy(&l2g);CHKERRQ(ierr); 58425d06dbeSstefano_zampini ierr = MatISSetLocalMat(T,BD2);CHKERRQ(ierr); 58525d06dbeSstefano_zampini ierr = MatDestroy(&BD2);CHKERRQ(ierr); 58625d06dbeSstefano_zampini ierr = MatISGetMPIXAIJ(T,MAT_INITIAL_MATRIX,&TA);CHKERRQ(ierr); 58725d06dbeSstefano_zampini ierr = MatDestroy(&T);CHKERRQ(ierr); 58825d06dbeSstefano_zampini ierr = ISCreateGeneral(comm,n_local_lambda,l2g_indices,PETSC_USE_POINTER,&is);CHKERRQ(ierr); 5897dae84e0SHong Zhang ierr = MatCreateSubMatrices(TA,1,&is,&is,MAT_INITIAL_MATRIX,&pT);CHKERRQ(ierr); 59025d06dbeSstefano_zampini ierr = MatDestroy(&TA);CHKERRQ(ierr); 59125d06dbeSstefano_zampini ierr = ISDestroy(&is);CHKERRQ(ierr); 59225d06dbeSstefano_zampini BD2 = pT[0]; 59325d06dbeSstefano_zampini ierr = PetscFree(pT);CHKERRQ(ierr); 59425d06dbeSstefano_zampini 59525d06dbeSstefano_zampini /* B_Ddelta for non-redundant multipliers with deluxe scaling */ 59625d06dbeSstefano_zampini ierr = PetscNew(&ctx);CHKERRQ(ierr); 59725d06dbeSstefano_zampini ierr = MatSetType(fetidpmat_ctx->B_Ddelta,MATSHELL);CHKERRQ(ierr); 59825d06dbeSstefano_zampini ierr = MatShellSetContext(fetidpmat_ctx->B_Ddelta,(void *)ctx);CHKERRQ(ierr); 59925d06dbeSstefano_zampini ierr = MatShellSetOperation(fetidpmat_ctx->B_Ddelta,MATOP_MULT,(void (*)(void))MatMult_BDdelta_deluxe_nonred);CHKERRQ(ierr); 60025d06dbeSstefano_zampini ierr = MatShellSetOperation(fetidpmat_ctx->B_Ddelta,MATOP_MULT_TRANSPOSE,(void (*)(void))MatMultTranspose_BDdelta_deluxe_nonred);CHKERRQ(ierr); 60125d06dbeSstefano_zampini ierr = MatShellSetOperation(fetidpmat_ctx->B_Ddelta,MATOP_DESTROY,(void (*)(void))MatDestroy_BDdelta_deluxe_nonred);CHKERRQ(ierr); 60225d06dbeSstefano_zampini ierr = MatSetUp(fetidpmat_ctx->B_Ddelta);CHKERRQ(ierr); 60325d06dbeSstefano_zampini 60425d06dbeSstefano_zampini ierr = PetscObjectReference((PetscObject)BD1);CHKERRQ(ierr); 60525d06dbeSstefano_zampini ctx->BD = BD1; 60625d06dbeSstefano_zampini ierr = KSPCreate(PETSC_COMM_SELF,&ctx->kBD);CHKERRQ(ierr); 60725d06dbeSstefano_zampini ierr = KSPSetOperators(ctx->kBD,BD2,BD2);CHKERRQ(ierr); 60825d06dbeSstefano_zampini ierr = VecDuplicate(fetidpmat_ctx->lambda_local,&ctx->work);CHKERRQ(ierr); 60925d06dbeSstefano_zampini fetidpmat_ctx->deluxe_nonred = PETSC_TRUE; 61025d06dbeSstefano_zampini } 61125d06dbeSstefano_zampini ierr = MatDestroy(&BD1);CHKERRQ(ierr); 61225d06dbeSstefano_zampini ierr = MatDestroy(&BD2);CHKERRQ(ierr); 61325d06dbeSstefano_zampini 61425d06dbeSstefano_zampini /* fetidpmat sizes */ 61525d06dbeSstefano_zampini fetidpmat_ctx->n += nPgl; 61625d06dbeSstefano_zampini fetidpmat_ctx->N = fetidpmat_ctx->n_lambda+nPg; 61725d06dbeSstefano_zampini 618457d4a33Sstefano_zampini /* Global vector for FETI-DP linear system */ 619457d4a33Sstefano_zampini ierr = VecCreate(comm,&fetidp_global);CHKERRQ(ierr); 620457d4a33Sstefano_zampini ierr = VecSetSizes(fetidp_global,fetidpmat_ctx->n,fetidpmat_ctx->N);CHKERRQ(ierr); 621457d4a33Sstefano_zampini ierr = VecSetType(fetidp_global,VECMPI);CHKERRQ(ierr); 622457d4a33Sstefano_zampini ierr = VecSetUp(fetidp_global);CHKERRQ(ierr); 623457d4a33Sstefano_zampini 6249eec4de8Sstefano_zampini /* Decide layout for fetidp dofs: if it is a saddle point problem 6259eec4de8Sstefano_zampini pressure is ordered first in the local part of the global vector 626457d4a33Sstefano_zampini of the FETI-DP linear system */ 627457d4a33Sstefano_zampini if (nPg) { 628af140850Sstefano_zampini IS IS_l2g_p,ais; 629457d4a33Sstefano_zampini PetscLayout alay; 630457d4a33Sstefano_zampini const PetscInt *idxs,*pranges,*aranges,*lranges; 631af140850Sstefano_zampini PetscInt *l2g_indices_p,rst; 632457d4a33Sstefano_zampini 633457d4a33Sstefano_zampini ierr = PetscMalloc1(nPl,&l2g_indices_p);CHKERRQ(ierr); 634457d4a33Sstefano_zampini ierr = VecGetLayout(fetidp_global,&alay);CHKERRQ(ierr); 635457d4a33Sstefano_zampini ierr = PetscLayoutGetRanges(alay,&aranges);CHKERRQ(ierr); 636457d4a33Sstefano_zampini ierr = PetscLayoutGetRanges(play,&pranges);CHKERRQ(ierr); 637457d4a33Sstefano_zampini ierr = PetscLayoutGetRanges(llay,&lranges);CHKERRQ(ierr); 638457d4a33Sstefano_zampini ierr = ISLocalToGlobalMappingGetIndices(l2gmap_p,&idxs);CHKERRQ(ierr); 639af140850Sstefano_zampini /* shift local to global indices for pressure */ 640457d4a33Sstefano_zampini for (i=0;i<nPl;i++) { 641457d4a33Sstefano_zampini PetscInt owner; 642457d4a33Sstefano_zampini 643457d4a33Sstefano_zampini ierr = PetscLayoutFindOwner(play,idxs[i],&owner);CHKERRQ(ierr); 644457d4a33Sstefano_zampini l2g_indices_p[i] = idxs[i]-pranges[owner]+aranges[owner]; 645457d4a33Sstefano_zampini } 646457d4a33Sstefano_zampini ierr = ISLocalToGlobalMappingRestoreIndices(l2gmap_p,&idxs);CHKERRQ(ierr); 647457d4a33Sstefano_zampini ierr = ISCreateGeneral(comm,nPl,l2g_indices_p,PETSC_OWN_POINTER,&IS_l2g_p);CHKERRQ(ierr); 648af140850Sstefano_zampini 649457d4a33Sstefano_zampini /* local to global scatter for interface pressure */ 650457d4a33Sstefano_zampini ierr = VecScatterCreate(fetidpmat_ctx->vP,NULL,fetidp_global,IS_l2g_p,&fetidpmat_ctx->l2g_p);CHKERRQ(ierr); 651457d4a33Sstefano_zampini ierr = ISDestroy(&IS_l2g_p);CHKERRQ(ierr); 652457d4a33Sstefano_zampini 653af140850Sstefano_zampini /* shift local to global indices for multipliers */ 654457d4a33Sstefano_zampini for (i=0;i<n_local_lambda;i++) { 655457d4a33Sstefano_zampini PetscInt owner,ps; 656457d4a33Sstefano_zampini 657457d4a33Sstefano_zampini ierr = PetscLayoutFindOwner(llay,l2g_indices[i],&owner);CHKERRQ(ierr); 658457d4a33Sstefano_zampini ps = pranges[owner+1]-pranges[owner]; 659457d4a33Sstefano_zampini l2g_indices[i] = l2g_indices[i]-lranges[owner]+aranges[owner]+ps; 660457d4a33Sstefano_zampini } 661457d4a33Sstefano_zampini 6629cc7774eSstefano_zampini /* scatter from alldofs to interface pressures global fetidp vector */ 6639cc7774eSstefano_zampini ierr = PetscLayoutGetRange(alay,&rst,NULL);CHKERRQ(ierr); 6649cc7774eSstefano_zampini ierr = ISCreateStride(comm,nPgl,rst,1,&ais);CHKERRQ(ierr); 665af140850Sstefano_zampini ierr = VecScatterCreate(pcis->vec1_global,pP,fetidp_global,ais,&fetidpmat_ctx->g2g_p);CHKERRQ(ierr); 6669cc7774eSstefano_zampini ierr = ISDestroy(&ais);CHKERRQ(ierr); 667457d4a33Sstefano_zampini } 668457d4a33Sstefano_zampini ierr = PetscLayoutDestroy(&llay);CHKERRQ(ierr); 669457d4a33Sstefano_zampini ierr = PetscLayoutDestroy(&play);CHKERRQ(ierr); 670457d4a33Sstefano_zampini ierr = ISCreateGeneral(comm,n_local_lambda,l2g_indices,PETSC_OWN_POINTER,&IS_l2g_lambda);CHKERRQ(ierr); 671a1c0d0daSstefano_zampini 6729cc7774eSstefano_zampini /* scatter from local to global multipliers */ 673457d4a33Sstefano_zampini ierr = VecScatterCreate(fetidpmat_ctx->lambda_local,NULL,fetidp_global,IS_l2g_lambda,&fetidpmat_ctx->l2g_lambda);CHKERRQ(ierr); 674457d4a33Sstefano_zampini ierr = ISDestroy(&IS_l2g_lambda);CHKERRQ(ierr); 675457d4a33Sstefano_zampini ierr = ISLocalToGlobalMappingDestroy(&l2gmap_p);CHKERRQ(ierr); 676a1c0d0daSstefano_zampini ierr = VecDestroy(&fetidp_global);CHKERRQ(ierr); 677457d4a33Sstefano_zampini 678a1c0d0daSstefano_zampini /* Create some work vectors needed by fetidp */ 679674ae819SStefano Zampini ierr = VecDuplicate(pcis->vec1_B,&fetidpmat_ctx->temp_solution_B);CHKERRQ(ierr); 680674ae819SStefano Zampini ierr = VecDuplicate(pcis->vec1_D,&fetidpmat_ctx->temp_solution_D);CHKERRQ(ierr); 681674ae819SStefano Zampini PetscFunctionReturn(0); 682674ae819SStefano Zampini } 683674ae819SStefano Zampini 68425d06dbeSstefano_zampini 685674ae819SStefano Zampini PetscErrorCode PCBDDCSetupFETIDPPCContext(Mat fetimat, FETIDPPC_ctx fetidppc_ctx) 686674ae819SStefano Zampini { 687674ae819SStefano Zampini FETIDPMat_ctx mat_ctx; 68825d06dbeSstefano_zampini PC_BDDC *pcbddc = (PC_BDDC*)fetidppc_ctx->pc->data; 68925d06dbeSstefano_zampini PC_IS *pcis = (PC_IS*)fetidppc_ctx->pc->data; 690f28b6018SStefano Zampini PetscBool lumped = PETSC_FALSE; 691674ae819SStefano Zampini PetscErrorCode ierr; 692674ae819SStefano Zampini 693674ae819SStefano Zampini PetscFunctionBegin; 694674ae819SStefano Zampini ierr = MatShellGetContext(fetimat,(void**)&mat_ctx);CHKERRQ(ierr); 695674ae819SStefano Zampini /* get references from objects created when setting up feti mat context */ 696674ae819SStefano Zampini ierr = PetscObjectReference((PetscObject)mat_ctx->lambda_local);CHKERRQ(ierr); 697674ae819SStefano Zampini fetidppc_ctx->lambda_local = mat_ctx->lambda_local; 698674ae819SStefano Zampini ierr = PetscObjectReference((PetscObject)mat_ctx->B_Ddelta);CHKERRQ(ierr); 699674ae819SStefano Zampini fetidppc_ctx->B_Ddelta = mat_ctx->B_Ddelta; 70025d06dbeSstefano_zampini if (mat_ctx->deluxe_nonred) { 70125d06dbeSstefano_zampini PC pc,mpc; 70225d06dbeSstefano_zampini BDdelta_DN ctx; 70325d06dbeSstefano_zampini MatSolverPackage solver; 70425d06dbeSstefano_zampini const char *prefix; 70525d06dbeSstefano_zampini 70625d06dbeSstefano_zampini ierr = MatShellGetContext(mat_ctx->B_Ddelta,&ctx);CHKERRQ(ierr); 70725d06dbeSstefano_zampini ierr = KSPSetType(ctx->kBD,KSPPREONLY);CHKERRQ(ierr); 70825d06dbeSstefano_zampini ierr = KSPGetPC(ctx->kBD,&mpc);CHKERRQ(ierr); 70925d06dbeSstefano_zampini ierr = KSPGetPC(pcbddc->ksp_D,&pc);CHKERRQ(ierr); 71025d06dbeSstefano_zampini ierr = PCSetType(mpc,PCLU);CHKERRQ(ierr); 71125d06dbeSstefano_zampini ierr = PCFactorGetMatSolverPackage(pc,(const MatSolverPackage*)&solver);CHKERRQ(ierr); 71225d06dbeSstefano_zampini if (solver) { 71325d06dbeSstefano_zampini ierr = PCFactorSetMatSolverPackage(mpc,solver);CHKERRQ(ierr); 71425d06dbeSstefano_zampini } 71525d06dbeSstefano_zampini ierr = MatGetOptionsPrefix(fetimat,&prefix);CHKERRQ(ierr); 71625d06dbeSstefano_zampini ierr = KSPSetOptionsPrefix(ctx->kBD,prefix);CHKERRQ(ierr); 71725d06dbeSstefano_zampini ierr = KSPAppendOptionsPrefix(ctx->kBD,"bddelta_");CHKERRQ(ierr); 71825d06dbeSstefano_zampini ierr = KSPSetFromOptions(ctx->kBD);CHKERRQ(ierr); 71925d06dbeSstefano_zampini } 72025d06dbeSstefano_zampini 721674ae819SStefano Zampini ierr = PetscObjectReference((PetscObject)mat_ctx->l2g_lambda);CHKERRQ(ierr); 722674ae819SStefano Zampini fetidppc_ctx->l2g_lambda = mat_ctx->l2g_lambda; 723f28b6018SStefano Zampini /* Dirichlet preconditioner */ 7249c2d02cdSstefano_zampini ierr = PetscOptionsGetBool(NULL,((PetscObject)fetimat)->prefix,"-pc_lumped",&lumped,NULL);CHKERRQ(ierr); 725f28b6018SStefano Zampini if (!lumped) { 726*9feb908dSstefano_zampini IS iV; 7279c2d02cdSstefano_zampini PetscBool discrete_harmonic = PETSC_FALSE; 7289c2d02cdSstefano_zampini 729*9feb908dSstefano_zampini ierr = PetscObjectQuery((PetscObject)fetidppc_ctx->pc,"__KSPFETIDP_iV",(PetscObject*)&iV);CHKERRQ(ierr); 730*9feb908dSstefano_zampini if (iV) { 7319c2d02cdSstefano_zampini ierr = PetscOptionsGetBool(NULL,((PetscObject)fetimat)->prefix,"-pc_discrete_harmonic",&discrete_harmonic,NULL);CHKERRQ(ierr); 7329c2d02cdSstefano_zampini } 7339c2d02cdSstefano_zampini if (discrete_harmonic) { 7349c2d02cdSstefano_zampini KSP sksp; 7359c2d02cdSstefano_zampini PC pc; 7369c2d02cdSstefano_zampini Mat A_II,A_IB,A_BI; 7379c2d02cdSstefano_zampini IS aB; 7389c2d02cdSstefano_zampini PetscInt nb; 7395334bea6Sstefano_zampini PetscBool isshell; 7409c2d02cdSstefano_zampini KSPType ksptype; 7419c2d02cdSstefano_zampini const char *prefix; 7429c2d02cdSstefano_zampini 7439c2d02cdSstefano_zampini /* 7449c2d02cdSstefano_zampini We constructs a Schur complement for 7459c2d02cdSstefano_zampini 7469c2d02cdSstefano_zampini | A_II A_ID | 7479c2d02cdSstefano_zampini | A_DI A_DD | 7489c2d02cdSstefano_zampini 7499c2d02cdSstefano_zampini instead of 7509c2d02cdSstefano_zampini 7519c2d02cdSstefano_zampini | A_II B^t_II A_ID | 7529c2d02cdSstefano_zampini | B_II -C_II B_ID | 7539c2d02cdSstefano_zampini | A_DI B^t_ID A_DD | 7549c2d02cdSstefano_zampini 7559c2d02cdSstefano_zampini */ 7569c2d02cdSstefano_zampini ierr = ISGetLocalSize(pcis->is_B_local,&nb);CHKERRQ(ierr); 7579c2d02cdSstefano_zampini ierr = ISCreateStride(PetscObjectComm((PetscObject)pcis->A_II),nb,0,1,&aB);CHKERRQ(ierr); 758*9feb908dSstefano_zampini ierr = MatCreateSubMatrix(pcis->A_II,iV,iV,MAT_INITIAL_MATRIX,&A_II);CHKERRQ(ierr); 759*9feb908dSstefano_zampini ierr = MatCreateSubMatrix(pcis->A_IB,iV,aB,MAT_INITIAL_MATRIX,&A_IB);CHKERRQ(ierr); 760*9feb908dSstefano_zampini ierr = MatCreateSubMatrix(pcis->A_BI,aB,iV,MAT_INITIAL_MATRIX,&A_BI);CHKERRQ(ierr); 7619c2d02cdSstefano_zampini ierr = MatCreateSchurComplement(A_II,A_II,A_IB,A_BI,pcis->A_BB,&fetidppc_ctx->S_j);CHKERRQ(ierr); 7629c2d02cdSstefano_zampini 7639c2d02cdSstefano_zampini /* propagate settings of solver */ 7649c2d02cdSstefano_zampini ierr = MatSchurComplementGetKSP(fetidppc_ctx->S_j,&sksp);CHKERRQ(ierr); 7659c2d02cdSstefano_zampini ierr = KSPGetType(pcis->ksp_D,&ksptype);CHKERRQ(ierr); 7669c2d02cdSstefano_zampini ierr = KSPSetType(sksp,ksptype);CHKERRQ(ierr); 7679c2d02cdSstefano_zampini ierr = KSPGetPC(pcis->ksp_D,&pc);CHKERRQ(ierr); 7685334bea6Sstefano_zampini ierr = PetscObjectTypeCompare((PetscObject)pc,PCSHELL,&isshell);CHKERRQ(ierr); 7695334bea6Sstefano_zampini if (!isshell) { 7705334bea6Sstefano_zampini MatSolverPackage solver; 7715334bea6Sstefano_zampini PCType pctype; 7725334bea6Sstefano_zampini 7739c2d02cdSstefano_zampini ierr = PCGetType(pc,&pctype);CHKERRQ(ierr); 7749c2d02cdSstefano_zampini ierr = PCFactorGetMatSolverPackage(pc,(const MatSolverPackage*)&solver);CHKERRQ(ierr); 7759c2d02cdSstefano_zampini ierr = KSPGetPC(sksp,&pc);CHKERRQ(ierr); 7769c2d02cdSstefano_zampini ierr = PCSetType(pc,pctype);CHKERRQ(ierr); 7779c2d02cdSstefano_zampini if (solver) { 7789c2d02cdSstefano_zampini ierr = PCFactorSetMatSolverPackage(pc,solver);CHKERRQ(ierr); 7799c2d02cdSstefano_zampini } 7805334bea6Sstefano_zampini } else { 7815334bea6Sstefano_zampini ierr = KSPGetPC(sksp,&pc);CHKERRQ(ierr); 7825334bea6Sstefano_zampini ierr = PCSetType(pc,PCLU);CHKERRQ(ierr); 7835334bea6Sstefano_zampini } 7849c2d02cdSstefano_zampini ierr = MatDestroy(&A_II);CHKERRQ(ierr); 7859c2d02cdSstefano_zampini ierr = MatDestroy(&A_IB);CHKERRQ(ierr); 7869c2d02cdSstefano_zampini ierr = MatDestroy(&A_BI);CHKERRQ(ierr); 7879c2d02cdSstefano_zampini ierr = ISDestroy(&aB);CHKERRQ(ierr); 7889c2d02cdSstefano_zampini ierr = MatGetOptionsPrefix(fetimat,&prefix);CHKERRQ(ierr); 7899c2d02cdSstefano_zampini ierr = KSPSetOptionsPrefix(sksp,prefix);CHKERRQ(ierr); 7909c2d02cdSstefano_zampini ierr = KSPAppendOptionsPrefix(sksp,"harmonic_");CHKERRQ(ierr); 7913016320fSstefano_zampini ierr = KSPSetFromOptions(sksp);CHKERRQ(ierr); 7929c2d02cdSstefano_zampini } else { /* default Dirichlet preconditioner is pde-harmonic */ 793ed6c3d69SStefano Zampini ierr = MatCreateSchurComplement(pcis->A_II,pcis->A_II,pcis->A_IB,pcis->A_BI,pcis->A_BB,&fetidppc_ctx->S_j);CHKERRQ(ierr); 794ed6c3d69SStefano Zampini ierr = MatSchurComplementSetKSP(fetidppc_ctx->S_j,pcis->ksp_D);CHKERRQ(ierr); 7959c2d02cdSstefano_zampini } 796f28b6018SStefano Zampini } else { 797f28b6018SStefano Zampini ierr = PetscObjectReference((PetscObject)pcis->A_BB);CHKERRQ(ierr); 798f28b6018SStefano Zampini fetidppc_ctx->S_j = pcis->A_BB; 799f28b6018SStefano Zampini } 800af140850Sstefano_zampini /* saddle-point */ 801af140850Sstefano_zampini if (mat_ctx->xPg) { 802af140850Sstefano_zampini ierr = PetscObjectReference((PetscObject)mat_ctx->xPg);CHKERRQ(ierr); 803af140850Sstefano_zampini fetidppc_ctx->xPg = mat_ctx->xPg; 804af140850Sstefano_zampini ierr = PetscObjectReference((PetscObject)mat_ctx->yPg);CHKERRQ(ierr); 805af140850Sstefano_zampini fetidppc_ctx->yPg = mat_ctx->yPg; 8066cc1294bSstefano_zampini ierr = PetscObjectQuery((PetscObject)fetidppc_ctx->pc,"__KSPFETIDP_PKSP",(PetscObject*)&fetidppc_ctx->kP);CHKERRQ(ierr); 8076cc1294bSstefano_zampini ierr = PetscObjectReference((PetscObject)fetidppc_ctx->kP);CHKERRQ(ierr); 808af140850Sstefano_zampini } 809674ae819SStefano Zampini PetscFunctionReturn(0); 810674ae819SStefano Zampini } 811674ae819SStefano Zampini 812674ae819SStefano Zampini PetscErrorCode FETIDPMatMult(Mat fetimat, Vec x, Vec y) 813674ae819SStefano Zampini { 814674ae819SStefano Zampini FETIDPMat_ctx mat_ctx; 815617d11aeSStefano Zampini PC_BDDC *pcbddc; 816674ae819SStefano Zampini PC_IS *pcis; 817674ae819SStefano Zampini PetscErrorCode ierr; 818674ae819SStefano Zampini 819674ae819SStefano Zampini PetscFunctionBegin; 820674ae819SStefano Zampini ierr = MatShellGetContext(fetimat,(void**)&mat_ctx);CHKERRQ(ierr); 821674ae819SStefano Zampini pcis = (PC_IS*)mat_ctx->pc->data; 822617d11aeSStefano Zampini pcbddc = (PC_BDDC*)mat_ctx->pc->data; 823674ae819SStefano Zampini /* Application of B_delta^T */ 824af140850Sstefano_zampini ierr = VecSet(pcis->vec1_B,0.);CHKERRQ(ierr); 825674ae819SStefano Zampini ierr = VecScatterBegin(mat_ctx->l2g_lambda,x,mat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 826674ae819SStefano Zampini ierr = VecScatterEnd(mat_ctx->l2g_lambda,x,mat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 827674ae819SStefano Zampini ierr = MatMultTranspose(mat_ctx->B_delta,mat_ctx->lambda_local,pcis->vec1_B);CHKERRQ(ierr); 828af140850Sstefano_zampini 829af140850Sstefano_zampini /* Add contribution from saddle point */ 830af140850Sstefano_zampini if (mat_ctx->l2g_p) { 831af140850Sstefano_zampini ierr = VecScatterBegin(mat_ctx->l2g_p,x,mat_ctx->vP,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 832af140850Sstefano_zampini ierr = VecScatterEnd(mat_ctx->l2g_p,x,mat_ctx->vP,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 833af140850Sstefano_zampini if (pcbddc->switch_static) { 834af140850Sstefano_zampini ierr = MatMult(mat_ctx->Bt_BI,mat_ctx->vP,pcis->vec1_D);CHKERRQ(ierr); 835af140850Sstefano_zampini } 836af140850Sstefano_zampini ierr = MatMultAdd(mat_ctx->Bt_BB,mat_ctx->vP,pcis->vec1_B,pcis->vec1_B);CHKERRQ(ierr); 837af140850Sstefano_zampini } else { 838af140850Sstefano_zampini if (pcbddc->switch_static) { 839674ae819SStefano Zampini ierr = VecSet(pcis->vec1_D,0.0);CHKERRQ(ierr); 840af140850Sstefano_zampini } 841af140850Sstefano_zampini } 842af140850Sstefano_zampini /* Application of \widetilde{S}^-1 */ 843617d11aeSStefano Zampini ierr = PetscMemzero(pcbddc->benign_p0,pcbddc->benign_n*sizeof(PetscScalar));CHKERRQ(ierr); 844dc359a40SStefano Zampini ierr = PCBDDCApplyInterfacePreconditioner(mat_ctx->pc,PETSC_FALSE);CHKERRQ(ierr); 845c7ffc8ceSStefano Zampini ierr = PetscMemzero(pcbddc->benign_p0,pcbddc->benign_n*sizeof(PetscScalar));CHKERRQ(ierr); 846af140850Sstefano_zampini ierr = VecSet(y,0.0);CHKERRQ(ierr); 847674ae819SStefano Zampini /* Application of B_delta */ 848674ae819SStefano Zampini ierr = MatMult(mat_ctx->B_delta,pcis->vec1_B,mat_ctx->lambda_local);CHKERRQ(ierr); 849af140850Sstefano_zampini /* Contribution from boundary pressures */ 850af140850Sstefano_zampini if (mat_ctx->C) { 851af140850Sstefano_zampini const PetscScalar *lx; 852af140850Sstefano_zampini PetscScalar *ly; 853af140850Sstefano_zampini 854af140850Sstefano_zampini /* pressures ordered first in x and y */ 855af140850Sstefano_zampini ierr = VecGetArrayRead(x,&lx);CHKERRQ(ierr); 856af140850Sstefano_zampini ierr = VecGetArray(y,&ly);CHKERRQ(ierr); 857af140850Sstefano_zampini ierr = VecPlaceArray(mat_ctx->xPg,lx);CHKERRQ(ierr); 858af140850Sstefano_zampini ierr = VecPlaceArray(mat_ctx->yPg,ly);CHKERRQ(ierr); 859af140850Sstefano_zampini ierr = MatMult(mat_ctx->C,mat_ctx->xPg,mat_ctx->yPg);CHKERRQ(ierr); 860af140850Sstefano_zampini ierr = VecResetArray(mat_ctx->xPg);CHKERRQ(ierr); 861af140850Sstefano_zampini ierr = VecResetArray(mat_ctx->yPg);CHKERRQ(ierr); 862af140850Sstefano_zampini ierr = VecRestoreArrayRead(x,&lx);CHKERRQ(ierr); 863af140850Sstefano_zampini ierr = VecRestoreArray(y,&ly);CHKERRQ(ierr); 864af140850Sstefano_zampini } 865af140850Sstefano_zampini /* Add contribution from saddle point */ 866af140850Sstefano_zampini if (mat_ctx->l2g_p) { 867af140850Sstefano_zampini ierr = MatMult(mat_ctx->B_BB,pcis->vec1_B,mat_ctx->vP);CHKERRQ(ierr); 868af140850Sstefano_zampini if (pcbddc->switch_static) { 869af140850Sstefano_zampini ierr = MatMultAdd(mat_ctx->B_BI,pcis->vec1_D,mat_ctx->vP,mat_ctx->vP);CHKERRQ(ierr); 870af140850Sstefano_zampini } 871af140850Sstefano_zampini ierr = VecScatterBegin(mat_ctx->l2g_p,mat_ctx->vP,y,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 872af140850Sstefano_zampini ierr = VecScatterEnd(mat_ctx->l2g_p,mat_ctx->vP,y,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 873af140850Sstefano_zampini } 874674ae819SStefano Zampini ierr = VecScatterBegin(mat_ctx->l2g_lambda,mat_ctx->lambda_local,y,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 875674ae819SStefano Zampini ierr = VecScatterEnd(mat_ctx->l2g_lambda,mat_ctx->lambda_local,y,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 876674ae819SStefano Zampini PetscFunctionReturn(0); 877674ae819SStefano Zampini } 878674ae819SStefano Zampini 879edf7251bSStefano Zampini PetscErrorCode FETIDPMatMultTranspose(Mat fetimat, Vec x, Vec y) 880edf7251bSStefano Zampini { 881edf7251bSStefano Zampini FETIDPMat_ctx mat_ctx; 882edf7251bSStefano Zampini PC_IS *pcis; 883edf7251bSStefano Zampini PetscErrorCode ierr; 884edf7251bSStefano Zampini 885edf7251bSStefano Zampini PetscFunctionBegin; 886edf7251bSStefano Zampini ierr = MatShellGetContext(fetimat,(void**)&mat_ctx);CHKERRQ(ierr); 887edf7251bSStefano Zampini pcis = (PC_IS*)mat_ctx->pc->data; 888edf7251bSStefano Zampini /* Application of B_delta^T */ 889edf7251bSStefano Zampini ierr = VecScatterBegin(mat_ctx->l2g_lambda,x,mat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 890edf7251bSStefano Zampini ierr = VecScatterEnd(mat_ctx->l2g_lambda,x,mat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 891edf7251bSStefano Zampini ierr = MatMultTranspose(mat_ctx->B_delta,mat_ctx->lambda_local,pcis->vec1_B);CHKERRQ(ierr); 892edf7251bSStefano Zampini /* Application of \widetilde{S}^-1 */ 893edf7251bSStefano Zampini ierr = VecSet(pcis->vec1_D,0.0);CHKERRQ(ierr); 894edf7251bSStefano Zampini ierr = PCBDDCApplyInterfacePreconditioner(mat_ctx->pc,PETSC_TRUE);CHKERRQ(ierr); 895edf7251bSStefano Zampini /* Application of B_delta */ 896edf7251bSStefano Zampini ierr = MatMult(mat_ctx->B_delta,pcis->vec1_B,mat_ctx->lambda_local);CHKERRQ(ierr); 897edf7251bSStefano Zampini ierr = VecSet(y,0.0);CHKERRQ(ierr); 898edf7251bSStefano Zampini ierr = VecScatterBegin(mat_ctx->l2g_lambda,mat_ctx->lambda_local,y,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 899edf7251bSStefano Zampini ierr = VecScatterEnd(mat_ctx->l2g_lambda,mat_ctx->lambda_local,y,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 900edf7251bSStefano Zampini PetscFunctionReturn(0); 901edf7251bSStefano Zampini } 902edf7251bSStefano Zampini 903674ae819SStefano Zampini PetscErrorCode FETIDPPCApply(PC fetipc, Vec x, Vec y) 904674ae819SStefano Zampini { 905674ae819SStefano Zampini FETIDPPC_ctx pc_ctx; 906674ae819SStefano Zampini PC_IS *pcis; 907674ae819SStefano Zampini PetscErrorCode ierr; 908674ae819SStefano Zampini 909674ae819SStefano Zampini PetscFunctionBegin; 910302440fdSBarry Smith ierr = PCShellGetContext(fetipc,(void**)&pc_ctx);CHKERRQ(ierr); 911674ae819SStefano Zampini pcis = (PC_IS*)pc_ctx->pc->data; 912674ae819SStefano Zampini /* Application of B_Ddelta^T */ 913674ae819SStefano Zampini ierr = VecScatterBegin(pc_ctx->l2g_lambda,x,pc_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 914674ae819SStefano Zampini ierr = VecScatterEnd(pc_ctx->l2g_lambda,x,pc_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 915674ae819SStefano Zampini ierr = VecSet(pcis->vec2_B,0.0);CHKERRQ(ierr); 916674ae819SStefano Zampini ierr = MatMultTranspose(pc_ctx->B_Ddelta,pc_ctx->lambda_local,pcis->vec2_B);CHKERRQ(ierr); 917ed6c3d69SStefano Zampini /* Application of local Schur complement */ 918ed6c3d69SStefano Zampini ierr = MatMult(pc_ctx->S_j,pcis->vec2_B,pcis->vec1_B);CHKERRQ(ierr); 919edf7251bSStefano Zampini /* Application of B_Ddelta */ 920edf7251bSStefano Zampini ierr = MatMult(pc_ctx->B_Ddelta,pcis->vec1_B,pc_ctx->lambda_local);CHKERRQ(ierr); 921edf7251bSStefano Zampini ierr = VecSet(y,0.0);CHKERRQ(ierr); 922edf7251bSStefano Zampini ierr = VecScatterBegin(pc_ctx->l2g_lambda,pc_ctx->lambda_local,y,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 923edf7251bSStefano Zampini ierr = VecScatterEnd(pc_ctx->l2g_lambda,pc_ctx->lambda_local,y,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 924b519380fSstefano_zampini /* interface pressure preconditioner */ 925b519380fSstefano_zampini if (pc_ctx->kP) { 926b519380fSstefano_zampini const PetscScalar *lx; 927b519380fSstefano_zampini PetscScalar *ly; 928b519380fSstefano_zampini 929b519380fSstefano_zampini /* pressures ordered first in x and y */ 930b519380fSstefano_zampini ierr = VecGetArrayRead(x,&lx);CHKERRQ(ierr); 931b519380fSstefano_zampini ierr = VecGetArray(y,&ly);CHKERRQ(ierr); 932b519380fSstefano_zampini ierr = VecPlaceArray(pc_ctx->xPg,lx);CHKERRQ(ierr); 933b519380fSstefano_zampini ierr = VecPlaceArray(pc_ctx->yPg,ly);CHKERRQ(ierr); 934b519380fSstefano_zampini ierr = KSPSolve(pc_ctx->kP,pc_ctx->xPg,pc_ctx->yPg);CHKERRQ(ierr); 935b519380fSstefano_zampini ierr = VecResetArray(pc_ctx->xPg);CHKERRQ(ierr); 936b519380fSstefano_zampini ierr = VecResetArray(pc_ctx->yPg);CHKERRQ(ierr); 937b519380fSstefano_zampini ierr = VecRestoreArrayRead(x,&lx);CHKERRQ(ierr); 938b519380fSstefano_zampini ierr = VecRestoreArray(y,&ly);CHKERRQ(ierr); 939b519380fSstefano_zampini } 940edf7251bSStefano Zampini PetscFunctionReturn(0); 941edf7251bSStefano Zampini } 942edf7251bSStefano Zampini 943edf7251bSStefano Zampini PetscErrorCode FETIDPPCApplyTranspose(PC fetipc, Vec x, Vec y) 944edf7251bSStefano Zampini { 945edf7251bSStefano Zampini FETIDPPC_ctx pc_ctx; 946edf7251bSStefano Zampini PC_IS *pcis; 947edf7251bSStefano Zampini PetscErrorCode ierr; 948edf7251bSStefano Zampini 949edf7251bSStefano Zampini PetscFunctionBegin; 950302440fdSBarry Smith ierr = PCShellGetContext(fetipc,(void**)&pc_ctx);CHKERRQ(ierr); 951edf7251bSStefano Zampini pcis = (PC_IS*)pc_ctx->pc->data; 952edf7251bSStefano Zampini /* Application of B_Ddelta^T */ 953edf7251bSStefano Zampini ierr = VecScatterBegin(pc_ctx->l2g_lambda,x,pc_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 954edf7251bSStefano Zampini ierr = VecScatterEnd(pc_ctx->l2g_lambda,x,pc_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 955edf7251bSStefano Zampini ierr = VecSet(pcis->vec2_B,0.0);CHKERRQ(ierr); 956edf7251bSStefano Zampini ierr = MatMultTranspose(pc_ctx->B_Ddelta,pc_ctx->lambda_local,pcis->vec2_B);CHKERRQ(ierr); 957ed6c3d69SStefano Zampini /* Application of local Schur complement */ 958ed6c3d69SStefano Zampini ierr = MatMultTranspose(pc_ctx->S_j,pcis->vec2_B,pcis->vec1_B);CHKERRQ(ierr); 959674ae819SStefano Zampini /* Application of B_Ddelta */ 960674ae819SStefano Zampini ierr = MatMult(pc_ctx->B_Ddelta,pcis->vec1_B,pc_ctx->lambda_local);CHKERRQ(ierr); 961674ae819SStefano Zampini ierr = VecSet(y,0.0);CHKERRQ(ierr); 962674ae819SStefano Zampini ierr = VecScatterBegin(pc_ctx->l2g_lambda,pc_ctx->lambda_local,y,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 963674ae819SStefano Zampini ierr = VecScatterEnd(pc_ctx->l2g_lambda,pc_ctx->lambda_local,y,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 964674ae819SStefano Zampini PetscFunctionReturn(0); 965674ae819SStefano Zampini } 966c45b8d2dSstefano_zampini 967c45b8d2dSstefano_zampini PetscErrorCode FETIDPPCView(PC pc, PetscViewer viewer) 968c45b8d2dSstefano_zampini { 969c45b8d2dSstefano_zampini FETIDPPC_ctx pc_ctx; 970c45b8d2dSstefano_zampini PetscBool iascii; 971c45b8d2dSstefano_zampini PetscViewer sviewer; 972c45b8d2dSstefano_zampini PetscErrorCode ierr; 973c45b8d2dSstefano_zampini 974c45b8d2dSstefano_zampini PetscFunctionBegin; 975c45b8d2dSstefano_zampini ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 976c45b8d2dSstefano_zampini if (iascii) { 977c45b8d2dSstefano_zampini PetscMPIInt rank; 97825d06dbeSstefano_zampini PetscBool isschur,isshell; 979c45b8d2dSstefano_zampini 980c45b8d2dSstefano_zampini ierr = PCShellGetContext(pc,(void**)&pc_ctx);CHKERRQ(ierr); 981c45b8d2dSstefano_zampini ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)pc),&rank);CHKERRQ(ierr); 982c45b8d2dSstefano_zampini ierr = PetscObjectTypeCompare((PetscObject)pc_ctx->S_j,MATSCHURCOMPLEMENT,&isschur);CHKERRQ(ierr); 983c45b8d2dSstefano_zampini if (isschur) { 984c45b8d2dSstefano_zampini ierr = PetscViewerASCIIPrintf(viewer," FETI-DP multipliers Dirichlet preconditioner (just from rank 0)\n");CHKERRQ(ierr); 985c45b8d2dSstefano_zampini } else { 986c45b8d2dSstefano_zampini ierr = PetscViewerASCIIPrintf(viewer," FETI-DP multipliers Lumped preconditioner (just from rank 0)\n");CHKERRQ(ierr); 987c45b8d2dSstefano_zampini } 988c45b8d2dSstefano_zampini ierr = PetscViewerGetSubViewer(viewer,PetscObjectComm((PetscObject)pc_ctx->S_j),&sviewer);CHKERRQ(ierr); 989c45b8d2dSstefano_zampini if (!rank) { 990c45b8d2dSstefano_zampini ierr = PetscViewerPushFormat(sviewer,PETSC_VIEWER_ASCII_INFO);CHKERRQ(ierr); 991c45b8d2dSstefano_zampini ierr = PetscViewerASCIIAddTab(sviewer,2);CHKERRQ(ierr); 992c45b8d2dSstefano_zampini ierr = MatView(pc_ctx->S_j,sviewer);CHKERRQ(ierr); 993c45b8d2dSstefano_zampini ierr = PetscViewerASCIISubtractTab(sviewer,2);CHKERRQ(ierr); 994c45b8d2dSstefano_zampini ierr = PetscViewerPopFormat(sviewer);CHKERRQ(ierr); 995c45b8d2dSstefano_zampini } 99625d06dbeSstefano_zampini ierr = PetscObjectTypeCompare((PetscObject)pc_ctx->B_Ddelta,MATSHELL,&isshell);CHKERRQ(ierr); 99725d06dbeSstefano_zampini if (isshell) { 99825d06dbeSstefano_zampini BDdelta_DN ctx; 99925d06dbeSstefano_zampini ierr = PetscViewerASCIIPrintf(viewer," FETI-DP BDdelta: DB^t * (B D^-1 B^t)^-1 for deluxe scaling (just from rank 0)\n");CHKERRQ(ierr); 100025d06dbeSstefano_zampini ierr = MatShellGetContext(pc_ctx->B_Ddelta,&ctx);CHKERRQ(ierr); 100125d06dbeSstefano_zampini if (!rank) { 100225d06dbeSstefano_zampini ierr = PetscViewerASCIIAddTab(sviewer,2);CHKERRQ(ierr); 100325d06dbeSstefano_zampini ierr = KSPView(ctx->kBD,sviewer);CHKERRQ(ierr); 100425d06dbeSstefano_zampini ierr = PetscViewerPushFormat(sviewer,PETSC_VIEWER_ASCII_INFO);CHKERRQ(ierr); 100525d06dbeSstefano_zampini ierr = MatView(ctx->BD,sviewer);CHKERRQ(ierr); 100625d06dbeSstefano_zampini ierr = PetscViewerASCIISubtractTab(sviewer,2);CHKERRQ(ierr); 100725d06dbeSstefano_zampini ierr = PetscViewerPopFormat(sviewer);CHKERRQ(ierr); 100825d06dbeSstefano_zampini } 100925d06dbeSstefano_zampini } 1010c45b8d2dSstefano_zampini ierr = PetscViewerRestoreSubViewer(viewer,PetscObjectComm((PetscObject)pc_ctx->S_j),&sviewer);CHKERRQ(ierr); 1011c45b8d2dSstefano_zampini ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 1012c45b8d2dSstefano_zampini if (pc_ctx->kP) { 1013c45b8d2dSstefano_zampini ierr = PetscViewerASCIIPrintf(viewer," FETI-DP pressure preconditioner\n");CHKERRQ(ierr); 1014c45b8d2dSstefano_zampini ierr = PetscViewerASCIIAddTab(viewer,2);CHKERRQ(ierr); 1015c45b8d2dSstefano_zampini ierr = KSPView(pc_ctx->kP,viewer);CHKERRQ(ierr); 1016c45b8d2dSstefano_zampini ierr = PetscViewerASCIISubtractTab(viewer,2);CHKERRQ(ierr); 1017c45b8d2dSstefano_zampini } 1018c45b8d2dSstefano_zampini } 1019c45b8d2dSstefano_zampini PetscFunctionReturn(0); 1020c45b8d2dSstefano_zampini } 1021