1ab5c6b0cSJed Brown #include <../src/ksp/pc/impls/bddc/bddc.h> 2ab5c6b0cSJed Brown #include <../src/ksp/pc/impls/bddc/bddcprivate.h> 3674ae819SStefano Zampini 4674ae819SStefano Zampini #undef __FUNCT__ 5674ae819SStefano Zampini #define __FUNCT__ "PCBDDCCreateFETIDPMatContext" 6674ae819SStefano Zampini PetscErrorCode PCBDDCCreateFETIDPMatContext(PC pc, FETIDPMat_ctx *fetidpmat_ctx) 7674ae819SStefano Zampini { 8674ae819SStefano Zampini FETIDPMat_ctx newctx; 9674ae819SStefano Zampini PetscErrorCode ierr; 10674ae819SStefano Zampini 11674ae819SStefano Zampini PetscFunctionBegin; 12854ce69bSBarry Smith ierr = PetscNew(&newctx);CHKERRQ(ierr); 13674ae819SStefano Zampini /* increase the reference count for BDDC preconditioner */ 14674ae819SStefano Zampini ierr = PetscObjectReference((PetscObject)pc);CHKERRQ(ierr); 15674ae819SStefano Zampini newctx->pc = pc; 16674ae819SStefano Zampini *fetidpmat_ctx = newctx; 17674ae819SStefano Zampini PetscFunctionReturn(0); 18674ae819SStefano Zampini } 19674ae819SStefano Zampini 20674ae819SStefano Zampini #undef __FUNCT__ 21674ae819SStefano Zampini #define __FUNCT__ "PCBDDCCreateFETIDPPCContext" 22674ae819SStefano Zampini PetscErrorCode PCBDDCCreateFETIDPPCContext(PC pc, FETIDPPC_ctx *fetidppc_ctx) 23674ae819SStefano Zampini { 24674ae819SStefano Zampini FETIDPPC_ctx newctx; 25674ae819SStefano Zampini PetscErrorCode ierr; 26674ae819SStefano Zampini 27674ae819SStefano Zampini PetscFunctionBegin; 28854ce69bSBarry Smith ierr = PetscNew(&newctx);CHKERRQ(ierr); 29674ae819SStefano Zampini /* increase the reference count for BDDC preconditioner */ 30674ae819SStefano Zampini ierr = PetscObjectReference((PetscObject)pc);CHKERRQ(ierr); 31674ae819SStefano Zampini newctx->pc = pc; 32674ae819SStefano Zampini *fetidppc_ctx = newctx; 33674ae819SStefano Zampini PetscFunctionReturn(0); 34674ae819SStefano Zampini } 35674ae819SStefano Zampini 36674ae819SStefano Zampini #undef __FUNCT__ 37674ae819SStefano Zampini #define __FUNCT__ "PCBDDCDestroyFETIDPMat" 38674ae819SStefano Zampini PetscErrorCode PCBDDCDestroyFETIDPMat(Mat A) 39674ae819SStefano Zampini { 40674ae819SStefano Zampini FETIDPMat_ctx mat_ctx; 41674ae819SStefano Zampini PetscErrorCode ierr; 42674ae819SStefano Zampini 43674ae819SStefano Zampini PetscFunctionBegin; 44674ae819SStefano Zampini ierr = MatShellGetContext(A,(void**)&mat_ctx);CHKERRQ(ierr); 45674ae819SStefano Zampini ierr = VecDestroy(&mat_ctx->lambda_local);CHKERRQ(ierr); 46674ae819SStefano Zampini ierr = VecDestroy(&mat_ctx->temp_solution_D);CHKERRQ(ierr); 47674ae819SStefano Zampini ierr = VecDestroy(&mat_ctx->temp_solution_B);CHKERRQ(ierr); 48674ae819SStefano Zampini ierr = MatDestroy(&mat_ctx->B_delta);CHKERRQ(ierr); 49674ae819SStefano Zampini ierr = MatDestroy(&mat_ctx->B_Ddelta);CHKERRQ(ierr); 50457d4a33Sstefano_zampini ierr = MatDestroy(&mat_ctx->B_BB);CHKERRQ(ierr); 51457d4a33Sstefano_zampini ierr = MatDestroy(&mat_ctx->B_BI);CHKERRQ(ierr); 52457d4a33Sstefano_zampini ierr = MatDestroy(&mat_ctx->Bt_BB);CHKERRQ(ierr); 53457d4a33Sstefano_zampini ierr = MatDestroy(&mat_ctx->Bt_BI);CHKERRQ(ierr); 54457d4a33Sstefano_zampini ierr = MatDestroy(&mat_ctx->C);CHKERRQ(ierr); 55*6cc1294bSstefano_zampini ierr = VecDestroy(&mat_ctx->rhs_flip);CHKERRQ(ierr); 56457d4a33Sstefano_zampini ierr = VecDestroy(&mat_ctx->vP);CHKERRQ(ierr); 57457d4a33Sstefano_zampini ierr = VecDestroy(&mat_ctx->xPg);CHKERRQ(ierr); 58457d4a33Sstefano_zampini ierr = VecDestroy(&mat_ctx->yPg);CHKERRQ(ierr); 59674ae819SStefano Zampini ierr = VecScatterDestroy(&mat_ctx->l2g_lambda);CHKERRQ(ierr); 60457d4a33Sstefano_zampini ierr = VecScatterDestroy(&mat_ctx->l2g_p);CHKERRQ(ierr); 619cc7774eSstefano_zampini ierr = VecScatterDestroy(&mat_ctx->g2g_p);CHKERRQ(ierr); 62674ae819SStefano Zampini ierr = PCDestroy(&mat_ctx->pc);CHKERRQ(ierr); /* decrease PCBDDC reference count */ 63674ae819SStefano Zampini ierr = PetscFree(mat_ctx);CHKERRQ(ierr); 64674ae819SStefano Zampini PetscFunctionReturn(0); 65674ae819SStefano Zampini } 66674ae819SStefano Zampini 67674ae819SStefano Zampini #undef __FUNCT__ 68674ae819SStefano Zampini #define __FUNCT__ "PCBDDCDestroyFETIDPPC" 69674ae819SStefano Zampini PetscErrorCode PCBDDCDestroyFETIDPPC(PC pc) 70674ae819SStefano Zampini { 71674ae819SStefano Zampini FETIDPPC_ctx pc_ctx; 72674ae819SStefano Zampini PetscErrorCode ierr; 73674ae819SStefano Zampini 74674ae819SStefano Zampini PetscFunctionBegin; 75674ae819SStefano Zampini ierr = PCShellGetContext(pc,(void**)&pc_ctx);CHKERRQ(ierr); 76674ae819SStefano Zampini ierr = VecDestroy(&pc_ctx->lambda_local);CHKERRQ(ierr); 77674ae819SStefano Zampini ierr = MatDestroy(&pc_ctx->B_Ddelta);CHKERRQ(ierr); 78674ae819SStefano Zampini ierr = VecScatterDestroy(&pc_ctx->l2g_lambda);CHKERRQ(ierr); 7968668abeSStefano Zampini ierr = MatDestroy(&pc_ctx->S_j);CHKERRQ(ierr); 80674ae819SStefano Zampini ierr = PCDestroy(&pc_ctx->pc);CHKERRQ(ierr); /* decrease PCBDDC reference count */ 81457d4a33Sstefano_zampini ierr = KSPDestroy(&pc_ctx->kP);CHKERRQ(ierr); 82457d4a33Sstefano_zampini ierr = VecDestroy(&pc_ctx->xPg);CHKERRQ(ierr); 83457d4a33Sstefano_zampini ierr = VecDestroy(&pc_ctx->yPg);CHKERRQ(ierr); 84674ae819SStefano Zampini ierr = PetscFree(pc_ctx);CHKERRQ(ierr); 85674ae819SStefano Zampini PetscFunctionReturn(0); 86674ae819SStefano Zampini } 87674ae819SStefano Zampini 88674ae819SStefano Zampini #undef __FUNCT__ 89674ae819SStefano Zampini #define __FUNCT__ "PCBDDCSetupFETIDPMatContext" 90674ae819SStefano Zampini PetscErrorCode PCBDDCSetupFETIDPMatContext(FETIDPMat_ctx fetidpmat_ctx ) 91674ae819SStefano Zampini { 92674ae819SStefano Zampini PetscErrorCode ierr; 93674ae819SStefano Zampini PC_IS *pcis=(PC_IS*)fetidpmat_ctx->pc->data; 94674ae819SStefano Zampini PC_BDDC *pcbddc=(PC_BDDC*)fetidpmat_ctx->pc->data; 95674ae819SStefano Zampini PCBDDCGraph mat_graph=pcbddc->mat_graph; 96674ae819SStefano Zampini Mat_IS *matis = (Mat_IS*)fetidpmat_ctx->pc->pmat->data; 97674ae819SStefano Zampini MPI_Comm comm; 98674ae819SStefano Zampini Mat ScalingMat; 99457d4a33Sstefano_zampini Vec fetidp_global; 100674ae819SStefano Zampini IS IS_l2g_lambda; 101dc456d91SStefano Zampini IS subset,subset_mult,subset_n; 102674ae819SStefano Zampini PetscBool skip_node,fully_redundant; 103674ae819SStefano Zampini PetscInt i,j,k,s,n_boundary_dofs,n_global_lambda,n_vertices,partial_sum; 104dc456d91SStefano Zampini PetscInt cum,n_local_lambda,n_lambda_for_dof,dual_size,n_neg_values,n_pos_values; 10576ec1555SBarry Smith PetscMPIInt rank,size,buf_size,neigh; 106674ae819SStefano Zampini PetscScalar scalar_value; 107674ae819SStefano Zampini PetscInt *vertex_indices; 108dc456d91SStefano Zampini PetscInt *dual_dofs_boundary_indices,*aux_local_numbering_1; 109dc456d91SStefano Zampini const PetscInt *aux_global_numbering; 110674ae819SStefano Zampini PetscInt *aux_sums,*cols_B_delta,*l2g_indices; 111674ae819SStefano Zampini PetscScalar *array,*scaling_factors,*vals_B_delta; 112674ae819SStefano Zampini PetscInt *aux_local_numbering_2; 113457d4a33Sstefano_zampini PetscLayout llay; 114674ae819SStefano Zampini /* For communication of scaling factors */ 115674ae819SStefano Zampini PetscInt *ptrs_buffer,neigh_position; 116674ae819SStefano Zampini PetscScalar **all_factors,*send_buffer,*recv_buffer; 117674ae819SStefano Zampini MPI_Request *send_reqs,*recv_reqs; 118674ae819SStefano Zampini /* tests */ 119674ae819SStefano Zampini PetscBool test_fetidp; 120674ae819SStefano Zampini PetscViewer viewer; 121457d4a33Sstefano_zampini /* saddlepoint */ 122457d4a33Sstefano_zampini ISLocalToGlobalMapping l2gmap_p; 123457d4a33Sstefano_zampini PetscLayout play; 124457d4a33Sstefano_zampini IS gP,pP; 125457d4a33Sstefano_zampini PetscInt nPl,nPg,nPgl; 126674ae819SStefano Zampini 127674ae819SStefano Zampini PetscFunctionBegin; 128674ae819SStefano Zampini ierr = PetscObjectGetComm((PetscObject)(fetidpmat_ctx->pc),&comm);CHKERRQ(ierr); 129674ae819SStefano Zampini ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 13076ec1555SBarry Smith ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 131674ae819SStefano Zampini 132457d4a33Sstefano_zampini /* saddlepoint */ 133457d4a33Sstefano_zampini nPl = 0; 134457d4a33Sstefano_zampini nPg = 0; 135457d4a33Sstefano_zampini nPgl = 0; 136457d4a33Sstefano_zampini gP = NULL; 137457d4a33Sstefano_zampini pP = NULL; 138457d4a33Sstefano_zampini l2gmap_p = NULL; 139457d4a33Sstefano_zampini play = NULL; 140457d4a33Sstefano_zampini ierr = PetscObjectQuery((PetscObject)fetidpmat_ctx->pc,"__KSPFETIDP_pP",(PetscObject*)&pP);CHKERRQ(ierr); 141022d8d2bSstefano_zampini if (pP) { /* saddle point */ 142457d4a33Sstefano_zampini /* subdomain pressures in global numbering */ 143457d4a33Sstefano_zampini ierr = PetscObjectQuery((PetscObject)fetidpmat_ctx->pc,"__KSPFETIDP_gP",(PetscObject*)&gP);CHKERRQ(ierr); 144457d4a33Sstefano_zampini if (!gP) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"gP not present"); 145457d4a33Sstefano_zampini ierr = ISGetLocalSize(gP,&nPl);CHKERRQ(ierr); 146457d4a33Sstefano_zampini ierr = VecCreate(PETSC_COMM_SELF,&fetidpmat_ctx->vP);CHKERRQ(ierr); 147457d4a33Sstefano_zampini ierr = VecSetSizes(fetidpmat_ctx->vP,nPl,nPl);CHKERRQ(ierr); 148457d4a33Sstefano_zampini ierr = VecSetType(fetidpmat_ctx->vP,VECSTANDARD);CHKERRQ(ierr); 149457d4a33Sstefano_zampini ierr = VecSetUp(fetidpmat_ctx->vP);CHKERRQ(ierr); 150457d4a33Sstefano_zampini 151457d4a33Sstefano_zampini /* interface pressure matrix */ 152457d4a33Sstefano_zampini ierr = PetscObjectQuery((PetscObject)fetidpmat_ctx->pc,"__KSPFETIDP_C",(PetscObject*)&fetidpmat_ctx->C);CHKERRQ(ierr); 153457d4a33Sstefano_zampini if (!fetidpmat_ctx->C) { /* null pressure block, compute layout and global numbering for interface pressures */ 154457d4a33Sstefano_zampini IS Pg; 155457d4a33Sstefano_zampini 156457d4a33Sstefano_zampini ierr = ISRenumber(gP,NULL,&nPg,&Pg);CHKERRQ(ierr); 157457d4a33Sstefano_zampini ierr = ISLocalToGlobalMappingCreateIS(Pg,&l2gmap_p);CHKERRQ(ierr); 158457d4a33Sstefano_zampini ierr = ISDestroy(&Pg);CHKERRQ(ierr); 159457d4a33Sstefano_zampini ierr = PetscLayoutCreate(comm,&play);CHKERRQ(ierr); 160457d4a33Sstefano_zampini ierr = PetscLayoutSetBlockSize(play,1);CHKERRQ(ierr); 161457d4a33Sstefano_zampini ierr = PetscLayoutSetSize(play,nPg);CHKERRQ(ierr); 162457d4a33Sstefano_zampini ierr = ISGetLocalSize(pP,&nPgl);CHKERRQ(ierr); 163457d4a33Sstefano_zampini ierr = PetscLayoutSetLocalSize(play,nPgl);CHKERRQ(ierr); 164457d4a33Sstefano_zampini ierr = PetscLayoutSetUp(play);CHKERRQ(ierr); 165457d4a33Sstefano_zampini } else { 166457d4a33Sstefano_zampini ierr = PetscObjectReference((PetscObject)fetidpmat_ctx->C);CHKERRQ(ierr); 167457d4a33Sstefano_zampini ierr = MatGetLocalToGlobalMapping(fetidpmat_ctx->C,&l2gmap_p,NULL);CHKERRQ(ierr); 168457d4a33Sstefano_zampini ierr = PetscObjectReference((PetscObject)l2gmap_p);CHKERRQ(ierr); 169457d4a33Sstefano_zampini ierr = MatGetSize(fetidpmat_ctx->C,&nPg,NULL);CHKERRQ(ierr); 170457d4a33Sstefano_zampini ierr = MatGetLocalSize(fetidpmat_ctx->C,NULL,&nPgl);CHKERRQ(ierr); 171457d4a33Sstefano_zampini ierr = MatGetLayouts(fetidpmat_ctx->C,NULL,&llay);CHKERRQ(ierr); 172457d4a33Sstefano_zampini ierr = PetscLayoutReference(llay,&play);CHKERRQ(ierr); 173457d4a33Sstefano_zampini } 174457d4a33Sstefano_zampini ierr = VecCreateMPIWithArray(comm,1,nPgl,nPg,NULL,&fetidpmat_ctx->xPg);CHKERRQ(ierr); 175457d4a33Sstefano_zampini ierr = VecCreateMPIWithArray(comm,1,nPgl,nPg,NULL,&fetidpmat_ctx->yPg);CHKERRQ(ierr); 176457d4a33Sstefano_zampini 177457d4a33Sstefano_zampini /* import matrices for interface pressures coupling */ 178457d4a33Sstefano_zampini ierr = PetscObjectQuery((PetscObject)fetidpmat_ctx->pc,"__KSPFETIDP_B_BI",(PetscObject*)&fetidpmat_ctx->B_BI);CHKERRQ(ierr); 179457d4a33Sstefano_zampini if (!fetidpmat_ctx->B_BI) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"B_BI not present"); 180457d4a33Sstefano_zampini ierr = PetscObjectReference((PetscObject)fetidpmat_ctx->B_BI);CHKERRQ(ierr); 181457d4a33Sstefano_zampini 182457d4a33Sstefano_zampini ierr = PetscObjectQuery((PetscObject)fetidpmat_ctx->pc,"__KSPFETIDP_B_BB",(PetscObject*)&fetidpmat_ctx->B_BB);CHKERRQ(ierr); 183457d4a33Sstefano_zampini if (!fetidpmat_ctx->B_BB) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"B_BB not present"); 184457d4a33Sstefano_zampini ierr = PetscObjectReference((PetscObject)fetidpmat_ctx->B_BB);CHKERRQ(ierr); 185457d4a33Sstefano_zampini 186457d4a33Sstefano_zampini ierr = PetscObjectQuery((PetscObject)fetidpmat_ctx->pc,"__KSPFETIDP_Bt_BI",(PetscObject*)&fetidpmat_ctx->Bt_BI);CHKERRQ(ierr); 187457d4a33Sstefano_zampini if (!fetidpmat_ctx->Bt_BI) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Bt_BI not present"); 188457d4a33Sstefano_zampini ierr = PetscObjectReference((PetscObject)fetidpmat_ctx->Bt_BI);CHKERRQ(ierr); 189457d4a33Sstefano_zampini 190457d4a33Sstefano_zampini ierr = PetscObjectQuery((PetscObject)fetidpmat_ctx->pc,"__KSPFETIDP_Bt_BB",(PetscObject*)&fetidpmat_ctx->Bt_BB);CHKERRQ(ierr); 191457d4a33Sstefano_zampini if (!fetidpmat_ctx->Bt_BB) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Bt_BB not present"); 192457d4a33Sstefano_zampini ierr = PetscObjectReference((PetscObject)fetidpmat_ctx->Bt_BB);CHKERRQ(ierr); 193*6cc1294bSstefano_zampini 194*6cc1294bSstefano_zampini ierr = PetscObjectQuery((PetscObject)fetidpmat_ctx->pc,"__KSPFETIDP_flip" ,(PetscObject*)&fetidpmat_ctx->rhs_flip);CHKERRQ(ierr); 195*6cc1294bSstefano_zampini if (fetidpmat_ctx->rhs_flip) { 196*6cc1294bSstefano_zampini ierr = PetscObjectReference((PetscObject)fetidpmat_ctx->rhs_flip);CHKERRQ(ierr); 197*6cc1294bSstefano_zampini } 198457d4a33Sstefano_zampini } 199457d4a33Sstefano_zampini 200674ae819SStefano Zampini /* Default type of lagrange multipliers is non-redundant */ 201329cd39dSStefano Zampini fully_redundant = fetidpmat_ctx->fully_redundant; 202674ae819SStefano Zampini 203674ae819SStefano Zampini /* Evaluate local and global number of lagrange multipliers */ 204674ae819SStefano Zampini ierr = VecSet(pcis->vec1_N,0.0);CHKERRQ(ierr); 205674ae819SStefano Zampini n_local_lambda = 0; 206674ae819SStefano Zampini partial_sum = 0; 207674ae819SStefano Zampini n_boundary_dofs = 0; 208674ae819SStefano Zampini s = 0; 209674ae819SStefano Zampini /* Get Vertices used to define the BDDC */ 210b371cd4fSStefano Zampini n_vertices = pcbddc->n_vertices; 2110e6343abSStefano Zampini vertex_indices = pcbddc->local_primal_ref_node; 212674ae819SStefano Zampini dual_size = pcis->n_B-n_vertices; 213785e854fSJed Brown ierr = PetscMalloc1(dual_size,&dual_dofs_boundary_indices);CHKERRQ(ierr); 214785e854fSJed Brown ierr = PetscMalloc1(dual_size,&aux_local_numbering_1);CHKERRQ(ierr); 215785e854fSJed Brown ierr = PetscMalloc1(dual_size,&aux_local_numbering_2);CHKERRQ(ierr); 216674ae819SStefano Zampini 217674ae819SStefano Zampini ierr = VecGetArray(pcis->vec1_N,&array);CHKERRQ(ierr); 218674ae819SStefano Zampini for (i=0;i<pcis->n;i++){ 219674ae819SStefano Zampini j = mat_graph->count[i]; /* RECALL: mat_graph->count[i] does not count myself */ 2202d456bbaSstefano_zampini if (j > 0) n_boundary_dofs++; 221674ae819SStefano Zampini skip_node = PETSC_FALSE; 222674ae819SStefano Zampini if (s < n_vertices && vertex_indices[s] == i) { /* it works for a sorted set of vertices */ 223674ae819SStefano Zampini skip_node = PETSC_TRUE; 224674ae819SStefano Zampini s++; 225674ae819SStefano Zampini } 2262d456bbaSstefano_zampini if (j < 1) skip_node = PETSC_TRUE; 2272d456bbaSstefano_zampini if (mat_graph->special_dof[i] == PCBDDCGRAPH_DIRICHLET_MARK) skip_node = PETSC_TRUE; 228674ae819SStefano Zampini if (!skip_node) { 229674ae819SStefano Zampini if (fully_redundant) { 230674ae819SStefano Zampini /* fully redundant set of lagrange multipliers */ 231674ae819SStefano Zampini n_lambda_for_dof = (j*(j+1))/2; 232674ae819SStefano Zampini } else { 233674ae819SStefano Zampini n_lambda_for_dof = j; 234674ae819SStefano Zampini } 235674ae819SStefano Zampini n_local_lambda += j; 236674ae819SStefano Zampini /* needed to evaluate global number of lagrange multipliers */ 237674ae819SStefano Zampini array[i]=(1.0*n_lambda_for_dof)/(j+1.0); /* already scaled for the next global sum */ 238674ae819SStefano Zampini /* store some data needed */ 239674ae819SStefano Zampini dual_dofs_boundary_indices[partial_sum] = n_boundary_dofs-1; 240674ae819SStefano Zampini aux_local_numbering_1[partial_sum] = i; 241674ae819SStefano Zampini aux_local_numbering_2[partial_sum] = n_lambda_for_dof; 242674ae819SStefano Zampini partial_sum++; 243674ae819SStefano Zampini } 244674ae819SStefano Zampini } 245674ae819SStefano Zampini ierr = VecRestoreArray(pcis->vec1_N,&array);CHKERRQ(ierr); 2462d456bbaSstefano_zampini dual_size = partial_sum; 247674ae819SStefano Zampini 248674ae819SStefano Zampini /* compute global ordering of lagrange multipliers and associate l2g map */ 249dc456d91SStefano Zampini ierr = ISCreateGeneral(comm,partial_sum,aux_local_numbering_1,PETSC_COPY_VALUES,&subset_n);CHKERRQ(ierr); 2503bbff08aSStefano Zampini ierr = ISLocalToGlobalMappingApplyIS(pcis->mapping,subset_n,&subset);CHKERRQ(ierr); 251dc456d91SStefano Zampini ierr = ISDestroy(&subset_n);CHKERRQ(ierr); 252dc456d91SStefano Zampini ierr = ISCreateGeneral(comm,partial_sum,aux_local_numbering_2,PETSC_OWN_POINTER,&subset_mult);CHKERRQ(ierr); 2533d996552SStefano Zampini ierr = ISRenumber(subset,subset_mult,&fetidpmat_ctx->n_lambda,&subset_n);CHKERRQ(ierr); 254dc456d91SStefano Zampini ierr = ISDestroy(&subset);CHKERRQ(ierr); 2553d996552SStefano Zampini 2564fe826edSStefano Zampini #if defined(PETSC_USE_DEBUG) 2574fe826edSStefano Zampini ierr = VecSet(pcis->vec1_global,0.0);CHKERRQ(ierr); 2584fe826edSStefano Zampini ierr = VecScatterBegin(matis->rctx,pcis->vec1_N,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 2594fe826edSStefano Zampini ierr = VecScatterEnd(matis->rctx,pcis->vec1_N,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 2604fe826edSStefano Zampini ierr = VecSum(pcis->vec1_global,&scalar_value);CHKERRQ(ierr); 2614fe826edSStefano Zampini i = (PetscInt)PetscRealPart(scalar_value); 2626c4ed002SBarry 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); 2634fe826edSStefano Zampini #endif 264674ae819SStefano Zampini 265674ae819SStefano Zampini /* init data for scaling factors exchange */ 266674ae819SStefano Zampini partial_sum = 0; 267785e854fSJed Brown ierr = PetscMalloc1(pcis->n_neigh,&ptrs_buffer);CHKERRQ(ierr); 2684b2aedd3SStefano Zampini ierr = PetscMalloc1(PetscMax(pcis->n_neigh-1,0),&send_reqs);CHKERRQ(ierr); 2694b2aedd3SStefano Zampini ierr = PetscMalloc1(PetscMax(pcis->n_neigh-1,0),&recv_reqs);CHKERRQ(ierr); 270785e854fSJed Brown ierr = PetscMalloc1(pcis->n,&all_factors);CHKERRQ(ierr); 2714b2aedd3SStefano Zampini if (pcis->n_neigh > 0) ptrs_buffer[0]=0; 272674ae819SStefano Zampini for (i=1;i<pcis->n_neigh;i++) { 273674ae819SStefano Zampini partial_sum += pcis->n_shared[i]; 274674ae819SStefano Zampini ptrs_buffer[i] = ptrs_buffer[i-1]+pcis->n_shared[i]; 275674ae819SStefano Zampini } 276785e854fSJed Brown ierr = PetscMalloc1(partial_sum,&send_buffer);CHKERRQ(ierr); 277785e854fSJed Brown ierr = PetscMalloc1(partial_sum,&recv_buffer);CHKERRQ(ierr); 278785e854fSJed Brown ierr = PetscMalloc1(partial_sum,&all_factors[0]);CHKERRQ(ierr); 279674ae819SStefano Zampini for (i=0;i<pcis->n-1;i++) { 280674ae819SStefano Zampini j = mat_graph->count[i]; 281674ae819SStefano Zampini all_factors[i+1]=all_factors[i]+j; 282674ae819SStefano Zampini } 283674ae819SStefano Zampini /* scatter B scaling to N vec */ 284674ae819SStefano Zampini ierr = VecScatterBegin(pcis->N_to_B,pcis->D,pcis->vec1_N,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 285674ae819SStefano Zampini ierr = VecScatterEnd(pcis->N_to_B,pcis->D,pcis->vec1_N,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 286674ae819SStefano Zampini /* communications */ 2872b095fd8SStefano Zampini ierr = VecGetArrayRead(pcis->vec1_N,(const PetscScalar**)&array);CHKERRQ(ierr); 288674ae819SStefano Zampini for (i=1;i<pcis->n_neigh;i++) { 289674ae819SStefano Zampini for (j=0;j<pcis->n_shared[i];j++) { 290674ae819SStefano Zampini send_buffer[ptrs_buffer[i-1]+j]=array[pcis->shared[i][j]]; 291674ae819SStefano Zampini } 292674ae819SStefano Zampini ierr = PetscMPIIntCast(ptrs_buffer[i]-ptrs_buffer[i-1],&buf_size);CHKERRQ(ierr); 293674ae819SStefano Zampini ierr = PetscMPIIntCast(pcis->neigh[i],&neigh);CHKERRQ(ierr); 294674ae819SStefano Zampini ierr = MPI_Isend(&send_buffer[ptrs_buffer[i-1]],buf_size,MPIU_SCALAR,neigh,0,comm,&send_reqs[i-1]);CHKERRQ(ierr); 295674ae819SStefano Zampini ierr = MPI_Irecv(&recv_buffer[ptrs_buffer[i-1]],buf_size,MPIU_SCALAR,neigh,0,comm,&recv_reqs[i-1]);CHKERRQ(ierr); 296674ae819SStefano Zampini } 2972b095fd8SStefano Zampini ierr = VecRestoreArrayRead(pcis->vec1_N,(const PetscScalar**)&array);CHKERRQ(ierr); 2984b2aedd3SStefano Zampini if (pcis->n_neigh > 0) { 2994b2aedd3SStefano Zampini ierr = MPI_Waitall(pcis->n_neigh-1,recv_reqs,MPI_STATUSES_IGNORE);CHKERRQ(ierr); 3004b2aedd3SStefano Zampini } 301674ae819SStefano Zampini /* put values in correct places */ 302674ae819SStefano Zampini for (i=1;i<pcis->n_neigh;i++) { 303674ae819SStefano Zampini for (j=0;j<pcis->n_shared[i];j++) { 304674ae819SStefano Zampini k = pcis->shared[i][j]; 305674ae819SStefano Zampini neigh_position = 0; 306674ae819SStefano Zampini while(mat_graph->neighbours_set[k][neigh_position] != pcis->neigh[i]) {neigh_position++;} 307674ae819SStefano Zampini all_factors[k][neigh_position]=recv_buffer[ptrs_buffer[i-1]+j]; 308674ae819SStefano Zampini } 309674ae819SStefano Zampini } 3104b2aedd3SStefano Zampini if (pcis->n_neigh > 0) { 3114b2aedd3SStefano Zampini ierr = MPI_Waitall(pcis->n_neigh-1,send_reqs,MPI_STATUSES_IGNORE);CHKERRQ(ierr); 3124b2aedd3SStefano Zampini } 313674ae819SStefano Zampini ierr = PetscFree(send_reqs);CHKERRQ(ierr); 314674ae819SStefano Zampini ierr = PetscFree(recv_reqs);CHKERRQ(ierr); 315674ae819SStefano Zampini ierr = PetscFree(send_buffer);CHKERRQ(ierr); 316674ae819SStefano Zampini ierr = PetscFree(recv_buffer);CHKERRQ(ierr); 317674ae819SStefano Zampini ierr = PetscFree(ptrs_buffer);CHKERRQ(ierr); 318674ae819SStefano Zampini 319674ae819SStefano Zampini /* Compute B and B_delta (local actions) */ 320785e854fSJed Brown ierr = PetscMalloc1(pcis->n_neigh,&aux_sums);CHKERRQ(ierr); 321785e854fSJed Brown ierr = PetscMalloc1(n_local_lambda,&l2g_indices);CHKERRQ(ierr); 322785e854fSJed Brown ierr = PetscMalloc1(n_local_lambda,&vals_B_delta);CHKERRQ(ierr); 323785e854fSJed Brown ierr = PetscMalloc1(n_local_lambda,&cols_B_delta);CHKERRQ(ierr); 324785e854fSJed Brown ierr = PetscMalloc1(n_local_lambda,&scaling_factors);CHKERRQ(ierr); 325dc456d91SStefano Zampini ierr = ISGetIndices(subset_n,&aux_global_numbering);CHKERRQ(ierr); 326674ae819SStefano Zampini partial_sum=0; 327dc456d91SStefano Zampini cum = 0; 328674ae819SStefano Zampini for (i=0;i<dual_size;i++) { 329dc456d91SStefano Zampini n_global_lambda = aux_global_numbering[cum]; 330674ae819SStefano Zampini j = mat_graph->count[aux_local_numbering_1[i]]; 331674ae819SStefano Zampini aux_sums[0]=0; 332674ae819SStefano Zampini for (s=1;s<j;s++) { 333674ae819SStefano Zampini aux_sums[s]=aux_sums[s-1]+j-s+1; 334674ae819SStefano Zampini } 335674ae819SStefano Zampini array = all_factors[aux_local_numbering_1[i]]; 336674ae819SStefano Zampini n_neg_values = 0; 3372a7da448SStefano Zampini while(n_neg_values < j && mat_graph->neighbours_set[aux_local_numbering_1[i]][n_neg_values] < rank) {n_neg_values++;} 338674ae819SStefano Zampini n_pos_values = j - n_neg_values; 339674ae819SStefano Zampini if (fully_redundant) { 340674ae819SStefano Zampini for (s=0;s<n_neg_values;s++) { 341674ae819SStefano Zampini l2g_indices [partial_sum+s]=aux_sums[s]+n_neg_values-s-1+n_global_lambda; 342674ae819SStefano Zampini cols_B_delta [partial_sum+s]=dual_dofs_boundary_indices[i]; 343674ae819SStefano Zampini vals_B_delta [partial_sum+s]=-1.0; 344674ae819SStefano Zampini scaling_factors[partial_sum+s]=array[s]; 345674ae819SStefano Zampini } 346674ae819SStefano Zampini for (s=0;s<n_pos_values;s++) { 347674ae819SStefano Zampini l2g_indices [partial_sum+s+n_neg_values]=aux_sums[n_neg_values]+s+n_global_lambda; 348674ae819SStefano Zampini cols_B_delta [partial_sum+s+n_neg_values]=dual_dofs_boundary_indices[i]; 349674ae819SStefano Zampini vals_B_delta [partial_sum+s+n_neg_values]=1.0; 350674ae819SStefano Zampini scaling_factors[partial_sum+s+n_neg_values]=array[s+n_neg_values]; 351674ae819SStefano Zampini } 352674ae819SStefano Zampini partial_sum += j; 353674ae819SStefano Zampini } else { 354674ae819SStefano Zampini /* l2g_indices and default cols and vals of B_delta */ 355674ae819SStefano Zampini for (s=0;s<j;s++) { 356674ae819SStefano Zampini l2g_indices [partial_sum+s]=n_global_lambda+s; 357674ae819SStefano Zampini cols_B_delta [partial_sum+s]=dual_dofs_boundary_indices[i]; 358674ae819SStefano Zampini vals_B_delta [partial_sum+s]=0.0; 359674ae819SStefano Zampini } 360674ae819SStefano Zampini /* B_delta */ 361674ae819SStefano Zampini if ( n_neg_values > 0 ) { /* there's a rank next to me to the left */ 362674ae819SStefano Zampini vals_B_delta [partial_sum+n_neg_values-1]=-1.0; 363674ae819SStefano Zampini } 364674ae819SStefano Zampini if ( n_neg_values < j ) { /* there's a rank next to me to the right */ 365674ae819SStefano Zampini vals_B_delta [partial_sum+n_neg_values]=1.0; 366674ae819SStefano Zampini } 367674ae819SStefano Zampini /* scaling as in Klawonn-Widlund 1999*/ 368674ae819SStefano Zampini for (s=0;s<n_neg_values;s++) { 369674ae819SStefano Zampini scalar_value = 0.0; 370674ae819SStefano Zampini for (k=0;k<s+1;k++) { 371674ae819SStefano Zampini scalar_value += array[k]; 372674ae819SStefano Zampini } 373674ae819SStefano Zampini scaling_factors[partial_sum+s] = -scalar_value; 374674ae819SStefano Zampini } 375674ae819SStefano Zampini for (s=0;s<n_pos_values;s++) { 376674ae819SStefano Zampini scalar_value = 0.0; 377674ae819SStefano Zampini for (k=s+n_neg_values;k<j;k++) { 378674ae819SStefano Zampini scalar_value += array[k]; 379674ae819SStefano Zampini } 380674ae819SStefano Zampini scaling_factors[partial_sum+s+n_neg_values] = scalar_value; 381674ae819SStefano Zampini } 382674ae819SStefano Zampini partial_sum += j; 383674ae819SStefano Zampini } 384dc456d91SStefano Zampini cum += aux_local_numbering_2[i]; 385674ae819SStefano Zampini } 386dc456d91SStefano Zampini ierr = ISRestoreIndices(subset_n,&aux_global_numbering);CHKERRQ(ierr); 387dc456d91SStefano Zampini ierr = ISDestroy(&subset_mult);CHKERRQ(ierr); 388dc456d91SStefano Zampini ierr = ISDestroy(&subset_n);CHKERRQ(ierr); 389674ae819SStefano Zampini ierr = PetscFree(aux_sums);CHKERRQ(ierr); 390674ae819SStefano Zampini ierr = PetscFree(aux_local_numbering_1);CHKERRQ(ierr); 391674ae819SStefano Zampini ierr = PetscFree(dual_dofs_boundary_indices);CHKERRQ(ierr); 392674ae819SStefano Zampini ierr = PetscFree(all_factors[0]);CHKERRQ(ierr); 393674ae819SStefano Zampini ierr = PetscFree(all_factors);CHKERRQ(ierr); 394674ae819SStefano Zampini 395674ae819SStefano Zampini /* Create local part of B_delta */ 396302440fdSBarry Smith ierr = MatCreate(PETSC_COMM_SELF,&fetidpmat_ctx->B_delta);CHKERRQ(ierr); 397674ae819SStefano Zampini ierr = MatSetSizes(fetidpmat_ctx->B_delta,n_local_lambda,pcis->n_B,n_local_lambda,pcis->n_B);CHKERRQ(ierr); 398674ae819SStefano Zampini ierr = MatSetType(fetidpmat_ctx->B_delta,MATSEQAIJ);CHKERRQ(ierr); 399674ae819SStefano Zampini ierr = MatSeqAIJSetPreallocation(fetidpmat_ctx->B_delta,1,NULL);CHKERRQ(ierr); 400674ae819SStefano Zampini ierr = MatSetOption(fetidpmat_ctx->B_delta,MAT_IGNORE_ZERO_ENTRIES,PETSC_TRUE);CHKERRQ(ierr); 401674ae819SStefano Zampini for (i=0;i<n_local_lambda;i++) { 402674ae819SStefano Zampini ierr = MatSetValue(fetidpmat_ctx->B_delta,i,cols_B_delta[i],vals_B_delta[i],INSERT_VALUES);CHKERRQ(ierr); 403674ae819SStefano Zampini } 404674ae819SStefano Zampini ierr = PetscFree(vals_B_delta);CHKERRQ(ierr); 405674ae819SStefano Zampini ierr = MatAssemblyBegin(fetidpmat_ctx->B_delta,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 406674ae819SStefano Zampini ierr = MatAssemblyEnd (fetidpmat_ctx->B_delta,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 407674ae819SStefano Zampini 408674ae819SStefano Zampini if (fully_redundant) { 409302440fdSBarry Smith ierr = MatCreate(PETSC_COMM_SELF,&ScalingMat);CHKERRQ(ierr); 410674ae819SStefano Zampini ierr = MatSetSizes(ScalingMat,n_local_lambda,n_local_lambda,n_local_lambda,n_local_lambda);CHKERRQ(ierr); 411674ae819SStefano Zampini ierr = MatSetType(ScalingMat,MATSEQAIJ);CHKERRQ(ierr); 412674ae819SStefano Zampini ierr = MatSeqAIJSetPreallocation(ScalingMat,1,NULL);CHKERRQ(ierr); 413674ae819SStefano Zampini for (i=0;i<n_local_lambda;i++) { 414674ae819SStefano Zampini ierr = MatSetValue(ScalingMat,i,i,scaling_factors[i],INSERT_VALUES);CHKERRQ(ierr); 415674ae819SStefano Zampini } 416674ae819SStefano Zampini ierr = MatAssemblyBegin(ScalingMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 417674ae819SStefano Zampini ierr = MatAssemblyEnd (ScalingMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 418674ae819SStefano Zampini ierr = MatMatMult(ScalingMat,fetidpmat_ctx->B_delta,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&fetidpmat_ctx->B_Ddelta);CHKERRQ(ierr); 419674ae819SStefano Zampini ierr = MatDestroy(&ScalingMat);CHKERRQ(ierr); 420674ae819SStefano Zampini } else { 421302440fdSBarry Smith ierr = MatCreate(PETSC_COMM_SELF,&fetidpmat_ctx->B_Ddelta);CHKERRQ(ierr); 422674ae819SStefano Zampini ierr = MatSetSizes(fetidpmat_ctx->B_Ddelta,n_local_lambda,pcis->n_B,n_local_lambda,pcis->n_B);CHKERRQ(ierr); 423674ae819SStefano Zampini ierr = MatSetType(fetidpmat_ctx->B_Ddelta,MATSEQAIJ);CHKERRQ(ierr); 424674ae819SStefano Zampini ierr = MatSeqAIJSetPreallocation(fetidpmat_ctx->B_Ddelta,1,NULL);CHKERRQ(ierr); 425674ae819SStefano Zampini for (i=0;i<n_local_lambda;i++) { 426674ae819SStefano Zampini ierr = MatSetValue(fetidpmat_ctx->B_Ddelta,i,cols_B_delta[i],scaling_factors[i],INSERT_VALUES);CHKERRQ(ierr); 427674ae819SStefano Zampini } 428674ae819SStefano Zampini ierr = MatAssemblyBegin(fetidpmat_ctx->B_Ddelta,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 429674ae819SStefano Zampini ierr = MatAssemblyEnd (fetidpmat_ctx->B_Ddelta,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 430674ae819SStefano Zampini } 431674ae819SStefano Zampini ierr = PetscFree(scaling_factors);CHKERRQ(ierr); 432674ae819SStefano Zampini ierr = PetscFree(cols_B_delta);CHKERRQ(ierr); 433674ae819SStefano Zampini 434457d4a33Sstefano_zampini /* Layout of multipliers */ 435457d4a33Sstefano_zampini ierr = PetscLayoutCreate(comm,&llay);CHKERRQ(ierr); 436457d4a33Sstefano_zampini ierr = PetscLayoutSetBlockSize(llay,1);CHKERRQ(ierr); 437457d4a33Sstefano_zampini ierr = PetscLayoutSetSize(llay,fetidpmat_ctx->n_lambda);CHKERRQ(ierr); 438457d4a33Sstefano_zampini ierr = PetscLayoutSetUp(llay);CHKERRQ(ierr); 439457d4a33Sstefano_zampini ierr = PetscLayoutGetLocalSize(llay,&fetidpmat_ctx->n);CHKERRQ(ierr); 440457d4a33Sstefano_zampini 441457d4a33Sstefano_zampini /* fetidpmat sizes */ 442457d4a33Sstefano_zampini fetidpmat_ctx->n += nPgl; 443457d4a33Sstefano_zampini fetidpmat_ctx->N = fetidpmat_ctx->n_lambda+nPg; 444457d4a33Sstefano_zampini 445457d4a33Sstefano_zampini /* Local work vector of multipliers */ 446457d4a33Sstefano_zampini ierr = VecCreate(PETSC_COMM_SELF,&fetidpmat_ctx->lambda_local);CHKERRQ(ierr); 447457d4a33Sstefano_zampini ierr = VecSetSizes(fetidpmat_ctx->lambda_local,n_local_lambda,n_local_lambda);CHKERRQ(ierr); 448457d4a33Sstefano_zampini ierr = VecSetType(fetidpmat_ctx->lambda_local,VECSEQ);CHKERRQ(ierr); 449457d4a33Sstefano_zampini 450457d4a33Sstefano_zampini /* Global vector for FETI-DP linear system */ 451457d4a33Sstefano_zampini ierr = VecCreate(comm,&fetidp_global);CHKERRQ(ierr); 452457d4a33Sstefano_zampini ierr = VecSetSizes(fetidp_global,fetidpmat_ctx->n,fetidpmat_ctx->N);CHKERRQ(ierr); 453457d4a33Sstefano_zampini ierr = VecSetType(fetidp_global,VECMPI);CHKERRQ(ierr); 454457d4a33Sstefano_zampini ierr = VecSetUp(fetidp_global);CHKERRQ(ierr); 455457d4a33Sstefano_zampini 4569eec4de8Sstefano_zampini /* Decide layout for fetidp dofs: if it is a saddle point problem 4579eec4de8Sstefano_zampini pressure is ordered first in the local part of the global vector 458457d4a33Sstefano_zampini of the FETI-DP linear system */ 459457d4a33Sstefano_zampini if (nPg) { 460af140850Sstefano_zampini IS IS_l2g_p,ais; 461457d4a33Sstefano_zampini PetscLayout alay; 462457d4a33Sstefano_zampini const PetscInt *idxs,*pranges,*aranges,*lranges; 463af140850Sstefano_zampini PetscInt *l2g_indices_p,rst; 464457d4a33Sstefano_zampini 465457d4a33Sstefano_zampini ierr = PetscMalloc1(nPl,&l2g_indices_p);CHKERRQ(ierr); 466457d4a33Sstefano_zampini ierr = VecGetLayout(fetidp_global,&alay);CHKERRQ(ierr); 467457d4a33Sstefano_zampini ierr = PetscLayoutGetRanges(alay,&aranges);CHKERRQ(ierr); 468457d4a33Sstefano_zampini ierr = PetscLayoutGetRanges(play,&pranges);CHKERRQ(ierr); 469457d4a33Sstefano_zampini ierr = PetscLayoutGetRanges(llay,&lranges);CHKERRQ(ierr); 470457d4a33Sstefano_zampini ierr = ISLocalToGlobalMappingGetIndices(l2gmap_p,&idxs);CHKERRQ(ierr); 471af140850Sstefano_zampini /* shift local to global indices for pressure */ 472457d4a33Sstefano_zampini for (i=0;i<nPl;i++) { 473457d4a33Sstefano_zampini PetscInt owner; 474457d4a33Sstefano_zampini 475457d4a33Sstefano_zampini ierr = PetscLayoutFindOwner(play,idxs[i],&owner);CHKERRQ(ierr); 476457d4a33Sstefano_zampini l2g_indices_p[i] = idxs[i]-pranges[owner]+aranges[owner]; 477457d4a33Sstefano_zampini } 478457d4a33Sstefano_zampini ierr = ISLocalToGlobalMappingRestoreIndices(l2gmap_p,&idxs);CHKERRQ(ierr); 479457d4a33Sstefano_zampini ierr = ISCreateGeneral(comm,nPl,l2g_indices_p,PETSC_OWN_POINTER,&IS_l2g_p);CHKERRQ(ierr); 480af140850Sstefano_zampini 481457d4a33Sstefano_zampini /* local to global scatter for interface pressure */ 482457d4a33Sstefano_zampini ierr = VecScatterCreate(fetidpmat_ctx->vP,NULL,fetidp_global,IS_l2g_p,&fetidpmat_ctx->l2g_p);CHKERRQ(ierr); 483457d4a33Sstefano_zampini ierr = ISDestroy(&IS_l2g_p);CHKERRQ(ierr); 484457d4a33Sstefano_zampini 485af140850Sstefano_zampini /* shift local to global indices for multipliers */ 486457d4a33Sstefano_zampini for (i=0;i<n_local_lambda;i++) { 487457d4a33Sstefano_zampini PetscInt owner,ps; 488457d4a33Sstefano_zampini 489457d4a33Sstefano_zampini ierr = PetscLayoutFindOwner(llay,l2g_indices[i],&owner);CHKERRQ(ierr); 490457d4a33Sstefano_zampini ps = pranges[owner+1]-pranges[owner]; 491457d4a33Sstefano_zampini l2g_indices[i] = l2g_indices[i]-lranges[owner]+aranges[owner]+ps; 492457d4a33Sstefano_zampini } 493457d4a33Sstefano_zampini 4949cc7774eSstefano_zampini /* scatter from alldofs to interface pressures global fetidp vector */ 4959cc7774eSstefano_zampini ierr = PetscLayoutGetRange(alay,&rst,NULL);CHKERRQ(ierr); 4969cc7774eSstefano_zampini ierr = ISCreateStride(comm,nPgl,rst,1,&ais);CHKERRQ(ierr); 497af140850Sstefano_zampini ierr = VecScatterCreate(pcis->vec1_global,pP,fetidp_global,ais,&fetidpmat_ctx->g2g_p);CHKERRQ(ierr); 4989cc7774eSstefano_zampini ierr = ISDestroy(&ais);CHKERRQ(ierr); 499457d4a33Sstefano_zampini } 500457d4a33Sstefano_zampini ierr = PetscLayoutDestroy(&llay);CHKERRQ(ierr); 501457d4a33Sstefano_zampini ierr = PetscLayoutDestroy(&play);CHKERRQ(ierr); 502457d4a33Sstefano_zampini ierr = ISCreateGeneral(comm,n_local_lambda,l2g_indices,PETSC_OWN_POINTER,&IS_l2g_lambda);CHKERRQ(ierr); 5039cc7774eSstefano_zampini /* scatter from local to global multipliers */ 504457d4a33Sstefano_zampini ierr = VecScatterCreate(fetidpmat_ctx->lambda_local,NULL,fetidp_global,IS_l2g_lambda,&fetidpmat_ctx->l2g_lambda);CHKERRQ(ierr); 505457d4a33Sstefano_zampini ierr = ISDestroy(&IS_l2g_lambda);CHKERRQ(ierr); 506457d4a33Sstefano_zampini ierr = ISLocalToGlobalMappingDestroy(&l2gmap_p);CHKERRQ(ierr); 507457d4a33Sstefano_zampini 508674ae819SStefano Zampini /* Create some vectors needed by fetidp */ 509674ae819SStefano Zampini ierr = VecDuplicate(pcis->vec1_B,&fetidpmat_ctx->temp_solution_B);CHKERRQ(ierr); 510674ae819SStefano Zampini ierr = VecDuplicate(pcis->vec1_D,&fetidpmat_ctx->temp_solution_D);CHKERRQ(ierr); 511674ae819SStefano Zampini 512674ae819SStefano Zampini test_fetidp = PETSC_FALSE; 513c5929fdfSBarry Smith ierr = PetscOptionsGetBool(NULL,NULL,"-fetidp_check",&test_fetidp,NULL);CHKERRQ(ierr); 514674ae819SStefano Zampini 515674ae819SStefano Zampini if (test_fetidp && !pcbddc->use_deluxe_scaling) { 5169eec4de8Sstefano_zampini Vec test_vec,test_vec_p = NULL; 5172d456bbaSstefano_zampini IS dirdofs; 5185b08dc53SStefano Zampini PetscReal real_value; 5195b08dc53SStefano Zampini 520674ae819SStefano Zampini ierr = PetscViewerASCIIGetStdout(comm,&viewer);CHKERRQ(ierr); 5211575c14dSBarry Smith ierr = PetscViewerASCIIPushSynchronized(viewer);CHKERRQ(ierr); 522457d4a33Sstefano_zampini ierr = PetscViewerASCIIPrintf(viewer,"----------FETI_DP SIZES--------------\n");CHKERRQ(ierr); 523457d4a33Sstefano_zampini ierr = PetscViewerASCIISynchronizedPrintf(viewer,"[%04d]: Sizes %D %D %D %D %D\n",rank,fetidpmat_ctx->n,fetidpmat_ctx->N,nPgl,nPg,nPl);CHKERRQ(ierr); 524457d4a33Sstefano_zampini ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 525674ae819SStefano Zampini ierr = PetscViewerASCIIPrintf(viewer,"----------FETI_DP TESTS--------------\n");CHKERRQ(ierr); 526674ae819SStefano Zampini ierr = PetscViewerASCIIPrintf(viewer,"All tests should return zero!\n");CHKERRQ(ierr); 527674ae819SStefano Zampini ierr = PetscViewerASCIIPrintf(viewer,"FETIDP MAT context in the ");CHKERRQ(ierr); 528674ae819SStefano Zampini if (fully_redundant) { 529674ae819SStefano Zampini ierr = PetscViewerASCIIPrintf(viewer,"fully redundant case for lagrange multipliers.\n");CHKERRQ(ierr); 530674ae819SStefano Zampini } else { 531674ae819SStefano Zampini ierr = PetscViewerASCIIPrintf(viewer,"Non-fully redundant case for lagrange multiplier.\n");CHKERRQ(ierr); 532674ae819SStefano Zampini } 533674ae819SStefano Zampini ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 534674ae819SStefano Zampini 535674ae819SStefano Zampini /******************************************************************/ 5369eec4de8Sstefano_zampini /* TEST A/B: Test numbering of global fetidp dofs */ 537674ae819SStefano Zampini /******************************************************************/ 538674ae819SStefano Zampini 539674ae819SStefano Zampini ierr = VecDuplicate(fetidpmat_ctx->lambda_local,&test_vec);CHKERRQ(ierr); 540457d4a33Sstefano_zampini ierr = VecSet(fetidp_global,1.0);CHKERRQ(ierr); 5419eec4de8Sstefano_zampini ierr = VecSet(test_vec,1.);CHKERRQ(ierr); 542457d4a33Sstefano_zampini ierr = VecScatterBegin(fetidpmat_ctx->l2g_lambda,fetidp_global,fetidpmat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 543457d4a33Sstefano_zampini ierr = VecScatterEnd(fetidpmat_ctx->l2g_lambda,fetidp_global,fetidpmat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 5449eec4de8Sstefano_zampini if (fetidpmat_ctx->l2g_p) { 5459eec4de8Sstefano_zampini ierr = VecDuplicate(fetidpmat_ctx->vP,&test_vec_p);CHKERRQ(ierr); 5469eec4de8Sstefano_zampini ierr = VecSet(test_vec_p,1.);CHKERRQ(ierr); 5479eec4de8Sstefano_zampini ierr = VecScatterBegin(fetidpmat_ctx->l2g_p,fetidp_global,fetidpmat_ctx->vP,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 5489eec4de8Sstefano_zampini ierr = VecScatterEnd(fetidpmat_ctx->l2g_p,fetidp_global,fetidpmat_ctx->vP,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 5499eec4de8Sstefano_zampini } 550674ae819SStefano Zampini scalar_value = -1.0; 551674ae819SStefano Zampini ierr = VecAXPY(test_vec,scalar_value,fetidpmat_ctx->lambda_local);CHKERRQ(ierr); 5525b08dc53SStefano Zampini ierr = VecNorm(test_vec,NORM_INFINITY,&real_value);CHKERRQ(ierr); 553674ae819SStefano Zampini ierr = VecDestroy(&test_vec);CHKERRQ(ierr); 5545b08dc53SStefano Zampini ierr = PetscViewerASCIISynchronizedPrintf(viewer,"A[%04d]: CHECK glob to loc: % 1.14e\n",rank,real_value);CHKERRQ(ierr); 5559eec4de8Sstefano_zampini if (fetidpmat_ctx->l2g_p) { 5569eec4de8Sstefano_zampini ierr = VecAXPY(test_vec_p,scalar_value,fetidpmat_ctx->vP);CHKERRQ(ierr); 5579eec4de8Sstefano_zampini ierr = VecNorm(test_vec_p,NORM_INFINITY,&real_value);CHKERRQ(ierr); 5589eec4de8Sstefano_zampini ierr = PetscViewerASCIISynchronizedPrintf(viewer,"A[%04d]: CHECK glob to loc (p): % 1.14e\n",rank,real_value);CHKERRQ(ierr); 5599eec4de8Sstefano_zampini } 560674ae819SStefano Zampini ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 561674ae819SStefano Zampini if (fully_redundant) { 562457d4a33Sstefano_zampini ierr = VecSet(fetidp_global,0.0);CHKERRQ(ierr); 563674ae819SStefano Zampini ierr = VecSet(fetidpmat_ctx->lambda_local,0.5);CHKERRQ(ierr); 564457d4a33Sstefano_zampini ierr = VecScatterBegin(fetidpmat_ctx->l2g_lambda,fetidpmat_ctx->lambda_local,fetidp_global,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 565457d4a33Sstefano_zampini ierr = VecScatterEnd(fetidpmat_ctx->l2g_lambda,fetidpmat_ctx->lambda_local,fetidp_global,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 566457d4a33Sstefano_zampini ierr = VecSum(fetidp_global,&scalar_value);CHKERRQ(ierr); 5675b08dc53SStefano Zampini ierr = PetscViewerASCIISynchronizedPrintf(viewer,"B[%04d]: CHECK loc to glob: % 1.14e\n",rank,PetscRealPart(scalar_value)-fetidpmat_ctx->n_lambda);CHKERRQ(ierr); 568674ae819SStefano Zampini ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 569674ae819SStefano Zampini } 5709eec4de8Sstefano_zampini if (fetidpmat_ctx->l2g_p) { 571af140850Sstefano_zampini ierr = VecSet(pcis->vec1_N,1.0);CHKERRQ(ierr); 572af140850Sstefano_zampini ierr = VecSet(pcis->vec1_global,0.0);CHKERRQ(ierr); 573af140850Sstefano_zampini ierr = VecScatterBegin(matis->rctx,pcis->vec1_N,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 574af140850Sstefano_zampini ierr = VecScatterEnd(matis->rctx,pcis->vec1_N,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 575af140850Sstefano_zampini 5769eec4de8Sstefano_zampini ierr = VecSet(fetidp_global,0.0);CHKERRQ(ierr); 577af140850Sstefano_zampini ierr = VecSet(fetidpmat_ctx->vP,-1.0);CHKERRQ(ierr); 578af140850Sstefano_zampini ierr = VecScatterBegin(fetidpmat_ctx->l2g_p,fetidpmat_ctx->vP,fetidp_global,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 579af140850Sstefano_zampini ierr = VecScatterEnd(fetidpmat_ctx->l2g_p,fetidpmat_ctx->vP,fetidp_global,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 580af140850Sstefano_zampini ierr = VecScatterBegin(fetidpmat_ctx->g2g_p,fetidp_global,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 581af140850Sstefano_zampini ierr = VecScatterEnd(fetidpmat_ctx->g2g_p,fetidp_global,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 582af140850Sstefano_zampini ierr = VecScatterBegin(fetidpmat_ctx->g2g_p,pcis->vec1_global,fetidp_global,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 583af140850Sstefano_zampini ierr = VecScatterEnd(fetidpmat_ctx->g2g_p,pcis->vec1_global,fetidp_global,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 5849eec4de8Sstefano_zampini ierr = VecSum(fetidp_global,&scalar_value);CHKERRQ(ierr); 585af140850Sstefano_zampini ierr = PetscViewerASCIISynchronizedPrintf(viewer,"B[%04d]: CHECK loc to glob (p): % 1.14e\n",rank,PetscRealPart(scalar_value));CHKERRQ(ierr); 5862d456bbaSstefano_zampini ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 5879eec4de8Sstefano_zampini } 588674ae819SStefano Zampini /******************************************************************/ 589674ae819SStefano Zampini /* TEST C: It should holds B_delta*w=0, w\in\widehat{W} */ 590674ae819SStefano Zampini /* This is the meaning of the B matrix */ 591674ae819SStefano Zampini /******************************************************************/ 592674ae819SStefano Zampini 593674ae819SStefano Zampini ierr = VecSetRandom(pcis->vec1_N,NULL);CHKERRQ(ierr); 594674ae819SStefano Zampini ierr = VecSet(pcis->vec1_global,0.0);CHKERRQ(ierr); 595e176bc59SStefano Zampini ierr = VecScatterBegin(matis->rctx,pcis->vec1_N,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 596e176bc59SStefano Zampini ierr = VecScatterEnd(matis->rctx,pcis->vec1_N,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 597e176bc59SStefano Zampini ierr = VecScatterBegin(matis->rctx,pcis->vec1_global,pcis->vec1_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 598e176bc59SStefano Zampini ierr = VecScatterEnd(matis->rctx,pcis->vec1_global,pcis->vec1_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 599674ae819SStefano Zampini ierr = VecScatterBegin(pcis->N_to_B,pcis->vec1_N,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 600674ae819SStefano Zampini ierr = VecScatterEnd(pcis->N_to_B,pcis->vec1_N,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 601674ae819SStefano Zampini /* Action of B_delta */ 602674ae819SStefano Zampini ierr = MatMult(fetidpmat_ctx->B_delta,pcis->vec1_B,fetidpmat_ctx->lambda_local);CHKERRQ(ierr); 603457d4a33Sstefano_zampini ierr = VecSet(fetidp_global,0.0);CHKERRQ(ierr); 604457d4a33Sstefano_zampini ierr = VecScatterBegin(fetidpmat_ctx->l2g_lambda,fetidpmat_ctx->lambda_local,fetidp_global,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 605457d4a33Sstefano_zampini ierr = VecScatterEnd(fetidpmat_ctx->l2g_lambda,fetidpmat_ctx->lambda_local,fetidp_global,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 606457d4a33Sstefano_zampini ierr = VecNorm(fetidp_global,NORM_INFINITY,&real_value);CHKERRQ(ierr); 6075b08dc53SStefano Zampini ierr = PetscViewerASCIIPrintf(viewer,"C[coll]: CHECK infty norm of B_delta*w (w continuous): % 1.14e\n",real_value);CHKERRQ(ierr); 608674ae819SStefano Zampini ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 609674ae819SStefano Zampini 610674ae819SStefano Zampini /******************************************************************/ 611674ae819SStefano Zampini /* TEST D: It should holds E_Dw = w - P_Dw w\in\widetilde{W} */ 612674ae819SStefano Zampini /* E_D = R_D^TR */ 613674ae819SStefano Zampini /* P_D = B_{D,delta}^T B_{delta} */ 614674ae819SStefano Zampini /* eq.44 Mandel Tezaur and Dohrmann 2005 */ 615674ae819SStefano Zampini /******************************************************************/ 616674ae819SStefano Zampini 617674ae819SStefano Zampini /* compute a random vector in \widetilde{W} */ 618674ae819SStefano Zampini ierr = VecSetRandom(pcis->vec1_N,NULL);CHKERRQ(ierr); 6192d456bbaSstefano_zampini /* set zero at vertices and essential dofs */ 6202d456bbaSstefano_zampini scalar_value = 0.0; 621674ae819SStefano Zampini ierr = VecGetArray(pcis->vec1_N,&array);CHKERRQ(ierr); 622674ae819SStefano Zampini for (i=0;i<n_vertices;i++) { array[vertex_indices[i]]=scalar_value; } 6232d456bbaSstefano_zampini ierr = PCBDDCGraphGetDirichletDofs(pcbddc->mat_graph,&dirdofs);CHKERRQ(ierr); 6242d456bbaSstefano_zampini if (dirdofs) { 6252d456bbaSstefano_zampini const PetscInt *idxs; 6262d456bbaSstefano_zampini PetscInt ndir; 6272d456bbaSstefano_zampini 6282d456bbaSstefano_zampini ierr = ISGetLocalSize(dirdofs,&ndir);CHKERRQ(ierr); 6292d456bbaSstefano_zampini ierr = ISGetIndices(dirdofs,&idxs);CHKERRQ(ierr); 6302d456bbaSstefano_zampini for (i=0;i<ndir;i++) { array[idxs[i]]=scalar_value; } 6312d456bbaSstefano_zampini ierr = ISRestoreIndices(dirdofs,&idxs);CHKERRQ(ierr); 6322d456bbaSstefano_zampini } 633674ae819SStefano Zampini ierr = VecRestoreArray(pcis->vec1_N,&array);CHKERRQ(ierr); 634674ae819SStefano Zampini /* store w for final comparison */ 635674ae819SStefano Zampini ierr = VecDuplicate(pcis->vec1_B,&test_vec);CHKERRQ(ierr); 636674ae819SStefano Zampini ierr = VecScatterBegin(pcis->N_to_B,pcis->vec1_N,test_vec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 637674ae819SStefano Zampini ierr = VecScatterEnd(pcis->N_to_B,pcis->vec1_N,test_vec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 638674ae819SStefano Zampini 639674ae819SStefano Zampini /* Jump operator P_D : results stored in pcis->vec1_B */ 640674ae819SStefano Zampini 641674ae819SStefano Zampini ierr = VecScatterBegin(pcis->N_to_B,pcis->vec1_N,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 642674ae819SStefano Zampini ierr = VecScatterEnd(pcis->N_to_B,pcis->vec1_N,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 643674ae819SStefano Zampini /* Action of B_delta */ 644674ae819SStefano Zampini ierr = MatMult(fetidpmat_ctx->B_delta,pcis->vec1_B,fetidpmat_ctx->lambda_local);CHKERRQ(ierr); 645457d4a33Sstefano_zampini ierr = VecSet(fetidp_global,0.0);CHKERRQ(ierr); 646457d4a33Sstefano_zampini ierr = VecScatterBegin(fetidpmat_ctx->l2g_lambda,fetidpmat_ctx->lambda_local,fetidp_global,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 647457d4a33Sstefano_zampini ierr = VecScatterEnd(fetidpmat_ctx->l2g_lambda,fetidpmat_ctx->lambda_local,fetidp_global,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 648674ae819SStefano Zampini /* Action of B_Ddelta^T */ 649457d4a33Sstefano_zampini ierr = VecScatterBegin(fetidpmat_ctx->l2g_lambda,fetidp_global,fetidpmat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 650457d4a33Sstefano_zampini ierr = VecScatterEnd(fetidpmat_ctx->l2g_lambda,fetidp_global,fetidpmat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 651674ae819SStefano Zampini ierr = MatMultTranspose(fetidpmat_ctx->B_Ddelta,fetidpmat_ctx->lambda_local,pcis->vec1_B);CHKERRQ(ierr); 652674ae819SStefano Zampini 653674ae819SStefano Zampini /* Average operator E_D : results stored in pcis->vec2_B */ 654674ae819SStefano Zampini ierr = VecScatterBegin(pcis->N_to_B,pcis->vec1_N,pcis->vec2_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 655674ae819SStefano Zampini ierr = VecScatterEnd(pcis->N_to_B,pcis->vec1_N,pcis->vec2_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 656674ae819SStefano Zampini ierr = PCBDDCScalingExtension(fetidpmat_ctx->pc,pcis->vec2_B,pcis->vec1_global);CHKERRQ(ierr); 657674ae819SStefano Zampini ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_global,pcis->vec2_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 658674ae819SStefano Zampini ierr = VecScatterEnd(pcis->global_to_B,pcis->vec1_global,pcis->vec2_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 659674ae819SStefano Zampini 660674ae819SStefano Zampini /* test E_D=I-P_D */ 661674ae819SStefano Zampini scalar_value = 1.0; 662674ae819SStefano Zampini ierr = VecAXPY(pcis->vec1_B,scalar_value,pcis->vec2_B);CHKERRQ(ierr); 663674ae819SStefano Zampini scalar_value = -1.0; 664674ae819SStefano Zampini ierr = VecAXPY(pcis->vec1_B,scalar_value,test_vec);CHKERRQ(ierr); 6655b08dc53SStefano Zampini ierr = VecNorm(pcis->vec1_B,NORM_INFINITY,&real_value);CHKERRQ(ierr); 666674ae819SStefano Zampini ierr = VecDestroy(&test_vec);CHKERRQ(ierr); 6675b08dc53SStefano Zampini ierr = PetscViewerASCIISynchronizedPrintf(viewer,"D[%04d] CHECK infty norm of E_D + P_D - I: % 1.14e\n",rank,real_value);CHKERRQ(ierr); 668674ae819SStefano Zampini ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 669674ae819SStefano Zampini 670674ae819SStefano Zampini /******************************************************************/ 6712d456bbaSstefano_zampini /* TEST E: It should hold R_D^TP_Dw=0 w\in\widetilde{W} */ 672674ae819SStefano Zampini /* eq.48 Mandel Tezaur and Dohrmann 2005 */ 673674ae819SStefano Zampini /******************************************************************/ 674674ae819SStefano Zampini 675674ae819SStefano Zampini ierr = VecSetRandom(pcis->vec1_N,NULL);CHKERRQ(ierr); 6762d456bbaSstefano_zampini /* set zero at vertices and essential dofs */ 6772d456bbaSstefano_zampini scalar_value = 0.0; 678674ae819SStefano Zampini ierr = VecGetArray(pcis->vec1_N,&array);CHKERRQ(ierr); 679674ae819SStefano Zampini for (i=0;i<n_vertices;i++) { array[vertex_indices[i]]=scalar_value; } 6802d456bbaSstefano_zampini if (dirdofs) { 6812d456bbaSstefano_zampini const PetscInt *idxs; 6822d456bbaSstefano_zampini PetscInt ndir; 6832d456bbaSstefano_zampini 6842d456bbaSstefano_zampini ierr = ISGetLocalSize(dirdofs,&ndir);CHKERRQ(ierr); 6852d456bbaSstefano_zampini ierr = ISGetIndices(dirdofs,&idxs);CHKERRQ(ierr); 6862d456bbaSstefano_zampini for (i=0;i<ndir;i++) { array[idxs[i]]=scalar_value; } 6872d456bbaSstefano_zampini ierr = ISRestoreIndices(dirdofs,&idxs);CHKERRQ(ierr); 6882d456bbaSstefano_zampini } 689674ae819SStefano Zampini ierr = VecRestoreArray(pcis->vec1_N,&array);CHKERRQ(ierr); 690674ae819SStefano Zampini 691674ae819SStefano Zampini /* Jump operator P_D : results stored in pcis->vec1_B */ 692674ae819SStefano Zampini 693674ae819SStefano Zampini ierr = VecScatterBegin(pcis->N_to_B,pcis->vec1_N,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 694674ae819SStefano Zampini ierr = VecScatterEnd(pcis->N_to_B,pcis->vec1_N,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 695674ae819SStefano Zampini /* Action of B_delta */ 696674ae819SStefano Zampini ierr = MatMult(fetidpmat_ctx->B_delta,pcis->vec1_B,fetidpmat_ctx->lambda_local);CHKERRQ(ierr); 697457d4a33Sstefano_zampini ierr = VecSet(fetidp_global,0.0);CHKERRQ(ierr); 698457d4a33Sstefano_zampini ierr = VecScatterBegin(fetidpmat_ctx->l2g_lambda,fetidpmat_ctx->lambda_local,fetidp_global,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 699457d4a33Sstefano_zampini ierr = VecScatterEnd(fetidpmat_ctx->l2g_lambda,fetidpmat_ctx->lambda_local,fetidp_global,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 700674ae819SStefano Zampini /* Action of B_Ddelta^T */ 701457d4a33Sstefano_zampini ierr = VecScatterBegin(fetidpmat_ctx->l2g_lambda,fetidp_global,fetidpmat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 702457d4a33Sstefano_zampini ierr = VecScatterEnd(fetidpmat_ctx->l2g_lambda,fetidp_global,fetidpmat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 703674ae819SStefano Zampini ierr = MatMultTranspose(fetidpmat_ctx->B_Ddelta,fetidpmat_ctx->lambda_local,pcis->vec1_B);CHKERRQ(ierr); 704674ae819SStefano Zampini /* scaling */ 705674ae819SStefano Zampini ierr = PCBDDCScalingExtension(fetidpmat_ctx->pc,pcis->vec1_B,pcis->vec1_global);CHKERRQ(ierr); 7065b08dc53SStefano Zampini ierr = VecNorm(pcis->vec1_global,NORM_INFINITY,&real_value);CHKERRQ(ierr); 7075b08dc53SStefano Zampini ierr = PetscViewerASCIIPrintf(viewer,"E[coll]: CHECK infty norm of R^T_D P_D: % 1.14e\n",real_value);CHKERRQ(ierr); 708674ae819SStefano Zampini ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 709674ae819SStefano Zampini 710674ae819SStefano Zampini if (!fully_redundant) { 711674ae819SStefano Zampini /******************************************************************/ 712674ae819SStefano Zampini /* TEST F: It should holds B_{delta}B^T_{D,delta}=I */ 713674ae819SStefano Zampini /* Corollary thm 14 Mandel Tezaur and Dohrmann 2005 */ 714674ae819SStefano Zampini /******************************************************************/ 715457d4a33Sstefano_zampini ierr = VecDuplicate(fetidp_global,&test_vec);CHKERRQ(ierr); 716457d4a33Sstefano_zampini ierr = VecSetRandom(fetidp_global,NULL);CHKERRQ(ierr); 717457d4a33Sstefano_zampini if (fetidpmat_ctx->l2g_p) { 718457d4a33Sstefano_zampini ierr = VecSet(fetidpmat_ctx->vP,0.);CHKERRQ(ierr); 719457d4a33Sstefano_zampini ierr = VecScatterBegin(fetidpmat_ctx->l2g_p,fetidpmat_ctx->vP,fetidp_global,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 720457d4a33Sstefano_zampini ierr = VecScatterEnd(fetidpmat_ctx->l2g_p,fetidpmat_ctx->vP,fetidp_global,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 721457d4a33Sstefano_zampini } 722674ae819SStefano Zampini /* Action of B_Ddelta^T */ 723457d4a33Sstefano_zampini ierr = VecScatterBegin(fetidpmat_ctx->l2g_lambda,fetidp_global,fetidpmat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 724457d4a33Sstefano_zampini ierr = VecScatterEnd(fetidpmat_ctx->l2g_lambda,fetidp_global,fetidpmat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 725674ae819SStefano Zampini ierr = MatMultTranspose(fetidpmat_ctx->B_Ddelta,fetidpmat_ctx->lambda_local,pcis->vec1_B);CHKERRQ(ierr); 726674ae819SStefano Zampini /* Action of B_delta */ 727674ae819SStefano Zampini ierr = MatMult(fetidpmat_ctx->B_delta,pcis->vec1_B,fetidpmat_ctx->lambda_local);CHKERRQ(ierr); 728674ae819SStefano Zampini ierr = VecSet(test_vec,0.0);CHKERRQ(ierr); 729674ae819SStefano Zampini ierr = VecScatterBegin(fetidpmat_ctx->l2g_lambda,fetidpmat_ctx->lambda_local,test_vec,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 730674ae819SStefano Zampini ierr = VecScatterEnd(fetidpmat_ctx->l2g_lambda,fetidpmat_ctx->lambda_local,test_vec,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 731674ae819SStefano Zampini scalar_value = -1.0; 732457d4a33Sstefano_zampini ierr = VecAXPY(fetidp_global,scalar_value,test_vec);CHKERRQ(ierr); 733457d4a33Sstefano_zampini ierr = VecNorm(fetidp_global,NORM_INFINITY,&real_value);CHKERRQ(ierr); 7345b08dc53SStefano Zampini ierr = PetscViewerASCIIPrintf(viewer,"E[coll]: CHECK infty norm of P^T_D - I: % 1.14e\n",real_value);CHKERRQ(ierr); 735674ae819SStefano Zampini ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 736674ae819SStefano Zampini ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 737674ae819SStefano Zampini ierr = VecDestroy(&test_vec);CHKERRQ(ierr); 738674ae819SStefano Zampini } 7399eec4de8Sstefano_zampini ierr = VecDestroy(&test_vec_p);CHKERRQ(ierr); 7402d456bbaSstefano_zampini ierr = ISDestroy(&dirdofs);CHKERRQ(ierr); 741674ae819SStefano Zampini } 742674ae819SStefano Zampini /* final cleanup */ 743457d4a33Sstefano_zampini ierr = VecDestroy(&fetidp_global);CHKERRQ(ierr); 744674ae819SStefano Zampini PetscFunctionReturn(0); 745674ae819SStefano Zampini } 746674ae819SStefano Zampini 747674ae819SStefano Zampini #undef __FUNCT__ 748674ae819SStefano Zampini #define __FUNCT__ "PCBDDCSetupFETIDPPCContext" 749674ae819SStefano Zampini PetscErrorCode PCBDDCSetupFETIDPPCContext(Mat fetimat, FETIDPPC_ctx fetidppc_ctx) 750674ae819SStefano Zampini { 751674ae819SStefano Zampini FETIDPMat_ctx mat_ctx; 752ed6c3d69SStefano Zampini PC_IS *pcis; 753f28b6018SStefano Zampini PetscBool lumped = PETSC_FALSE; 754674ae819SStefano Zampini PetscErrorCode ierr; 755674ae819SStefano Zampini 756674ae819SStefano Zampini PetscFunctionBegin; 757674ae819SStefano Zampini ierr = MatShellGetContext(fetimat,(void**)&mat_ctx);CHKERRQ(ierr); 758674ae819SStefano Zampini /* get references from objects created when setting up feti mat context */ 759674ae819SStefano Zampini ierr = PetscObjectReference((PetscObject)mat_ctx->lambda_local);CHKERRQ(ierr); 760674ae819SStefano Zampini fetidppc_ctx->lambda_local = mat_ctx->lambda_local; 761674ae819SStefano Zampini ierr = PetscObjectReference((PetscObject)mat_ctx->B_Ddelta);CHKERRQ(ierr); 762674ae819SStefano Zampini fetidppc_ctx->B_Ddelta = mat_ctx->B_Ddelta; 763674ae819SStefano Zampini ierr = PetscObjectReference((PetscObject)mat_ctx->l2g_lambda);CHKERRQ(ierr); 764674ae819SStefano Zampini fetidppc_ctx->l2g_lambda = mat_ctx->l2g_lambda; 765f28b6018SStefano Zampini /* create preconditioner */ 766ed6c3d69SStefano Zampini pcis = (PC_IS*)fetidppc_ctx->pc->data; 767f28b6018SStefano Zampini /* Dirichlet preconditioner */ 768f28b6018SStefano Zampini ierr = PetscOptionsGetBool(NULL,NULL,"-fetidp_lumped",&lumped,NULL);CHKERRQ(ierr); 769f28b6018SStefano Zampini if (!lumped) { 770ed6c3d69SStefano Zampini ierr = MatCreateSchurComplement(pcis->A_II,pcis->A_II,pcis->A_IB,pcis->A_BI,pcis->A_BB,&fetidppc_ctx->S_j);CHKERRQ(ierr); 771ed6c3d69SStefano Zampini ierr = MatSchurComplementSetKSP(fetidppc_ctx->S_j,pcis->ksp_D);CHKERRQ(ierr); 772f28b6018SStefano Zampini } else { 773f28b6018SStefano Zampini ierr = PetscObjectReference((PetscObject)pcis->A_BB);CHKERRQ(ierr); 774f28b6018SStefano Zampini fetidppc_ctx->S_j = pcis->A_BB; 775f28b6018SStefano Zampini } 776af140850Sstefano_zampini /* saddle-point */ 777af140850Sstefano_zampini if (mat_ctx->xPg) { 778af140850Sstefano_zampini ierr = PetscObjectReference((PetscObject)mat_ctx->xPg);CHKERRQ(ierr); 779af140850Sstefano_zampini fetidppc_ctx->xPg = mat_ctx->xPg; 780af140850Sstefano_zampini ierr = PetscObjectReference((PetscObject)mat_ctx->yPg);CHKERRQ(ierr); 781af140850Sstefano_zampini fetidppc_ctx->yPg = mat_ctx->yPg; 782*6cc1294bSstefano_zampini ierr = PetscObjectQuery((PetscObject)fetidppc_ctx->pc,"__KSPFETIDP_PKSP",(PetscObject*)&fetidppc_ctx->kP);CHKERRQ(ierr); 783*6cc1294bSstefano_zampini ierr = PetscObjectReference((PetscObject)fetidppc_ctx->kP);CHKERRQ(ierr); 784af140850Sstefano_zampini } 785674ae819SStefano Zampini PetscFunctionReturn(0); 786674ae819SStefano Zampini } 787674ae819SStefano Zampini 788674ae819SStefano Zampini #undef __FUNCT__ 789674ae819SStefano Zampini #define __FUNCT__ "FETIDPMatMult" 790674ae819SStefano Zampini PetscErrorCode FETIDPMatMult(Mat fetimat, Vec x, Vec y) 791674ae819SStefano Zampini { 792674ae819SStefano Zampini FETIDPMat_ctx mat_ctx; 793617d11aeSStefano Zampini PC_BDDC *pcbddc; 794674ae819SStefano Zampini PC_IS *pcis; 795674ae819SStefano Zampini PetscErrorCode ierr; 796674ae819SStefano Zampini 797674ae819SStefano Zampini PetscFunctionBegin; 798674ae819SStefano Zampini ierr = MatShellGetContext(fetimat,(void**)&mat_ctx);CHKERRQ(ierr); 799674ae819SStefano Zampini pcis = (PC_IS*)mat_ctx->pc->data; 800617d11aeSStefano Zampini pcbddc = (PC_BDDC*)mat_ctx->pc->data; 801674ae819SStefano Zampini /* Application of B_delta^T */ 802af140850Sstefano_zampini ierr = VecSet(pcis->vec1_B,0.);CHKERRQ(ierr); 803674ae819SStefano Zampini ierr = VecScatterBegin(mat_ctx->l2g_lambda,x,mat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 804674ae819SStefano Zampini ierr = VecScatterEnd(mat_ctx->l2g_lambda,x,mat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 805674ae819SStefano Zampini ierr = MatMultTranspose(mat_ctx->B_delta,mat_ctx->lambda_local,pcis->vec1_B);CHKERRQ(ierr); 806af140850Sstefano_zampini 807af140850Sstefano_zampini /* Add contribution from saddle point */ 808af140850Sstefano_zampini if (mat_ctx->l2g_p) { 809af140850Sstefano_zampini ierr = VecScatterBegin(mat_ctx->l2g_p,x,mat_ctx->vP,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 810af140850Sstefano_zampini ierr = VecScatterEnd(mat_ctx->l2g_p,x,mat_ctx->vP,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 811af140850Sstefano_zampini if (pcbddc->switch_static) { 812af140850Sstefano_zampini ierr = MatMult(mat_ctx->Bt_BI,mat_ctx->vP,pcis->vec1_D);CHKERRQ(ierr); 813af140850Sstefano_zampini } 814af140850Sstefano_zampini ierr = MatMultAdd(mat_ctx->Bt_BB,mat_ctx->vP,pcis->vec1_B,pcis->vec1_B);CHKERRQ(ierr); 815af140850Sstefano_zampini } else { 816af140850Sstefano_zampini if (pcbddc->switch_static) { 817674ae819SStefano Zampini ierr = VecSet(pcis->vec1_D,0.0);CHKERRQ(ierr); 818af140850Sstefano_zampini } 819af140850Sstefano_zampini } 820af140850Sstefano_zampini /* Application of \widetilde{S}^-1 */ 821617d11aeSStefano Zampini ierr = PetscMemzero(pcbddc->benign_p0,pcbddc->benign_n*sizeof(PetscScalar));CHKERRQ(ierr); 822dc359a40SStefano Zampini ierr = PCBDDCApplyInterfacePreconditioner(mat_ctx->pc,PETSC_FALSE);CHKERRQ(ierr); 823c7ffc8ceSStefano Zampini ierr = PetscMemzero(pcbddc->benign_p0,pcbddc->benign_n*sizeof(PetscScalar));CHKERRQ(ierr); 824af140850Sstefano_zampini ierr = VecSet(y,0.0);CHKERRQ(ierr); 825674ae819SStefano Zampini /* Application of B_delta */ 826674ae819SStefano Zampini ierr = MatMult(mat_ctx->B_delta,pcis->vec1_B,mat_ctx->lambda_local);CHKERRQ(ierr); 827af140850Sstefano_zampini /* Contribution from boundary pressures */ 828af140850Sstefano_zampini if (mat_ctx->C) { 829af140850Sstefano_zampini const PetscScalar *lx; 830af140850Sstefano_zampini PetscScalar *ly; 831af140850Sstefano_zampini 832af140850Sstefano_zampini /* pressures ordered first in x and y */ 833af140850Sstefano_zampini ierr = VecGetArrayRead(x,&lx);CHKERRQ(ierr); 834af140850Sstefano_zampini ierr = VecGetArray(y,&ly);CHKERRQ(ierr); 835af140850Sstefano_zampini ierr = VecPlaceArray(mat_ctx->xPg,lx);CHKERRQ(ierr); 836af140850Sstefano_zampini ierr = VecPlaceArray(mat_ctx->yPg,ly);CHKERRQ(ierr); 837af140850Sstefano_zampini ierr = MatMult(mat_ctx->C,mat_ctx->xPg,mat_ctx->yPg);CHKERRQ(ierr); 838af140850Sstefano_zampini ierr = VecResetArray(mat_ctx->xPg);CHKERRQ(ierr); 839af140850Sstefano_zampini ierr = VecResetArray(mat_ctx->yPg);CHKERRQ(ierr); 840af140850Sstefano_zampini ierr = VecRestoreArrayRead(x,&lx);CHKERRQ(ierr); 841af140850Sstefano_zampini ierr = VecRestoreArray(y,&ly);CHKERRQ(ierr); 842af140850Sstefano_zampini } 843af140850Sstefano_zampini /* Add contribution from saddle point */ 844af140850Sstefano_zampini if (mat_ctx->l2g_p) { 845af140850Sstefano_zampini ierr = MatMult(mat_ctx->B_BB,pcis->vec1_B,mat_ctx->vP);CHKERRQ(ierr); 846af140850Sstefano_zampini if (pcbddc->switch_static) { 847af140850Sstefano_zampini ierr = MatMultAdd(mat_ctx->B_BI,pcis->vec1_D,mat_ctx->vP,mat_ctx->vP);CHKERRQ(ierr); 848af140850Sstefano_zampini } 849af140850Sstefano_zampini ierr = VecScatterBegin(mat_ctx->l2g_p,mat_ctx->vP,y,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 850af140850Sstefano_zampini ierr = VecScatterEnd(mat_ctx->l2g_p,mat_ctx->vP,y,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 851af140850Sstefano_zampini } 852674ae819SStefano Zampini ierr = VecScatterBegin(mat_ctx->l2g_lambda,mat_ctx->lambda_local,y,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 853674ae819SStefano Zampini ierr = VecScatterEnd(mat_ctx->l2g_lambda,mat_ctx->lambda_local,y,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 854674ae819SStefano Zampini PetscFunctionReturn(0); 855674ae819SStefano Zampini } 856674ae819SStefano Zampini 857674ae819SStefano Zampini #undef __FUNCT__ 858edf7251bSStefano Zampini #define __FUNCT__ "FETIDPMatMultTranspose" 859edf7251bSStefano Zampini PetscErrorCode FETIDPMatMultTranspose(Mat fetimat, Vec x, Vec y) 860edf7251bSStefano Zampini { 861edf7251bSStefano Zampini FETIDPMat_ctx mat_ctx; 862edf7251bSStefano Zampini PC_IS *pcis; 863edf7251bSStefano Zampini PetscErrorCode ierr; 864edf7251bSStefano Zampini 865edf7251bSStefano Zampini PetscFunctionBegin; 866edf7251bSStefano Zampini ierr = MatShellGetContext(fetimat,(void**)&mat_ctx);CHKERRQ(ierr); 867edf7251bSStefano Zampini pcis = (PC_IS*)mat_ctx->pc->data; 868edf7251bSStefano Zampini /* Application of B_delta^T */ 869edf7251bSStefano Zampini ierr = VecScatterBegin(mat_ctx->l2g_lambda,x,mat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 870edf7251bSStefano Zampini ierr = VecScatterEnd(mat_ctx->l2g_lambda,x,mat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 871edf7251bSStefano Zampini ierr = MatMultTranspose(mat_ctx->B_delta,mat_ctx->lambda_local,pcis->vec1_B);CHKERRQ(ierr); 872edf7251bSStefano Zampini /* Application of \widetilde{S}^-1 */ 873edf7251bSStefano Zampini ierr = VecSet(pcis->vec1_D,0.0);CHKERRQ(ierr); 874edf7251bSStefano Zampini ierr = PCBDDCApplyInterfacePreconditioner(mat_ctx->pc,PETSC_TRUE);CHKERRQ(ierr); 875edf7251bSStefano Zampini /* Application of B_delta */ 876edf7251bSStefano Zampini ierr = MatMult(mat_ctx->B_delta,pcis->vec1_B,mat_ctx->lambda_local);CHKERRQ(ierr); 877edf7251bSStefano Zampini ierr = VecSet(y,0.0);CHKERRQ(ierr); 878edf7251bSStefano Zampini ierr = VecScatterBegin(mat_ctx->l2g_lambda,mat_ctx->lambda_local,y,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 879edf7251bSStefano Zampini ierr = VecScatterEnd(mat_ctx->l2g_lambda,mat_ctx->lambda_local,y,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 880edf7251bSStefano Zampini PetscFunctionReturn(0); 881edf7251bSStefano Zampini } 882edf7251bSStefano Zampini 883edf7251bSStefano Zampini #undef __FUNCT__ 884674ae819SStefano Zampini #define __FUNCT__ "FETIDPPCApply" 885674ae819SStefano Zampini PetscErrorCode FETIDPPCApply(PC fetipc, Vec x, Vec y) 886674ae819SStefano Zampini { 887674ae819SStefano Zampini FETIDPPC_ctx pc_ctx; 888674ae819SStefano Zampini PC_IS *pcis; 889674ae819SStefano Zampini PetscErrorCode ierr; 890674ae819SStefano Zampini 891674ae819SStefano Zampini PetscFunctionBegin; 892302440fdSBarry Smith ierr = PCShellGetContext(fetipc,(void**)&pc_ctx);CHKERRQ(ierr); 893674ae819SStefano Zampini pcis = (PC_IS*)pc_ctx->pc->data; 894674ae819SStefano Zampini /* Application of B_Ddelta^T */ 895674ae819SStefano Zampini ierr = VecScatterBegin(pc_ctx->l2g_lambda,x,pc_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 896674ae819SStefano Zampini ierr = VecScatterEnd(pc_ctx->l2g_lambda,x,pc_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 897674ae819SStefano Zampini ierr = VecSet(pcis->vec2_B,0.0);CHKERRQ(ierr); 898674ae819SStefano Zampini ierr = MatMultTranspose(pc_ctx->B_Ddelta,pc_ctx->lambda_local,pcis->vec2_B);CHKERRQ(ierr); 899ed6c3d69SStefano Zampini /* Application of local Schur complement */ 900ed6c3d69SStefano Zampini ierr = MatMult(pc_ctx->S_j,pcis->vec2_B,pcis->vec1_B);CHKERRQ(ierr); 901edf7251bSStefano Zampini /* Application of B_Ddelta */ 902edf7251bSStefano Zampini ierr = MatMult(pc_ctx->B_Ddelta,pcis->vec1_B,pc_ctx->lambda_local);CHKERRQ(ierr); 903edf7251bSStefano Zampini ierr = VecSet(y,0.0);CHKERRQ(ierr); 904edf7251bSStefano Zampini ierr = VecScatterBegin(pc_ctx->l2g_lambda,pc_ctx->lambda_local,y,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 905edf7251bSStefano Zampini ierr = VecScatterEnd(pc_ctx->l2g_lambda,pc_ctx->lambda_local,y,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 906b519380fSstefano_zampini /* interface pressure preconditioner */ 907b519380fSstefano_zampini if (pc_ctx->kP) { 908b519380fSstefano_zampini const PetscScalar *lx; 909b519380fSstefano_zampini PetscScalar *ly; 910b519380fSstefano_zampini 911b519380fSstefano_zampini /* pressures ordered first in x and y */ 912b519380fSstefano_zampini ierr = VecGetArrayRead(x,&lx);CHKERRQ(ierr); 913b519380fSstefano_zampini ierr = VecGetArray(y,&ly);CHKERRQ(ierr); 914b519380fSstefano_zampini ierr = VecPlaceArray(pc_ctx->xPg,lx);CHKERRQ(ierr); 915b519380fSstefano_zampini ierr = VecPlaceArray(pc_ctx->yPg,ly);CHKERRQ(ierr); 916b519380fSstefano_zampini ierr = KSPSolve(pc_ctx->kP,pc_ctx->xPg,pc_ctx->yPg);CHKERRQ(ierr); 917b519380fSstefano_zampini ierr = VecResetArray(pc_ctx->xPg);CHKERRQ(ierr); 918b519380fSstefano_zampini ierr = VecResetArray(pc_ctx->yPg);CHKERRQ(ierr); 919b519380fSstefano_zampini ierr = VecRestoreArrayRead(x,&lx);CHKERRQ(ierr); 920b519380fSstefano_zampini ierr = VecRestoreArray(y,&ly);CHKERRQ(ierr); 921b519380fSstefano_zampini } 922edf7251bSStefano Zampini PetscFunctionReturn(0); 923edf7251bSStefano Zampini } 924edf7251bSStefano Zampini 925edf7251bSStefano Zampini #undef __FUNCT__ 926edf7251bSStefano Zampini #define __FUNCT__ "FETIDPPCApplyTranspose" 927edf7251bSStefano Zampini PetscErrorCode FETIDPPCApplyTranspose(PC fetipc, Vec x, Vec y) 928edf7251bSStefano Zampini { 929edf7251bSStefano Zampini FETIDPPC_ctx pc_ctx; 930edf7251bSStefano Zampini PC_IS *pcis; 931edf7251bSStefano Zampini PetscErrorCode ierr; 932edf7251bSStefano Zampini 933edf7251bSStefano Zampini PetscFunctionBegin; 934302440fdSBarry Smith ierr = PCShellGetContext(fetipc,(void**)&pc_ctx);CHKERRQ(ierr); 935edf7251bSStefano Zampini pcis = (PC_IS*)pc_ctx->pc->data; 936edf7251bSStefano Zampini /* Application of B_Ddelta^T */ 937edf7251bSStefano Zampini ierr = VecScatterBegin(pc_ctx->l2g_lambda,x,pc_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 938edf7251bSStefano Zampini ierr = VecScatterEnd(pc_ctx->l2g_lambda,x,pc_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 939edf7251bSStefano Zampini ierr = VecSet(pcis->vec2_B,0.0);CHKERRQ(ierr); 940edf7251bSStefano Zampini ierr = MatMultTranspose(pc_ctx->B_Ddelta,pc_ctx->lambda_local,pcis->vec2_B);CHKERRQ(ierr); 941ed6c3d69SStefano Zampini /* Application of local Schur complement */ 942ed6c3d69SStefano Zampini ierr = MatMultTranspose(pc_ctx->S_j,pcis->vec2_B,pcis->vec1_B);CHKERRQ(ierr); 943674ae819SStefano Zampini /* Application of B_Ddelta */ 944674ae819SStefano Zampini ierr = MatMult(pc_ctx->B_Ddelta,pcis->vec1_B,pc_ctx->lambda_local);CHKERRQ(ierr); 945674ae819SStefano Zampini ierr = VecSet(y,0.0);CHKERRQ(ierr); 946674ae819SStefano Zampini ierr = VecScatterBegin(pc_ctx->l2g_lambda,pc_ctx->lambda_local,y,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 947674ae819SStefano Zampini ierr = VecScatterEnd(pc_ctx->l2g_lambda,pc_ctx->lambda_local,y,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 948674ae819SStefano Zampini PetscFunctionReturn(0); 949674ae819SStefano Zampini } 950