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); 55457d4a33Sstefano_zampini ierr = VecDestroy(&mat_ctx->vP);CHKERRQ(ierr); 56457d4a33Sstefano_zampini ierr = VecDestroy(&mat_ctx->xPg);CHKERRQ(ierr); 57457d4a33Sstefano_zampini ierr = VecDestroy(&mat_ctx->yPg);CHKERRQ(ierr); 58674ae819SStefano Zampini ierr = VecScatterDestroy(&mat_ctx->l2g_lambda);CHKERRQ(ierr); 59457d4a33Sstefano_zampini ierr = VecScatterDestroy(&mat_ctx->l2g_p);CHKERRQ(ierr); 609cc7774eSstefano_zampini ierr = VecScatterDestroy(&mat_ctx->g2g_p);CHKERRQ(ierr); 61674ae819SStefano Zampini ierr = PCDestroy(&mat_ctx->pc);CHKERRQ(ierr); /* decrease PCBDDC reference count */ 62674ae819SStefano Zampini ierr = PetscFree(mat_ctx);CHKERRQ(ierr); 63674ae819SStefano Zampini PetscFunctionReturn(0); 64674ae819SStefano Zampini } 65674ae819SStefano Zampini 66674ae819SStefano Zampini #undef __FUNCT__ 67674ae819SStefano Zampini #define __FUNCT__ "PCBDDCDestroyFETIDPPC" 68674ae819SStefano Zampini PetscErrorCode PCBDDCDestroyFETIDPPC(PC pc) 69674ae819SStefano Zampini { 70674ae819SStefano Zampini FETIDPPC_ctx pc_ctx; 71674ae819SStefano Zampini PetscErrorCode ierr; 72674ae819SStefano Zampini 73674ae819SStefano Zampini PetscFunctionBegin; 74674ae819SStefano Zampini ierr = PCShellGetContext(pc,(void**)&pc_ctx);CHKERRQ(ierr); 75674ae819SStefano Zampini ierr = VecDestroy(&pc_ctx->lambda_local);CHKERRQ(ierr); 76674ae819SStefano Zampini ierr = MatDestroy(&pc_ctx->B_Ddelta);CHKERRQ(ierr); 77674ae819SStefano Zampini ierr = VecScatterDestroy(&pc_ctx->l2g_lambda);CHKERRQ(ierr); 7868668abeSStefano Zampini ierr = MatDestroy(&pc_ctx->S_j);CHKERRQ(ierr); 79674ae819SStefano Zampini ierr = PCDestroy(&pc_ctx->pc);CHKERRQ(ierr); /* decrease PCBDDC reference count */ 80457d4a33Sstefano_zampini ierr = KSPDestroy(&pc_ctx->kP);CHKERRQ(ierr); 81457d4a33Sstefano_zampini ierr = VecDestroy(&pc_ctx->xPg);CHKERRQ(ierr); 82457d4a33Sstefano_zampini ierr = VecDestroy(&pc_ctx->yPg);CHKERRQ(ierr); 83674ae819SStefano Zampini ierr = PetscFree(pc_ctx);CHKERRQ(ierr); 84674ae819SStefano Zampini PetscFunctionReturn(0); 85674ae819SStefano Zampini } 86674ae819SStefano Zampini 87674ae819SStefano Zampini #undef __FUNCT__ 88674ae819SStefano Zampini #define __FUNCT__ "PCBDDCSetupFETIDPMatContext" 89674ae819SStefano Zampini PetscErrorCode PCBDDCSetupFETIDPMatContext(FETIDPMat_ctx fetidpmat_ctx ) 90674ae819SStefano Zampini { 91674ae819SStefano Zampini PetscErrorCode ierr; 92674ae819SStefano Zampini PC_IS *pcis=(PC_IS*)fetidpmat_ctx->pc->data; 93674ae819SStefano Zampini PC_BDDC *pcbddc=(PC_BDDC*)fetidpmat_ctx->pc->data; 94674ae819SStefano Zampini PCBDDCGraph mat_graph=pcbddc->mat_graph; 95674ae819SStefano Zampini Mat_IS *matis = (Mat_IS*)fetidpmat_ctx->pc->pmat->data; 96674ae819SStefano Zampini MPI_Comm comm; 97674ae819SStefano Zampini Mat ScalingMat; 98457d4a33Sstefano_zampini Vec fetidp_global; 99674ae819SStefano Zampini IS IS_l2g_lambda; 100dc456d91SStefano Zampini IS subset,subset_mult,subset_n; 101674ae819SStefano Zampini PetscBool skip_node,fully_redundant; 102674ae819SStefano Zampini PetscInt i,j,k,s,n_boundary_dofs,n_global_lambda,n_vertices,partial_sum; 103dc456d91SStefano Zampini PetscInt cum,n_local_lambda,n_lambda_for_dof,dual_size,n_neg_values,n_pos_values; 10476ec1555SBarry Smith PetscMPIInt rank,size,buf_size,neigh; 105674ae819SStefano Zampini PetscScalar scalar_value; 106674ae819SStefano Zampini PetscInt *vertex_indices; 107dc456d91SStefano Zampini PetscInt *dual_dofs_boundary_indices,*aux_local_numbering_1; 108dc456d91SStefano Zampini const PetscInt *aux_global_numbering; 109674ae819SStefano Zampini PetscInt *aux_sums,*cols_B_delta,*l2g_indices; 110674ae819SStefano Zampini PetscScalar *array,*scaling_factors,*vals_B_delta; 111674ae819SStefano Zampini PetscInt *aux_local_numbering_2; 112457d4a33Sstefano_zampini PetscLayout llay; 113674ae819SStefano Zampini /* For communication of scaling factors */ 114674ae819SStefano Zampini PetscInt *ptrs_buffer,neigh_position; 115674ae819SStefano Zampini PetscScalar **all_factors,*send_buffer,*recv_buffer; 116674ae819SStefano Zampini MPI_Request *send_reqs,*recv_reqs; 117674ae819SStefano Zampini /* tests */ 118674ae819SStefano Zampini PetscBool test_fetidp; 119674ae819SStefano Zampini PetscViewer viewer; 120457d4a33Sstefano_zampini /* saddlepoint */ 121457d4a33Sstefano_zampini ISLocalToGlobalMapping l2gmap_p; 122457d4a33Sstefano_zampini PetscLayout play; 123457d4a33Sstefano_zampini IS gP,pP; 124457d4a33Sstefano_zampini PetscInt nPl,nPg,nPgl; 125674ae819SStefano Zampini 126674ae819SStefano Zampini PetscFunctionBegin; 127674ae819SStefano Zampini ierr = PetscObjectGetComm((PetscObject)(fetidpmat_ctx->pc),&comm);CHKERRQ(ierr); 128674ae819SStefano Zampini ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 12976ec1555SBarry Smith ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 130674ae819SStefano Zampini 131457d4a33Sstefano_zampini /* saddlepoint */ 132457d4a33Sstefano_zampini nPl = 0; 133457d4a33Sstefano_zampini nPg = 0; 134457d4a33Sstefano_zampini nPgl = 0; 135457d4a33Sstefano_zampini gP = NULL; 136457d4a33Sstefano_zampini pP = NULL; 137457d4a33Sstefano_zampini l2gmap_p = NULL; 138457d4a33Sstefano_zampini play = NULL; 139457d4a33Sstefano_zampini ierr = PetscObjectQuery((PetscObject)fetidpmat_ctx->pc,"__KSPFETIDP_pP",(PetscObject*)&pP);CHKERRQ(ierr); 140*022d8d2bSstefano_zampini if (pP) { /* saddle point */ 141457d4a33Sstefano_zampini /* subdomain pressures in global numbering */ 142457d4a33Sstefano_zampini ierr = PetscObjectQuery((PetscObject)fetidpmat_ctx->pc,"__KSPFETIDP_gP",(PetscObject*)&gP);CHKERRQ(ierr); 143457d4a33Sstefano_zampini if (!gP) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"gP not present"); 144457d4a33Sstefano_zampini ierr = ISGetLocalSize(gP,&nPl);CHKERRQ(ierr); 145457d4a33Sstefano_zampini ierr = VecCreate(PETSC_COMM_SELF,&fetidpmat_ctx->vP);CHKERRQ(ierr); 146457d4a33Sstefano_zampini ierr = VecSetSizes(fetidpmat_ctx->vP,nPl,nPl);CHKERRQ(ierr); 147457d4a33Sstefano_zampini ierr = VecSetType(fetidpmat_ctx->vP,VECSTANDARD);CHKERRQ(ierr); 148457d4a33Sstefano_zampini ierr = VecSetUp(fetidpmat_ctx->vP);CHKERRQ(ierr); 149457d4a33Sstefano_zampini 150457d4a33Sstefano_zampini /* interface pressure matrix */ 151457d4a33Sstefano_zampini ierr = PetscObjectQuery((PetscObject)fetidpmat_ctx->pc,"__KSPFETIDP_C",(PetscObject*)&fetidpmat_ctx->C);CHKERRQ(ierr); 152457d4a33Sstefano_zampini if (!fetidpmat_ctx->C) { /* null pressure block, compute layout and global numbering for interface pressures */ 153457d4a33Sstefano_zampini IS Pg; 154457d4a33Sstefano_zampini 155457d4a33Sstefano_zampini ierr = ISRenumber(gP,NULL,&nPg,&Pg);CHKERRQ(ierr); 156457d4a33Sstefano_zampini ierr = ISLocalToGlobalMappingCreateIS(Pg,&l2gmap_p);CHKERRQ(ierr); 157457d4a33Sstefano_zampini ierr = ISDestroy(&Pg);CHKERRQ(ierr); 158457d4a33Sstefano_zampini ierr = PetscLayoutCreate(comm,&play);CHKERRQ(ierr); 159457d4a33Sstefano_zampini ierr = PetscLayoutSetBlockSize(play,1);CHKERRQ(ierr); 160457d4a33Sstefano_zampini ierr = PetscLayoutSetSize(play,nPg);CHKERRQ(ierr); 161457d4a33Sstefano_zampini ierr = ISGetLocalSize(pP,&nPgl);CHKERRQ(ierr); 162457d4a33Sstefano_zampini ierr = PetscLayoutSetLocalSize(play,nPgl);CHKERRQ(ierr); 163457d4a33Sstefano_zampini ierr = PetscLayoutSetUp(play);CHKERRQ(ierr); 164457d4a33Sstefano_zampini } else { 165457d4a33Sstefano_zampini ierr = PetscObjectReference((PetscObject)fetidpmat_ctx->C);CHKERRQ(ierr); 166457d4a33Sstefano_zampini ierr = MatGetLocalToGlobalMapping(fetidpmat_ctx->C,&l2gmap_p,NULL);CHKERRQ(ierr); 167457d4a33Sstefano_zampini ierr = PetscObjectReference((PetscObject)l2gmap_p);CHKERRQ(ierr); 168457d4a33Sstefano_zampini ierr = MatGetSize(fetidpmat_ctx->C,&nPg,NULL);CHKERRQ(ierr); 169457d4a33Sstefano_zampini ierr = MatGetLocalSize(fetidpmat_ctx->C,NULL,&nPgl);CHKERRQ(ierr); 170457d4a33Sstefano_zampini ierr = MatGetLayouts(fetidpmat_ctx->C,NULL,&llay);CHKERRQ(ierr); 171457d4a33Sstefano_zampini ierr = PetscLayoutReference(llay,&play);CHKERRQ(ierr); 172457d4a33Sstefano_zampini } 173457d4a33Sstefano_zampini ierr = VecCreateMPIWithArray(comm,1,nPgl,nPg,NULL,&fetidpmat_ctx->xPg);CHKERRQ(ierr); 174457d4a33Sstefano_zampini ierr = VecCreateMPIWithArray(comm,1,nPgl,nPg,NULL,&fetidpmat_ctx->yPg);CHKERRQ(ierr); 175457d4a33Sstefano_zampini 176457d4a33Sstefano_zampini /* import matrices for interface pressures coupling */ 177457d4a33Sstefano_zampini ierr = PetscObjectQuery((PetscObject)fetidpmat_ctx->pc,"__KSPFETIDP_B_BI",(PetscObject*)&fetidpmat_ctx->B_BI);CHKERRQ(ierr); 178457d4a33Sstefano_zampini if (!fetidpmat_ctx->B_BI) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"B_BI not present"); 179457d4a33Sstefano_zampini ierr = PetscObjectReference((PetscObject)fetidpmat_ctx->B_BI);CHKERRQ(ierr); 180457d4a33Sstefano_zampini 181457d4a33Sstefano_zampini ierr = PetscObjectQuery((PetscObject)fetidpmat_ctx->pc,"__KSPFETIDP_B_BB",(PetscObject*)&fetidpmat_ctx->B_BB);CHKERRQ(ierr); 182457d4a33Sstefano_zampini if (!fetidpmat_ctx->B_BB) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"B_BB not present"); 183457d4a33Sstefano_zampini ierr = PetscObjectReference((PetscObject)fetidpmat_ctx->B_BB);CHKERRQ(ierr); 184457d4a33Sstefano_zampini 185457d4a33Sstefano_zampini ierr = PetscObjectQuery((PetscObject)fetidpmat_ctx->pc,"__KSPFETIDP_Bt_BI",(PetscObject*)&fetidpmat_ctx->Bt_BI);CHKERRQ(ierr); 186457d4a33Sstefano_zampini if (!fetidpmat_ctx->Bt_BI) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Bt_BI not present"); 187457d4a33Sstefano_zampini ierr = PetscObjectReference((PetscObject)fetidpmat_ctx->Bt_BI);CHKERRQ(ierr); 188457d4a33Sstefano_zampini 189457d4a33Sstefano_zampini ierr = PetscObjectQuery((PetscObject)fetidpmat_ctx->pc,"__KSPFETIDP_Bt_BB",(PetscObject*)&fetidpmat_ctx->Bt_BB);CHKERRQ(ierr); 190457d4a33Sstefano_zampini if (!fetidpmat_ctx->Bt_BB) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Bt_BB not present"); 191457d4a33Sstefano_zampini ierr = PetscObjectReference((PetscObject)fetidpmat_ctx->Bt_BB);CHKERRQ(ierr); 192457d4a33Sstefano_zampini } 193457d4a33Sstefano_zampini 194674ae819SStefano Zampini /* Default type of lagrange multipliers is non-redundant */ 195329cd39dSStefano Zampini fully_redundant = fetidpmat_ctx->fully_redundant; 196674ae819SStefano Zampini 197674ae819SStefano Zampini /* Evaluate local and global number of lagrange multipliers */ 198674ae819SStefano Zampini ierr = VecSet(pcis->vec1_N,0.0);CHKERRQ(ierr); 199674ae819SStefano Zampini n_local_lambda = 0; 200674ae819SStefano Zampini partial_sum = 0; 201674ae819SStefano Zampini n_boundary_dofs = 0; 202674ae819SStefano Zampini s = 0; 203674ae819SStefano Zampini /* Get Vertices used to define the BDDC */ 204b371cd4fSStefano Zampini n_vertices = pcbddc->n_vertices; 2050e6343abSStefano Zampini vertex_indices = pcbddc->local_primal_ref_node; 206674ae819SStefano Zampini dual_size = pcis->n_B-n_vertices; 207785e854fSJed Brown ierr = PetscMalloc1(dual_size,&dual_dofs_boundary_indices);CHKERRQ(ierr); 208785e854fSJed Brown ierr = PetscMalloc1(dual_size,&aux_local_numbering_1);CHKERRQ(ierr); 209785e854fSJed Brown ierr = PetscMalloc1(dual_size,&aux_local_numbering_2);CHKERRQ(ierr); 210674ae819SStefano Zampini 211674ae819SStefano Zampini ierr = VecGetArray(pcis->vec1_N,&array);CHKERRQ(ierr); 212674ae819SStefano Zampini for (i=0;i<pcis->n;i++){ 213674ae819SStefano Zampini j = mat_graph->count[i]; /* RECALL: mat_graph->count[i] does not count myself */ 214674ae819SStefano Zampini if ( j > 0 ) { 215674ae819SStefano Zampini n_boundary_dofs++; 216674ae819SStefano Zampini } 217674ae819SStefano Zampini skip_node = PETSC_FALSE; 218674ae819SStefano Zampini if ( s < n_vertices && vertex_indices[s]==i) { /* it works for a sorted set of vertices */ 219674ae819SStefano Zampini skip_node = PETSC_TRUE; 220674ae819SStefano Zampini s++; 221674ae819SStefano Zampini } 222674ae819SStefano Zampini if (j < 1) { 223674ae819SStefano Zampini skip_node = PETSC_TRUE; 224674ae819SStefano Zampini } 225674ae819SStefano Zampini if ( !skip_node ) { 226674ae819SStefano Zampini if (fully_redundant) { 227674ae819SStefano Zampini /* fully redundant set of lagrange multipliers */ 228674ae819SStefano Zampini n_lambda_for_dof = (j*(j+1))/2; 229674ae819SStefano Zampini } else { 230674ae819SStefano Zampini n_lambda_for_dof = j; 231674ae819SStefano Zampini } 232674ae819SStefano Zampini n_local_lambda += j; 233674ae819SStefano Zampini /* needed to evaluate global number of lagrange multipliers */ 234674ae819SStefano Zampini array[i]=(1.0*n_lambda_for_dof)/(j+1.0); /* already scaled for the next global sum */ 235674ae819SStefano Zampini /* store some data needed */ 236674ae819SStefano Zampini dual_dofs_boundary_indices[partial_sum] = n_boundary_dofs-1; 237674ae819SStefano Zampini aux_local_numbering_1[partial_sum] = i; 238674ae819SStefano Zampini aux_local_numbering_2[partial_sum] = n_lambda_for_dof; 239674ae819SStefano Zampini partial_sum++; 240674ae819SStefano Zampini } 241674ae819SStefano Zampini } 242674ae819SStefano Zampini ierr = VecRestoreArray(pcis->vec1_N,&array);CHKERRQ(ierr); 243674ae819SStefano Zampini 244674ae819SStefano Zampini /* compute global ordering of lagrange multipliers and associate l2g map */ 245dc456d91SStefano Zampini ierr = ISCreateGeneral(comm,partial_sum,aux_local_numbering_1,PETSC_COPY_VALUES,&subset_n);CHKERRQ(ierr); 2463bbff08aSStefano Zampini ierr = ISLocalToGlobalMappingApplyIS(pcis->mapping,subset_n,&subset);CHKERRQ(ierr); 247dc456d91SStefano Zampini ierr = ISDestroy(&subset_n);CHKERRQ(ierr); 248dc456d91SStefano Zampini ierr = ISCreateGeneral(comm,partial_sum,aux_local_numbering_2,PETSC_OWN_POINTER,&subset_mult);CHKERRQ(ierr); 2493d996552SStefano Zampini ierr = ISRenumber(subset,subset_mult,&fetidpmat_ctx->n_lambda,&subset_n);CHKERRQ(ierr); 250dc456d91SStefano Zampini ierr = ISDestroy(&subset);CHKERRQ(ierr); 2513d996552SStefano Zampini 2524fe826edSStefano Zampini #if defined(PETSC_USE_DEBUG) 2534fe826edSStefano Zampini ierr = VecSet(pcis->vec1_global,0.0);CHKERRQ(ierr); 2544fe826edSStefano Zampini ierr = VecScatterBegin(matis->rctx,pcis->vec1_N,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 2554fe826edSStefano Zampini ierr = VecScatterEnd(matis->rctx,pcis->vec1_N,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 2564fe826edSStefano Zampini ierr = VecSum(pcis->vec1_global,&scalar_value);CHKERRQ(ierr); 2574fe826edSStefano Zampini i = (PetscInt)PetscRealPart(scalar_value); 2586c4ed002SBarry 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); 2594fe826edSStefano Zampini #endif 260674ae819SStefano Zampini 261674ae819SStefano Zampini /* init data for scaling factors exchange */ 262674ae819SStefano Zampini partial_sum = 0; 263785e854fSJed Brown ierr = PetscMalloc1(pcis->n_neigh,&ptrs_buffer);CHKERRQ(ierr); 2644b2aedd3SStefano Zampini ierr = PetscMalloc1(PetscMax(pcis->n_neigh-1,0),&send_reqs);CHKERRQ(ierr); 2654b2aedd3SStefano Zampini ierr = PetscMalloc1(PetscMax(pcis->n_neigh-1,0),&recv_reqs);CHKERRQ(ierr); 266785e854fSJed Brown ierr = PetscMalloc1(pcis->n,&all_factors);CHKERRQ(ierr); 2674b2aedd3SStefano Zampini if (pcis->n_neigh > 0) ptrs_buffer[0]=0; 268674ae819SStefano Zampini for (i=1;i<pcis->n_neigh;i++) { 269674ae819SStefano Zampini partial_sum += pcis->n_shared[i]; 270674ae819SStefano Zampini ptrs_buffer[i] = ptrs_buffer[i-1]+pcis->n_shared[i]; 271674ae819SStefano Zampini } 272785e854fSJed Brown ierr = PetscMalloc1(partial_sum,&send_buffer);CHKERRQ(ierr); 273785e854fSJed Brown ierr = PetscMalloc1(partial_sum,&recv_buffer);CHKERRQ(ierr); 274785e854fSJed Brown ierr = PetscMalloc1(partial_sum,&all_factors[0]);CHKERRQ(ierr); 275674ae819SStefano Zampini for (i=0;i<pcis->n-1;i++) { 276674ae819SStefano Zampini j = mat_graph->count[i]; 277674ae819SStefano Zampini all_factors[i+1]=all_factors[i]+j; 278674ae819SStefano Zampini } 279674ae819SStefano Zampini /* scatter B scaling to N vec */ 280674ae819SStefano Zampini ierr = VecScatterBegin(pcis->N_to_B,pcis->D,pcis->vec1_N,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 281674ae819SStefano Zampini ierr = VecScatterEnd(pcis->N_to_B,pcis->D,pcis->vec1_N,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 282674ae819SStefano Zampini /* communications */ 2832b095fd8SStefano Zampini ierr = VecGetArrayRead(pcis->vec1_N,(const PetscScalar**)&array);CHKERRQ(ierr); 284674ae819SStefano Zampini for (i=1;i<pcis->n_neigh;i++) { 285674ae819SStefano Zampini for (j=0;j<pcis->n_shared[i];j++) { 286674ae819SStefano Zampini send_buffer[ptrs_buffer[i-1]+j]=array[pcis->shared[i][j]]; 287674ae819SStefano Zampini } 288674ae819SStefano Zampini ierr = PetscMPIIntCast(ptrs_buffer[i]-ptrs_buffer[i-1],&buf_size);CHKERRQ(ierr); 289674ae819SStefano Zampini ierr = PetscMPIIntCast(pcis->neigh[i],&neigh);CHKERRQ(ierr); 290674ae819SStefano Zampini ierr = MPI_Isend(&send_buffer[ptrs_buffer[i-1]],buf_size,MPIU_SCALAR,neigh,0,comm,&send_reqs[i-1]);CHKERRQ(ierr); 291674ae819SStefano Zampini ierr = MPI_Irecv(&recv_buffer[ptrs_buffer[i-1]],buf_size,MPIU_SCALAR,neigh,0,comm,&recv_reqs[i-1]);CHKERRQ(ierr); 292674ae819SStefano Zampini } 2932b095fd8SStefano Zampini ierr = VecRestoreArrayRead(pcis->vec1_N,(const PetscScalar**)&array);CHKERRQ(ierr); 2944b2aedd3SStefano Zampini if (pcis->n_neigh > 0) { 2954b2aedd3SStefano Zampini ierr = MPI_Waitall(pcis->n_neigh-1,recv_reqs,MPI_STATUSES_IGNORE);CHKERRQ(ierr); 2964b2aedd3SStefano Zampini } 297674ae819SStefano Zampini /* put values in correct places */ 298674ae819SStefano Zampini for (i=1;i<pcis->n_neigh;i++) { 299674ae819SStefano Zampini for (j=0;j<pcis->n_shared[i];j++) { 300674ae819SStefano Zampini k = pcis->shared[i][j]; 301674ae819SStefano Zampini neigh_position = 0; 302674ae819SStefano Zampini while(mat_graph->neighbours_set[k][neigh_position] != pcis->neigh[i]) {neigh_position++;} 303674ae819SStefano Zampini all_factors[k][neigh_position]=recv_buffer[ptrs_buffer[i-1]+j]; 304674ae819SStefano Zampini } 305674ae819SStefano Zampini } 3064b2aedd3SStefano Zampini if (pcis->n_neigh > 0) { 3074b2aedd3SStefano Zampini ierr = MPI_Waitall(pcis->n_neigh-1,send_reqs,MPI_STATUSES_IGNORE);CHKERRQ(ierr); 3084b2aedd3SStefano Zampini } 309674ae819SStefano Zampini ierr = PetscFree(send_reqs);CHKERRQ(ierr); 310674ae819SStefano Zampini ierr = PetscFree(recv_reqs);CHKERRQ(ierr); 311674ae819SStefano Zampini ierr = PetscFree(send_buffer);CHKERRQ(ierr); 312674ae819SStefano Zampini ierr = PetscFree(recv_buffer);CHKERRQ(ierr); 313674ae819SStefano Zampini ierr = PetscFree(ptrs_buffer);CHKERRQ(ierr); 314674ae819SStefano Zampini 315674ae819SStefano Zampini /* Compute B and B_delta (local actions) */ 316785e854fSJed Brown ierr = PetscMalloc1(pcis->n_neigh,&aux_sums);CHKERRQ(ierr); 317785e854fSJed Brown ierr = PetscMalloc1(n_local_lambda,&l2g_indices);CHKERRQ(ierr); 318785e854fSJed Brown ierr = PetscMalloc1(n_local_lambda,&vals_B_delta);CHKERRQ(ierr); 319785e854fSJed Brown ierr = PetscMalloc1(n_local_lambda,&cols_B_delta);CHKERRQ(ierr); 320785e854fSJed Brown ierr = PetscMalloc1(n_local_lambda,&scaling_factors);CHKERRQ(ierr); 321dc456d91SStefano Zampini ierr = ISGetIndices(subset_n,&aux_global_numbering);CHKERRQ(ierr); 322674ae819SStefano Zampini partial_sum=0; 323dc456d91SStefano Zampini cum = 0; 324674ae819SStefano Zampini for (i=0;i<dual_size;i++) { 325dc456d91SStefano Zampini n_global_lambda = aux_global_numbering[cum]; 326674ae819SStefano Zampini j = mat_graph->count[aux_local_numbering_1[i]]; 327674ae819SStefano Zampini aux_sums[0]=0; 328674ae819SStefano Zampini for (s=1;s<j;s++) { 329674ae819SStefano Zampini aux_sums[s]=aux_sums[s-1]+j-s+1; 330674ae819SStefano Zampini } 331674ae819SStefano Zampini array = all_factors[aux_local_numbering_1[i]]; 332674ae819SStefano Zampini n_neg_values = 0; 3332a7da448SStefano Zampini while(n_neg_values < j && mat_graph->neighbours_set[aux_local_numbering_1[i]][n_neg_values] < rank) {n_neg_values++;} 334674ae819SStefano Zampini n_pos_values = j - n_neg_values; 335674ae819SStefano Zampini if (fully_redundant) { 336674ae819SStefano Zampini for (s=0;s<n_neg_values;s++) { 337674ae819SStefano Zampini l2g_indices [partial_sum+s]=aux_sums[s]+n_neg_values-s-1+n_global_lambda; 338674ae819SStefano Zampini cols_B_delta [partial_sum+s]=dual_dofs_boundary_indices[i]; 339674ae819SStefano Zampini vals_B_delta [partial_sum+s]=-1.0; 340674ae819SStefano Zampini scaling_factors[partial_sum+s]=array[s]; 341674ae819SStefano Zampini } 342674ae819SStefano Zampini for (s=0;s<n_pos_values;s++) { 343674ae819SStefano Zampini l2g_indices [partial_sum+s+n_neg_values]=aux_sums[n_neg_values]+s+n_global_lambda; 344674ae819SStefano Zampini cols_B_delta [partial_sum+s+n_neg_values]=dual_dofs_boundary_indices[i]; 345674ae819SStefano Zampini vals_B_delta [partial_sum+s+n_neg_values]=1.0; 346674ae819SStefano Zampini scaling_factors[partial_sum+s+n_neg_values]=array[s+n_neg_values]; 347674ae819SStefano Zampini } 348674ae819SStefano Zampini partial_sum += j; 349674ae819SStefano Zampini } else { 350674ae819SStefano Zampini /* l2g_indices and default cols and vals of B_delta */ 351674ae819SStefano Zampini for (s=0;s<j;s++) { 352674ae819SStefano Zampini l2g_indices [partial_sum+s]=n_global_lambda+s; 353674ae819SStefano Zampini cols_B_delta [partial_sum+s]=dual_dofs_boundary_indices[i]; 354674ae819SStefano Zampini vals_B_delta [partial_sum+s]=0.0; 355674ae819SStefano Zampini } 356674ae819SStefano Zampini /* B_delta */ 357674ae819SStefano Zampini if ( n_neg_values > 0 ) { /* there's a rank next to me to the left */ 358674ae819SStefano Zampini vals_B_delta [partial_sum+n_neg_values-1]=-1.0; 359674ae819SStefano Zampini } 360674ae819SStefano Zampini if ( n_neg_values < j ) { /* there's a rank next to me to the right */ 361674ae819SStefano Zampini vals_B_delta [partial_sum+n_neg_values]=1.0; 362674ae819SStefano Zampini } 363674ae819SStefano Zampini /* scaling as in Klawonn-Widlund 1999*/ 364674ae819SStefano Zampini for (s=0;s<n_neg_values;s++) { 365674ae819SStefano Zampini scalar_value = 0.0; 366674ae819SStefano Zampini for (k=0;k<s+1;k++) { 367674ae819SStefano Zampini scalar_value += array[k]; 368674ae819SStefano Zampini } 369674ae819SStefano Zampini scaling_factors[partial_sum+s] = -scalar_value; 370674ae819SStefano Zampini } 371674ae819SStefano Zampini for (s=0;s<n_pos_values;s++) { 372674ae819SStefano Zampini scalar_value = 0.0; 373674ae819SStefano Zampini for (k=s+n_neg_values;k<j;k++) { 374674ae819SStefano Zampini scalar_value += array[k]; 375674ae819SStefano Zampini } 376674ae819SStefano Zampini scaling_factors[partial_sum+s+n_neg_values] = scalar_value; 377674ae819SStefano Zampini } 378674ae819SStefano Zampini partial_sum += j; 379674ae819SStefano Zampini } 380dc456d91SStefano Zampini cum += aux_local_numbering_2[i]; 381674ae819SStefano Zampini } 382dc456d91SStefano Zampini ierr = ISRestoreIndices(subset_n,&aux_global_numbering);CHKERRQ(ierr); 383dc456d91SStefano Zampini ierr = ISDestroy(&subset_mult);CHKERRQ(ierr); 384dc456d91SStefano Zampini ierr = ISDestroy(&subset_n);CHKERRQ(ierr); 385674ae819SStefano Zampini ierr = PetscFree(aux_sums);CHKERRQ(ierr); 386674ae819SStefano Zampini ierr = PetscFree(aux_local_numbering_1);CHKERRQ(ierr); 387674ae819SStefano Zampini ierr = PetscFree(dual_dofs_boundary_indices);CHKERRQ(ierr); 388674ae819SStefano Zampini ierr = PetscFree(all_factors[0]);CHKERRQ(ierr); 389674ae819SStefano Zampini ierr = PetscFree(all_factors);CHKERRQ(ierr); 390674ae819SStefano Zampini 391674ae819SStefano Zampini /* Create local part of B_delta */ 392302440fdSBarry Smith ierr = MatCreate(PETSC_COMM_SELF,&fetidpmat_ctx->B_delta);CHKERRQ(ierr); 393674ae819SStefano Zampini ierr = MatSetSizes(fetidpmat_ctx->B_delta,n_local_lambda,pcis->n_B,n_local_lambda,pcis->n_B);CHKERRQ(ierr); 394674ae819SStefano Zampini ierr = MatSetType(fetidpmat_ctx->B_delta,MATSEQAIJ);CHKERRQ(ierr); 395674ae819SStefano Zampini ierr = MatSeqAIJSetPreallocation(fetidpmat_ctx->B_delta,1,NULL);CHKERRQ(ierr); 396674ae819SStefano Zampini ierr = MatSetOption(fetidpmat_ctx->B_delta,MAT_IGNORE_ZERO_ENTRIES,PETSC_TRUE);CHKERRQ(ierr); 397674ae819SStefano Zampini for (i=0;i<n_local_lambda;i++) { 398674ae819SStefano Zampini ierr = MatSetValue(fetidpmat_ctx->B_delta,i,cols_B_delta[i],vals_B_delta[i],INSERT_VALUES);CHKERRQ(ierr); 399674ae819SStefano Zampini } 400674ae819SStefano Zampini ierr = PetscFree(vals_B_delta);CHKERRQ(ierr); 401674ae819SStefano Zampini ierr = MatAssemblyBegin(fetidpmat_ctx->B_delta,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 402674ae819SStefano Zampini ierr = MatAssemblyEnd (fetidpmat_ctx->B_delta,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 403674ae819SStefano Zampini 404674ae819SStefano Zampini if (fully_redundant) { 405302440fdSBarry Smith ierr = MatCreate(PETSC_COMM_SELF,&ScalingMat);CHKERRQ(ierr); 406674ae819SStefano Zampini ierr = MatSetSizes(ScalingMat,n_local_lambda,n_local_lambda,n_local_lambda,n_local_lambda);CHKERRQ(ierr); 407674ae819SStefano Zampini ierr = MatSetType(ScalingMat,MATSEQAIJ);CHKERRQ(ierr); 408674ae819SStefano Zampini ierr = MatSeqAIJSetPreallocation(ScalingMat,1,NULL);CHKERRQ(ierr); 409674ae819SStefano Zampini for (i=0;i<n_local_lambda;i++) { 410674ae819SStefano Zampini ierr = MatSetValue(ScalingMat,i,i,scaling_factors[i],INSERT_VALUES);CHKERRQ(ierr); 411674ae819SStefano Zampini } 412674ae819SStefano Zampini ierr = MatAssemblyBegin(ScalingMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 413674ae819SStefano Zampini ierr = MatAssemblyEnd (ScalingMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 414674ae819SStefano Zampini ierr = MatMatMult(ScalingMat,fetidpmat_ctx->B_delta,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&fetidpmat_ctx->B_Ddelta);CHKERRQ(ierr); 415674ae819SStefano Zampini ierr = MatDestroy(&ScalingMat);CHKERRQ(ierr); 416674ae819SStefano Zampini } else { 417302440fdSBarry Smith ierr = MatCreate(PETSC_COMM_SELF,&fetidpmat_ctx->B_Ddelta);CHKERRQ(ierr); 418674ae819SStefano Zampini ierr = MatSetSizes(fetidpmat_ctx->B_Ddelta,n_local_lambda,pcis->n_B,n_local_lambda,pcis->n_B);CHKERRQ(ierr); 419674ae819SStefano Zampini ierr = MatSetType(fetidpmat_ctx->B_Ddelta,MATSEQAIJ);CHKERRQ(ierr); 420674ae819SStefano Zampini ierr = MatSeqAIJSetPreallocation(fetidpmat_ctx->B_Ddelta,1,NULL);CHKERRQ(ierr); 421674ae819SStefano Zampini for (i=0;i<n_local_lambda;i++) { 422674ae819SStefano Zampini ierr = MatSetValue(fetidpmat_ctx->B_Ddelta,i,cols_B_delta[i],scaling_factors[i],INSERT_VALUES);CHKERRQ(ierr); 423674ae819SStefano Zampini } 424674ae819SStefano Zampini ierr = MatAssemblyBegin(fetidpmat_ctx->B_Ddelta,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 425674ae819SStefano Zampini ierr = MatAssemblyEnd (fetidpmat_ctx->B_Ddelta,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 426674ae819SStefano Zampini } 427674ae819SStefano Zampini ierr = PetscFree(scaling_factors);CHKERRQ(ierr); 428674ae819SStefano Zampini ierr = PetscFree(cols_B_delta);CHKERRQ(ierr); 429674ae819SStefano Zampini 430457d4a33Sstefano_zampini /* Layout of multipliers */ 431457d4a33Sstefano_zampini ierr = PetscLayoutCreate(comm,&llay);CHKERRQ(ierr); 432457d4a33Sstefano_zampini ierr = PetscLayoutSetBlockSize(llay,1);CHKERRQ(ierr); 433457d4a33Sstefano_zampini ierr = PetscLayoutSetSize(llay,fetidpmat_ctx->n_lambda);CHKERRQ(ierr); 434457d4a33Sstefano_zampini ierr = PetscLayoutSetUp(llay);CHKERRQ(ierr); 435457d4a33Sstefano_zampini ierr = PetscLayoutGetLocalSize(llay,&fetidpmat_ctx->n);CHKERRQ(ierr); 436457d4a33Sstefano_zampini 437457d4a33Sstefano_zampini /* fetidpmat sizes */ 438457d4a33Sstefano_zampini fetidpmat_ctx->n += nPgl; 439457d4a33Sstefano_zampini fetidpmat_ctx->N = fetidpmat_ctx->n_lambda+nPg; 440457d4a33Sstefano_zampini 441457d4a33Sstefano_zampini /* Local work vector of multipliers */ 442457d4a33Sstefano_zampini ierr = VecCreate(PETSC_COMM_SELF,&fetidpmat_ctx->lambda_local);CHKERRQ(ierr); 443457d4a33Sstefano_zampini ierr = VecSetSizes(fetidpmat_ctx->lambda_local,n_local_lambda,n_local_lambda);CHKERRQ(ierr); 444457d4a33Sstefano_zampini ierr = VecSetType(fetidpmat_ctx->lambda_local,VECSEQ);CHKERRQ(ierr); 445457d4a33Sstefano_zampini 446457d4a33Sstefano_zampini /* Global vector for FETI-DP linear system */ 447457d4a33Sstefano_zampini ierr = VecCreate(comm,&fetidp_global);CHKERRQ(ierr); 448457d4a33Sstefano_zampini ierr = VecSetSizes(fetidp_global,fetidpmat_ctx->n,fetidpmat_ctx->N);CHKERRQ(ierr); 449457d4a33Sstefano_zampini ierr = VecSetType(fetidp_global,VECMPI);CHKERRQ(ierr); 450457d4a33Sstefano_zampini ierr = VecSetUp(fetidp_global);CHKERRQ(ierr); 451457d4a33Sstefano_zampini 4529eec4de8Sstefano_zampini /* Decide layout for fetidp dofs: if it is a saddle point problem 4539eec4de8Sstefano_zampini pressure is ordered first in the local part of the global vector 454457d4a33Sstefano_zampini of the FETI-DP linear system */ 455457d4a33Sstefano_zampini if (nPg) { 456af140850Sstefano_zampini IS IS_l2g_p,ais; 457457d4a33Sstefano_zampini PetscLayout alay; 458457d4a33Sstefano_zampini const PetscInt *idxs,*pranges,*aranges,*lranges; 459af140850Sstefano_zampini PetscInt *l2g_indices_p,rst; 460457d4a33Sstefano_zampini 461457d4a33Sstefano_zampini ierr = PetscMalloc1(nPl,&l2g_indices_p);CHKERRQ(ierr); 462457d4a33Sstefano_zampini ierr = VecGetLayout(fetidp_global,&alay);CHKERRQ(ierr); 463457d4a33Sstefano_zampini ierr = PetscLayoutGetRanges(alay,&aranges);CHKERRQ(ierr); 464457d4a33Sstefano_zampini ierr = PetscLayoutGetRanges(play,&pranges);CHKERRQ(ierr); 465457d4a33Sstefano_zampini ierr = PetscLayoutGetRanges(llay,&lranges);CHKERRQ(ierr); 466457d4a33Sstefano_zampini ierr = ISLocalToGlobalMappingGetIndices(l2gmap_p,&idxs);CHKERRQ(ierr); 467af140850Sstefano_zampini /* shift local to global indices for pressure */ 468457d4a33Sstefano_zampini for (i=0;i<nPl;i++) { 469457d4a33Sstefano_zampini PetscInt owner; 470457d4a33Sstefano_zampini 471457d4a33Sstefano_zampini ierr = PetscLayoutFindOwner(play,idxs[i],&owner);CHKERRQ(ierr); 472457d4a33Sstefano_zampini l2g_indices_p[i] = idxs[i]-pranges[owner]+aranges[owner]; 473457d4a33Sstefano_zampini } 474457d4a33Sstefano_zampini ierr = ISLocalToGlobalMappingRestoreIndices(l2gmap_p,&idxs);CHKERRQ(ierr); 475457d4a33Sstefano_zampini ierr = ISCreateGeneral(comm,nPl,l2g_indices_p,PETSC_OWN_POINTER,&IS_l2g_p);CHKERRQ(ierr); 476af140850Sstefano_zampini 477457d4a33Sstefano_zampini /* local to global scatter for interface pressure */ 478457d4a33Sstefano_zampini ierr = VecScatterCreate(fetidpmat_ctx->vP,NULL,fetidp_global,IS_l2g_p,&fetidpmat_ctx->l2g_p);CHKERRQ(ierr); 479457d4a33Sstefano_zampini ierr = ISDestroy(&IS_l2g_p);CHKERRQ(ierr); 480457d4a33Sstefano_zampini 481af140850Sstefano_zampini /* shift local to global indices for multipliers */ 482457d4a33Sstefano_zampini for (i=0;i<n_local_lambda;i++) { 483457d4a33Sstefano_zampini PetscInt owner,ps; 484457d4a33Sstefano_zampini 485457d4a33Sstefano_zampini ierr = PetscLayoutFindOwner(llay,l2g_indices[i],&owner);CHKERRQ(ierr); 486457d4a33Sstefano_zampini ps = pranges[owner+1]-pranges[owner]; 487457d4a33Sstefano_zampini l2g_indices[i] = l2g_indices[i]-lranges[owner]+aranges[owner]+ps; 488457d4a33Sstefano_zampini } 489457d4a33Sstefano_zampini 4909cc7774eSstefano_zampini /* scatter from alldofs to interface pressures global fetidp vector */ 4919cc7774eSstefano_zampini ierr = PetscLayoutGetRange(alay,&rst,NULL);CHKERRQ(ierr); 4929cc7774eSstefano_zampini ierr = ISCreateStride(comm,nPgl,rst,1,&ais);CHKERRQ(ierr); 493af140850Sstefano_zampini ierr = VecScatterCreate(pcis->vec1_global,pP,fetidp_global,ais,&fetidpmat_ctx->g2g_p);CHKERRQ(ierr); 4949cc7774eSstefano_zampini ierr = ISDestroy(&ais);CHKERRQ(ierr); 495457d4a33Sstefano_zampini } 496457d4a33Sstefano_zampini ierr = PetscLayoutDestroy(&llay);CHKERRQ(ierr); 497457d4a33Sstefano_zampini ierr = PetscLayoutDestroy(&play);CHKERRQ(ierr); 498457d4a33Sstefano_zampini ierr = ISCreateGeneral(comm,n_local_lambda,l2g_indices,PETSC_OWN_POINTER,&IS_l2g_lambda);CHKERRQ(ierr); 4999cc7774eSstefano_zampini /* scatter from local to global multipliers */ 500457d4a33Sstefano_zampini ierr = VecScatterCreate(fetidpmat_ctx->lambda_local,NULL,fetidp_global,IS_l2g_lambda,&fetidpmat_ctx->l2g_lambda);CHKERRQ(ierr); 501457d4a33Sstefano_zampini ierr = ISDestroy(&IS_l2g_lambda);CHKERRQ(ierr); 502457d4a33Sstefano_zampini ierr = ISLocalToGlobalMappingDestroy(&l2gmap_p);CHKERRQ(ierr); 503457d4a33Sstefano_zampini 504674ae819SStefano Zampini /* Create some vectors needed by fetidp */ 505674ae819SStefano Zampini ierr = VecDuplicate(pcis->vec1_B,&fetidpmat_ctx->temp_solution_B);CHKERRQ(ierr); 506674ae819SStefano Zampini ierr = VecDuplicate(pcis->vec1_D,&fetidpmat_ctx->temp_solution_D);CHKERRQ(ierr); 507674ae819SStefano Zampini 508674ae819SStefano Zampini test_fetidp = PETSC_FALSE; 509c5929fdfSBarry Smith ierr = PetscOptionsGetBool(NULL,NULL,"-fetidp_check",&test_fetidp,NULL);CHKERRQ(ierr); 510674ae819SStefano Zampini 511674ae819SStefano Zampini if (test_fetidp && !pcbddc->use_deluxe_scaling) { 5129eec4de8Sstefano_zampini Vec test_vec,test_vec_p = NULL; 5135b08dc53SStefano Zampini PetscReal real_value; 5145b08dc53SStefano Zampini 515674ae819SStefano Zampini ierr = PetscViewerASCIIGetStdout(comm,&viewer);CHKERRQ(ierr); 5161575c14dSBarry Smith ierr = PetscViewerASCIIPushSynchronized(viewer);CHKERRQ(ierr); 517457d4a33Sstefano_zampini ierr = PetscViewerASCIIPrintf(viewer,"----------FETI_DP SIZES--------------\n");CHKERRQ(ierr); 518457d4a33Sstefano_zampini ierr = PetscViewerASCIISynchronizedPrintf(viewer,"[%04d]: Sizes %D %D %D %D %D\n",rank,fetidpmat_ctx->n,fetidpmat_ctx->N,nPgl,nPg,nPl);CHKERRQ(ierr); 519457d4a33Sstefano_zampini ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 520674ae819SStefano Zampini ierr = PetscViewerASCIIPrintf(viewer,"----------FETI_DP TESTS--------------\n");CHKERRQ(ierr); 521674ae819SStefano Zampini ierr = PetscViewerASCIIPrintf(viewer,"All tests should return zero!\n");CHKERRQ(ierr); 522674ae819SStefano Zampini ierr = PetscViewerASCIIPrintf(viewer,"FETIDP MAT context in the ");CHKERRQ(ierr); 523674ae819SStefano Zampini if (fully_redundant) { 524674ae819SStefano Zampini ierr = PetscViewerASCIIPrintf(viewer,"fully redundant case for lagrange multipliers.\n");CHKERRQ(ierr); 525674ae819SStefano Zampini } else { 526674ae819SStefano Zampini ierr = PetscViewerASCIIPrintf(viewer,"Non-fully redundant case for lagrange multiplier.\n");CHKERRQ(ierr); 527674ae819SStefano Zampini } 528674ae819SStefano Zampini ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 529674ae819SStefano Zampini 530674ae819SStefano Zampini /******************************************************************/ 5319eec4de8Sstefano_zampini /* TEST A/B: Test numbering of global fetidp dofs */ 532674ae819SStefano Zampini /******************************************************************/ 533674ae819SStefano Zampini 534674ae819SStefano Zampini ierr = VecDuplicate(fetidpmat_ctx->lambda_local,&test_vec);CHKERRQ(ierr); 535457d4a33Sstefano_zampini ierr = VecSet(fetidp_global,1.0);CHKERRQ(ierr); 5369eec4de8Sstefano_zampini ierr = VecSet(test_vec,1.);CHKERRQ(ierr); 537457d4a33Sstefano_zampini ierr = VecScatterBegin(fetidpmat_ctx->l2g_lambda,fetidp_global,fetidpmat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 538457d4a33Sstefano_zampini ierr = VecScatterEnd(fetidpmat_ctx->l2g_lambda,fetidp_global,fetidpmat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 5399eec4de8Sstefano_zampini if (fetidpmat_ctx->l2g_p) { 5409eec4de8Sstefano_zampini ierr = VecDuplicate(fetidpmat_ctx->vP,&test_vec_p);CHKERRQ(ierr); 5419eec4de8Sstefano_zampini ierr = VecSet(test_vec_p,1.);CHKERRQ(ierr); 5429eec4de8Sstefano_zampini ierr = VecScatterBegin(fetidpmat_ctx->l2g_p,fetidp_global,fetidpmat_ctx->vP,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 5439eec4de8Sstefano_zampini ierr = VecScatterEnd(fetidpmat_ctx->l2g_p,fetidp_global,fetidpmat_ctx->vP,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 5449eec4de8Sstefano_zampini } 545674ae819SStefano Zampini scalar_value = -1.0; 546674ae819SStefano Zampini ierr = VecAXPY(test_vec,scalar_value,fetidpmat_ctx->lambda_local);CHKERRQ(ierr); 5475b08dc53SStefano Zampini ierr = VecNorm(test_vec,NORM_INFINITY,&real_value);CHKERRQ(ierr); 548674ae819SStefano Zampini ierr = VecDestroy(&test_vec);CHKERRQ(ierr); 5495b08dc53SStefano Zampini ierr = PetscViewerASCIISynchronizedPrintf(viewer,"A[%04d]: CHECK glob to loc: % 1.14e\n",rank,real_value);CHKERRQ(ierr); 5509eec4de8Sstefano_zampini if (fetidpmat_ctx->l2g_p) { 5519eec4de8Sstefano_zampini ierr = VecAXPY(test_vec_p,scalar_value,fetidpmat_ctx->vP);CHKERRQ(ierr); 5529eec4de8Sstefano_zampini ierr = VecNorm(test_vec_p,NORM_INFINITY,&real_value);CHKERRQ(ierr); 5539eec4de8Sstefano_zampini ierr = PetscViewerASCIISynchronizedPrintf(viewer,"A[%04d]: CHECK glob to loc (p): % 1.14e\n",rank,real_value);CHKERRQ(ierr); 5549eec4de8Sstefano_zampini } 555674ae819SStefano Zampini ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 556674ae819SStefano Zampini if (fully_redundant) { 557457d4a33Sstefano_zampini ierr = VecSet(fetidp_global,0.0);CHKERRQ(ierr); 558674ae819SStefano Zampini ierr = VecSet(fetidpmat_ctx->lambda_local,0.5);CHKERRQ(ierr); 559457d4a33Sstefano_zampini ierr = VecScatterBegin(fetidpmat_ctx->l2g_lambda,fetidpmat_ctx->lambda_local,fetidp_global,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 560457d4a33Sstefano_zampini ierr = VecScatterEnd(fetidpmat_ctx->l2g_lambda,fetidpmat_ctx->lambda_local,fetidp_global,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 561457d4a33Sstefano_zampini ierr = VecSum(fetidp_global,&scalar_value);CHKERRQ(ierr); 5625b08dc53SStefano Zampini ierr = PetscViewerASCIISynchronizedPrintf(viewer,"B[%04d]: CHECK loc to glob: % 1.14e\n",rank,PetscRealPart(scalar_value)-fetidpmat_ctx->n_lambda);CHKERRQ(ierr); 563674ae819SStefano Zampini ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 564674ae819SStefano Zampini } 5659eec4de8Sstefano_zampini if (fetidpmat_ctx->l2g_p) { 566af140850Sstefano_zampini ierr = VecSet(pcis->vec1_N,1.0);CHKERRQ(ierr); 567af140850Sstefano_zampini ierr = VecSet(pcis->vec1_global,0.0);CHKERRQ(ierr); 568af140850Sstefano_zampini ierr = VecScatterBegin(matis->rctx,pcis->vec1_N,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 569af140850Sstefano_zampini ierr = VecScatterEnd(matis->rctx,pcis->vec1_N,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 570af140850Sstefano_zampini 5719eec4de8Sstefano_zampini ierr = VecSet(fetidp_global,0.0);CHKERRQ(ierr); 572af140850Sstefano_zampini ierr = VecSet(fetidpmat_ctx->vP,-1.0);CHKERRQ(ierr); 573af140850Sstefano_zampini ierr = VecScatterBegin(fetidpmat_ctx->l2g_p,fetidpmat_ctx->vP,fetidp_global,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 574af140850Sstefano_zampini ierr = VecScatterEnd(fetidpmat_ctx->l2g_p,fetidpmat_ctx->vP,fetidp_global,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 575af140850Sstefano_zampini ierr = VecScatterBegin(fetidpmat_ctx->g2g_p,fetidp_global,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 576af140850Sstefano_zampini ierr = VecScatterEnd(fetidpmat_ctx->g2g_p,fetidp_global,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 577af140850Sstefano_zampini ierr = VecScatterBegin(fetidpmat_ctx->g2g_p,pcis->vec1_global,fetidp_global,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 578af140850Sstefano_zampini ierr = VecScatterEnd(fetidpmat_ctx->g2g_p,pcis->vec1_global,fetidp_global,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 5799eec4de8Sstefano_zampini ierr = VecSum(fetidp_global,&scalar_value);CHKERRQ(ierr); 580af140850Sstefano_zampini ierr = PetscViewerASCIISynchronizedPrintf(viewer,"B[%04d]: CHECK loc to glob (p): % 1.14e\n",rank,PetscRealPart(scalar_value));CHKERRQ(ierr); 5819eec4de8Sstefano_zampini } 582674ae819SStefano Zampini /******************************************************************/ 583674ae819SStefano Zampini /* TEST C: It should holds B_delta*w=0, w\in\widehat{W} */ 584674ae819SStefano Zampini /* This is the meaning of the B matrix */ 585674ae819SStefano Zampini /******************************************************************/ 586674ae819SStefano Zampini 587674ae819SStefano Zampini ierr = VecSetRandom(pcis->vec1_N,NULL);CHKERRQ(ierr); 588674ae819SStefano Zampini ierr = VecSet(pcis->vec1_global,0.0);CHKERRQ(ierr); 589e176bc59SStefano Zampini ierr = VecScatterBegin(matis->rctx,pcis->vec1_N,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 590e176bc59SStefano Zampini ierr = VecScatterEnd(matis->rctx,pcis->vec1_N,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 591e176bc59SStefano Zampini ierr = VecScatterBegin(matis->rctx,pcis->vec1_global,pcis->vec1_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 592e176bc59SStefano Zampini ierr = VecScatterEnd(matis->rctx,pcis->vec1_global,pcis->vec1_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 593674ae819SStefano Zampini ierr = VecScatterBegin(pcis->N_to_B,pcis->vec1_N,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 594674ae819SStefano Zampini ierr = VecScatterEnd(pcis->N_to_B,pcis->vec1_N,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 595674ae819SStefano Zampini /* Action of B_delta */ 596674ae819SStefano Zampini ierr = MatMult(fetidpmat_ctx->B_delta,pcis->vec1_B,fetidpmat_ctx->lambda_local);CHKERRQ(ierr); 597457d4a33Sstefano_zampini ierr = VecSet(fetidp_global,0.0);CHKERRQ(ierr); 598457d4a33Sstefano_zampini ierr = VecScatterBegin(fetidpmat_ctx->l2g_lambda,fetidpmat_ctx->lambda_local,fetidp_global,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 599457d4a33Sstefano_zampini ierr = VecScatterEnd(fetidpmat_ctx->l2g_lambda,fetidpmat_ctx->lambda_local,fetidp_global,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 600457d4a33Sstefano_zampini ierr = VecNorm(fetidp_global,NORM_INFINITY,&real_value);CHKERRQ(ierr); 6015b08dc53SStefano Zampini ierr = PetscViewerASCIIPrintf(viewer,"C[coll]: CHECK infty norm of B_delta*w (w continuous): % 1.14e\n",real_value);CHKERRQ(ierr); 602674ae819SStefano Zampini ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 603674ae819SStefano Zampini 604674ae819SStefano Zampini /******************************************************************/ 605674ae819SStefano Zampini /* TEST D: It should holds E_Dw = w - P_Dw w\in\widetilde{W} */ 606674ae819SStefano Zampini /* E_D = R_D^TR */ 607674ae819SStefano Zampini /* P_D = B_{D,delta}^T B_{delta} */ 608674ae819SStefano Zampini /* eq.44 Mandel Tezaur and Dohrmann 2005 */ 609674ae819SStefano Zampini /******************************************************************/ 610674ae819SStefano Zampini 611674ae819SStefano Zampini /* compute a random vector in \widetilde{W} */ 612674ae819SStefano Zampini ierr = VecSetRandom(pcis->vec1_N,NULL);CHKERRQ(ierr); 613674ae819SStefano Zampini scalar_value = 0.0; /* set zero at vertices */ 614674ae819SStefano Zampini ierr = VecGetArray(pcis->vec1_N,&array);CHKERRQ(ierr); 615674ae819SStefano Zampini for (i=0;i<n_vertices;i++) { array[vertex_indices[i]]=scalar_value; } 616674ae819SStefano Zampini ierr = VecRestoreArray(pcis->vec1_N,&array);CHKERRQ(ierr); 617674ae819SStefano Zampini /* store w for final comparison */ 618674ae819SStefano Zampini ierr = VecDuplicate(pcis->vec1_B,&test_vec);CHKERRQ(ierr); 619674ae819SStefano Zampini ierr = VecScatterBegin(pcis->N_to_B,pcis->vec1_N,test_vec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 620674ae819SStefano Zampini ierr = VecScatterEnd(pcis->N_to_B,pcis->vec1_N,test_vec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 621674ae819SStefano Zampini 622674ae819SStefano Zampini /* Jump operator P_D : results stored in pcis->vec1_B */ 623674ae819SStefano Zampini 624674ae819SStefano Zampini ierr = VecScatterBegin(pcis->N_to_B,pcis->vec1_N,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 625674ae819SStefano Zampini ierr = VecScatterEnd(pcis->N_to_B,pcis->vec1_N,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 626674ae819SStefano Zampini /* Action of B_delta */ 627674ae819SStefano Zampini ierr = MatMult(fetidpmat_ctx->B_delta,pcis->vec1_B,fetidpmat_ctx->lambda_local);CHKERRQ(ierr); 628457d4a33Sstefano_zampini ierr = VecSet(fetidp_global,0.0);CHKERRQ(ierr); 629457d4a33Sstefano_zampini ierr = VecScatterBegin(fetidpmat_ctx->l2g_lambda,fetidpmat_ctx->lambda_local,fetidp_global,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 630457d4a33Sstefano_zampini ierr = VecScatterEnd(fetidpmat_ctx->l2g_lambda,fetidpmat_ctx->lambda_local,fetidp_global,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 631674ae819SStefano Zampini /* Action of B_Ddelta^T */ 632457d4a33Sstefano_zampini ierr = VecScatterBegin(fetidpmat_ctx->l2g_lambda,fetidp_global,fetidpmat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 633457d4a33Sstefano_zampini ierr = VecScatterEnd(fetidpmat_ctx->l2g_lambda,fetidp_global,fetidpmat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 634674ae819SStefano Zampini ierr = MatMultTranspose(fetidpmat_ctx->B_Ddelta,fetidpmat_ctx->lambda_local,pcis->vec1_B);CHKERRQ(ierr); 635674ae819SStefano Zampini 636674ae819SStefano Zampini /* Average operator E_D : results stored in pcis->vec2_B */ 637674ae819SStefano Zampini ierr = VecScatterBegin(pcis->N_to_B,pcis->vec1_N,pcis->vec2_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 638674ae819SStefano Zampini ierr = VecScatterEnd(pcis->N_to_B,pcis->vec1_N,pcis->vec2_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 639674ae819SStefano Zampini ierr = PCBDDCScalingExtension(fetidpmat_ctx->pc,pcis->vec2_B,pcis->vec1_global);CHKERRQ(ierr); 640674ae819SStefano Zampini ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_global,pcis->vec2_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 641674ae819SStefano Zampini ierr = VecScatterEnd(pcis->global_to_B,pcis->vec1_global,pcis->vec2_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 642674ae819SStefano Zampini 643674ae819SStefano Zampini /* test E_D=I-P_D */ 644674ae819SStefano Zampini scalar_value = 1.0; 645674ae819SStefano Zampini ierr = VecAXPY(pcis->vec1_B,scalar_value,pcis->vec2_B);CHKERRQ(ierr); 646674ae819SStefano Zampini scalar_value = -1.0; 647674ae819SStefano Zampini ierr = VecAXPY(pcis->vec1_B,scalar_value,test_vec);CHKERRQ(ierr); 6485b08dc53SStefano Zampini ierr = VecNorm(pcis->vec1_B,NORM_INFINITY,&real_value);CHKERRQ(ierr); 649674ae819SStefano Zampini ierr = VecDestroy(&test_vec);CHKERRQ(ierr); 6505b08dc53SStefano Zampini ierr = PetscViewerASCIISynchronizedPrintf(viewer,"D[%04d] CHECK infty norm of E_D + P_D - I: % 1.14e\n",rank,real_value);CHKERRQ(ierr); 651674ae819SStefano Zampini ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 652674ae819SStefano Zampini 653674ae819SStefano Zampini /******************************************************************/ 654674ae819SStefano Zampini /* TEST E: It should holds R_D^TP_Dw=0 w\in\widetilde{W} */ 655674ae819SStefano Zampini /* eq.48 Mandel Tezaur and Dohrmann 2005 */ 656674ae819SStefano Zampini /******************************************************************/ 657674ae819SStefano Zampini 658674ae819SStefano Zampini ierr = VecSetRandom(pcis->vec1_N,NULL);CHKERRQ(ierr); 659674ae819SStefano Zampini ierr = VecGetArray(pcis->vec1_N,&array);CHKERRQ(ierr); 660674ae819SStefano Zampini scalar_value = 0.0; /* set zero at vertices */ 661674ae819SStefano Zampini for (i=0;i<n_vertices;i++) { array[vertex_indices[i]]=scalar_value; } 662674ae819SStefano Zampini ierr = VecRestoreArray(pcis->vec1_N,&array);CHKERRQ(ierr); 663674ae819SStefano Zampini 664674ae819SStefano Zampini /* Jump operator P_D : results stored in pcis->vec1_B */ 665674ae819SStefano Zampini 666674ae819SStefano Zampini ierr = VecScatterBegin(pcis->N_to_B,pcis->vec1_N,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 667674ae819SStefano Zampini ierr = VecScatterEnd(pcis->N_to_B,pcis->vec1_N,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 668674ae819SStefano Zampini /* Action of B_delta */ 669674ae819SStefano Zampini ierr = MatMult(fetidpmat_ctx->B_delta,pcis->vec1_B,fetidpmat_ctx->lambda_local);CHKERRQ(ierr); 670457d4a33Sstefano_zampini ierr = VecSet(fetidp_global,0.0);CHKERRQ(ierr); 671457d4a33Sstefano_zampini ierr = VecScatterBegin(fetidpmat_ctx->l2g_lambda,fetidpmat_ctx->lambda_local,fetidp_global,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 672457d4a33Sstefano_zampini ierr = VecScatterEnd(fetidpmat_ctx->l2g_lambda,fetidpmat_ctx->lambda_local,fetidp_global,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 673674ae819SStefano Zampini /* Action of B_Ddelta^T */ 674457d4a33Sstefano_zampini ierr = VecScatterBegin(fetidpmat_ctx->l2g_lambda,fetidp_global,fetidpmat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 675457d4a33Sstefano_zampini ierr = VecScatterEnd(fetidpmat_ctx->l2g_lambda,fetidp_global,fetidpmat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 676674ae819SStefano Zampini ierr = MatMultTranspose(fetidpmat_ctx->B_Ddelta,fetidpmat_ctx->lambda_local,pcis->vec1_B);CHKERRQ(ierr); 677674ae819SStefano Zampini /* scaling */ 678674ae819SStefano Zampini ierr = PCBDDCScalingExtension(fetidpmat_ctx->pc,pcis->vec1_B,pcis->vec1_global);CHKERRQ(ierr); 6795b08dc53SStefano Zampini ierr = VecNorm(pcis->vec1_global,NORM_INFINITY,&real_value);CHKERRQ(ierr); 6805b08dc53SStefano Zampini ierr = PetscViewerASCIIPrintf(viewer,"E[coll]: CHECK infty norm of R^T_D P_D: % 1.14e\n",real_value);CHKERRQ(ierr); 681674ae819SStefano Zampini ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 682674ae819SStefano Zampini 683674ae819SStefano Zampini if (!fully_redundant) { 684674ae819SStefano Zampini /******************************************************************/ 685674ae819SStefano Zampini /* TEST F: It should holds B_{delta}B^T_{D,delta}=I */ 686674ae819SStefano Zampini /* Corollary thm 14 Mandel Tezaur and Dohrmann 2005 */ 687674ae819SStefano Zampini /******************************************************************/ 688457d4a33Sstefano_zampini ierr = VecDuplicate(fetidp_global,&test_vec);CHKERRQ(ierr); 689457d4a33Sstefano_zampini ierr = VecSetRandom(fetidp_global,NULL);CHKERRQ(ierr); 690457d4a33Sstefano_zampini if (fetidpmat_ctx->l2g_p) { 691457d4a33Sstefano_zampini ierr = VecSet(fetidpmat_ctx->vP,0.);CHKERRQ(ierr); 692457d4a33Sstefano_zampini ierr = VecScatterBegin(fetidpmat_ctx->l2g_p,fetidpmat_ctx->vP,fetidp_global,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 693457d4a33Sstefano_zampini ierr = VecScatterEnd(fetidpmat_ctx->l2g_p,fetidpmat_ctx->vP,fetidp_global,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 694457d4a33Sstefano_zampini } 695674ae819SStefano Zampini /* Action of B_Ddelta^T */ 696457d4a33Sstefano_zampini ierr = VecScatterBegin(fetidpmat_ctx->l2g_lambda,fetidp_global,fetidpmat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 697457d4a33Sstefano_zampini ierr = VecScatterEnd(fetidpmat_ctx->l2g_lambda,fetidp_global,fetidpmat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 698674ae819SStefano Zampini ierr = MatMultTranspose(fetidpmat_ctx->B_Ddelta,fetidpmat_ctx->lambda_local,pcis->vec1_B);CHKERRQ(ierr); 699674ae819SStefano Zampini /* Action of B_delta */ 700674ae819SStefano Zampini ierr = MatMult(fetidpmat_ctx->B_delta,pcis->vec1_B,fetidpmat_ctx->lambda_local);CHKERRQ(ierr); 701674ae819SStefano Zampini ierr = VecSet(test_vec,0.0);CHKERRQ(ierr); 702674ae819SStefano Zampini ierr = VecScatterBegin(fetidpmat_ctx->l2g_lambda,fetidpmat_ctx->lambda_local,test_vec,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 703674ae819SStefano Zampini ierr = VecScatterEnd(fetidpmat_ctx->l2g_lambda,fetidpmat_ctx->lambda_local,test_vec,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 704674ae819SStefano Zampini scalar_value = -1.0; 705457d4a33Sstefano_zampini ierr = VecAXPY(fetidp_global,scalar_value,test_vec);CHKERRQ(ierr); 706457d4a33Sstefano_zampini ierr = VecNorm(fetidp_global,NORM_INFINITY,&real_value);CHKERRQ(ierr); 7075b08dc53SStefano Zampini ierr = PetscViewerASCIIPrintf(viewer,"E[coll]: CHECK infty norm of P^T_D - I: % 1.14e\n",real_value);CHKERRQ(ierr); 708674ae819SStefano Zampini ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 709674ae819SStefano Zampini ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 710674ae819SStefano Zampini ierr = VecDestroy(&test_vec);CHKERRQ(ierr); 711674ae819SStefano Zampini } 7129eec4de8Sstefano_zampini ierr = VecDestroy(&test_vec_p);CHKERRQ(ierr); 713674ae819SStefano Zampini } 714674ae819SStefano Zampini /* final cleanup */ 715457d4a33Sstefano_zampini ierr = VecDestroy(&fetidp_global);CHKERRQ(ierr); 716674ae819SStefano Zampini PetscFunctionReturn(0); 717674ae819SStefano Zampini } 718674ae819SStefano Zampini 719674ae819SStefano Zampini #undef __FUNCT__ 720674ae819SStefano Zampini #define __FUNCT__ "PCBDDCSetupFETIDPPCContext" 721674ae819SStefano Zampini PetscErrorCode PCBDDCSetupFETIDPPCContext(Mat fetimat, FETIDPPC_ctx fetidppc_ctx) 722674ae819SStefano Zampini { 723674ae819SStefano Zampini FETIDPMat_ctx mat_ctx; 724ed6c3d69SStefano Zampini PC_IS *pcis; 725f28b6018SStefano Zampini PetscBool lumped = PETSC_FALSE; 726674ae819SStefano Zampini PetscErrorCode ierr; 727674ae819SStefano Zampini 728674ae819SStefano Zampini PetscFunctionBegin; 729674ae819SStefano Zampini ierr = MatShellGetContext(fetimat,(void**)&mat_ctx);CHKERRQ(ierr); 730674ae819SStefano Zampini /* get references from objects created when setting up feti mat context */ 731674ae819SStefano Zampini ierr = PetscObjectReference((PetscObject)mat_ctx->lambda_local);CHKERRQ(ierr); 732674ae819SStefano Zampini fetidppc_ctx->lambda_local = mat_ctx->lambda_local; 733674ae819SStefano Zampini ierr = PetscObjectReference((PetscObject)mat_ctx->B_Ddelta);CHKERRQ(ierr); 734674ae819SStefano Zampini fetidppc_ctx->B_Ddelta = mat_ctx->B_Ddelta; 735674ae819SStefano Zampini ierr = PetscObjectReference((PetscObject)mat_ctx->l2g_lambda);CHKERRQ(ierr); 736674ae819SStefano Zampini fetidppc_ctx->l2g_lambda = mat_ctx->l2g_lambda; 737f28b6018SStefano Zampini /* create preconditioner */ 738ed6c3d69SStefano Zampini pcis = (PC_IS*)fetidppc_ctx->pc->data; 739f28b6018SStefano Zampini /* Dirichlet preconditioner */ 740f28b6018SStefano Zampini ierr = PetscOptionsGetBool(NULL,NULL,"-fetidp_lumped",&lumped,NULL);CHKERRQ(ierr); 741f28b6018SStefano Zampini if (!lumped) { 742ed6c3d69SStefano Zampini ierr = MatCreateSchurComplement(pcis->A_II,pcis->A_II,pcis->A_IB,pcis->A_BI,pcis->A_BB,&fetidppc_ctx->S_j);CHKERRQ(ierr); 743ed6c3d69SStefano Zampini ierr = MatSchurComplementSetKSP(fetidppc_ctx->S_j,pcis->ksp_D);CHKERRQ(ierr); 744f28b6018SStefano Zampini } else { 745f28b6018SStefano Zampini ierr = PetscObjectReference((PetscObject)pcis->A_BB);CHKERRQ(ierr); 746f28b6018SStefano Zampini fetidppc_ctx->S_j = pcis->A_BB; 747f28b6018SStefano Zampini } 748af140850Sstefano_zampini /* saddle-point */ 749af140850Sstefano_zampini if (mat_ctx->xPg) { 750af140850Sstefano_zampini Mat PAmat = NULL,PPmat = NULL; 751af140850Sstefano_zampini 752af140850Sstefano_zampini ierr = PetscObjectReference((PetscObject)mat_ctx->xPg);CHKERRQ(ierr); 753af140850Sstefano_zampini fetidppc_ctx->xPg = mat_ctx->xPg; 754af140850Sstefano_zampini ierr = PetscObjectReference((PetscObject)mat_ctx->yPg);CHKERRQ(ierr); 755af140850Sstefano_zampini fetidppc_ctx->yPg = mat_ctx->yPg; 756af140850Sstefano_zampini ierr = PetscObjectQuery((PetscObject)fetidppc_ctx->pc,"__KSPFETIDP_PPmat",(PetscObject*)&PPmat);CHKERRQ(ierr); 757af140850Sstefano_zampini ierr = PetscObjectQuery((PetscObject)fetidppc_ctx->pc,"__KSPFETIDP_PAmat",(PetscObject*)&PAmat);CHKERRQ(ierr); 758af140850Sstefano_zampini ierr = KSPCreate(PetscObjectComm((PetscObject)fetidppc_ctx->pc),&fetidppc_ctx->kP);CHKERRQ(ierr); 759af140850Sstefano_zampini ierr = KSPSetOperators(fetidppc_ctx->kP,PAmat,PPmat);CHKERRQ(ierr); 760af140850Sstefano_zampini ierr = KSPSetOptionsPrefix(fetidppc_ctx->kP,"fetidp_pmass_");CHKERRQ(ierr); 761af140850Sstefano_zampini ierr = KSPSetFromOptions(fetidppc_ctx->kP);CHKERRQ(ierr); 762af140850Sstefano_zampini } 763674ae819SStefano Zampini PetscFunctionReturn(0); 764674ae819SStefano Zampini } 765674ae819SStefano Zampini 766674ae819SStefano Zampini #undef __FUNCT__ 767674ae819SStefano Zampini #define __FUNCT__ "FETIDPMatMult" 768674ae819SStefano Zampini PetscErrorCode FETIDPMatMult(Mat fetimat, Vec x, Vec y) 769674ae819SStefano Zampini { 770674ae819SStefano Zampini FETIDPMat_ctx mat_ctx; 771617d11aeSStefano Zampini PC_BDDC *pcbddc; 772674ae819SStefano Zampini PC_IS *pcis; 773674ae819SStefano Zampini PetscErrorCode ierr; 774674ae819SStefano Zampini 775674ae819SStefano Zampini PetscFunctionBegin; 776674ae819SStefano Zampini ierr = MatShellGetContext(fetimat,(void**)&mat_ctx);CHKERRQ(ierr); 777674ae819SStefano Zampini pcis = (PC_IS*)mat_ctx->pc->data; 778617d11aeSStefano Zampini pcbddc = (PC_BDDC*)mat_ctx->pc->data; 779674ae819SStefano Zampini /* Application of B_delta^T */ 780af140850Sstefano_zampini ierr = VecSet(pcis->vec1_B,0.);CHKERRQ(ierr); 781674ae819SStefano Zampini ierr = VecScatterBegin(mat_ctx->l2g_lambda,x,mat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 782674ae819SStefano Zampini ierr = VecScatterEnd(mat_ctx->l2g_lambda,x,mat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 783674ae819SStefano Zampini ierr = MatMultTranspose(mat_ctx->B_delta,mat_ctx->lambda_local,pcis->vec1_B);CHKERRQ(ierr); 784af140850Sstefano_zampini 785af140850Sstefano_zampini /* Add contribution from saddle point */ 786af140850Sstefano_zampini if (mat_ctx->l2g_p) { 787af140850Sstefano_zampini ierr = VecScatterBegin(mat_ctx->l2g_p,x,mat_ctx->vP,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 788af140850Sstefano_zampini ierr = VecScatterEnd(mat_ctx->l2g_p,x,mat_ctx->vP,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 789af140850Sstefano_zampini if (pcbddc->switch_static) { 790af140850Sstefano_zampini ierr = MatMult(mat_ctx->Bt_BI,mat_ctx->vP,pcis->vec1_D);CHKERRQ(ierr); 791af140850Sstefano_zampini } 792af140850Sstefano_zampini ierr = MatMultAdd(mat_ctx->Bt_BB,mat_ctx->vP,pcis->vec1_B,pcis->vec1_B);CHKERRQ(ierr); 793af140850Sstefano_zampini } else { 794af140850Sstefano_zampini if (pcbddc->switch_static) { 795674ae819SStefano Zampini ierr = VecSet(pcis->vec1_D,0.0);CHKERRQ(ierr); 796af140850Sstefano_zampini } 797af140850Sstefano_zampini } 798af140850Sstefano_zampini /* Application of \widetilde{S}^-1 */ 799617d11aeSStefano Zampini ierr = PetscMemzero(pcbddc->benign_p0,pcbddc->benign_n*sizeof(PetscScalar));CHKERRQ(ierr); 800dc359a40SStefano Zampini ierr = PCBDDCApplyInterfacePreconditioner(mat_ctx->pc,PETSC_FALSE);CHKERRQ(ierr); 801c7ffc8ceSStefano Zampini ierr = PetscMemzero(pcbddc->benign_p0,pcbddc->benign_n*sizeof(PetscScalar));CHKERRQ(ierr); 802af140850Sstefano_zampini ierr = VecSet(y,0.0);CHKERRQ(ierr); 803674ae819SStefano Zampini /* Application of B_delta */ 804674ae819SStefano Zampini ierr = MatMult(mat_ctx->B_delta,pcis->vec1_B,mat_ctx->lambda_local);CHKERRQ(ierr); 805af140850Sstefano_zampini /* Contribution from boundary pressures */ 806af140850Sstefano_zampini if (mat_ctx->C) { 807af140850Sstefano_zampini const PetscScalar *lx; 808af140850Sstefano_zampini PetscScalar *ly; 809af140850Sstefano_zampini 810af140850Sstefano_zampini /* pressures ordered first in x and y */ 811af140850Sstefano_zampini ierr = VecGetArrayRead(x,&lx);CHKERRQ(ierr); 812af140850Sstefano_zampini ierr = VecGetArray(y,&ly);CHKERRQ(ierr); 813af140850Sstefano_zampini ierr = VecPlaceArray(mat_ctx->xPg,lx);CHKERRQ(ierr); 814af140850Sstefano_zampini ierr = VecPlaceArray(mat_ctx->yPg,ly);CHKERRQ(ierr); 815af140850Sstefano_zampini ierr = MatMult(mat_ctx->C,mat_ctx->xPg,mat_ctx->yPg);CHKERRQ(ierr); 816af140850Sstefano_zampini ierr = VecResetArray(mat_ctx->xPg);CHKERRQ(ierr); 817af140850Sstefano_zampini ierr = VecResetArray(mat_ctx->yPg);CHKERRQ(ierr); 818af140850Sstefano_zampini ierr = VecRestoreArrayRead(x,&lx);CHKERRQ(ierr); 819af140850Sstefano_zampini ierr = VecRestoreArray(y,&ly);CHKERRQ(ierr); 820af140850Sstefano_zampini } 821af140850Sstefano_zampini /* Add contribution from saddle point */ 822af140850Sstefano_zampini if (mat_ctx->l2g_p) { 823af140850Sstefano_zampini ierr = MatMult(mat_ctx->B_BB,pcis->vec1_B,mat_ctx->vP);CHKERRQ(ierr); 824af140850Sstefano_zampini if (pcbddc->switch_static) { 825af140850Sstefano_zampini ierr = MatMultAdd(mat_ctx->B_BI,pcis->vec1_D,mat_ctx->vP,mat_ctx->vP);CHKERRQ(ierr); 826af140850Sstefano_zampini } 827af140850Sstefano_zampini ierr = VecScatterBegin(mat_ctx->l2g_p,mat_ctx->vP,y,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 828af140850Sstefano_zampini ierr = VecScatterEnd(mat_ctx->l2g_p,mat_ctx->vP,y,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 829af140850Sstefano_zampini } 830674ae819SStefano Zampini ierr = VecScatterBegin(mat_ctx->l2g_lambda,mat_ctx->lambda_local,y,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 831674ae819SStefano Zampini ierr = VecScatterEnd(mat_ctx->l2g_lambda,mat_ctx->lambda_local,y,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 832674ae819SStefano Zampini PetscFunctionReturn(0); 833674ae819SStefano Zampini } 834674ae819SStefano Zampini 835674ae819SStefano Zampini #undef __FUNCT__ 836edf7251bSStefano Zampini #define __FUNCT__ "FETIDPMatMultTranspose" 837edf7251bSStefano Zampini PetscErrorCode FETIDPMatMultTranspose(Mat fetimat, Vec x, Vec y) 838edf7251bSStefano Zampini { 839edf7251bSStefano Zampini FETIDPMat_ctx mat_ctx; 840edf7251bSStefano Zampini PC_IS *pcis; 841edf7251bSStefano Zampini PetscErrorCode ierr; 842edf7251bSStefano Zampini 843edf7251bSStefano Zampini PetscFunctionBegin; 844edf7251bSStefano Zampini ierr = MatShellGetContext(fetimat,(void**)&mat_ctx);CHKERRQ(ierr); 845edf7251bSStefano Zampini pcis = (PC_IS*)mat_ctx->pc->data; 846edf7251bSStefano Zampini /* Application of B_delta^T */ 847edf7251bSStefano Zampini ierr = VecScatterBegin(mat_ctx->l2g_lambda,x,mat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 848edf7251bSStefano Zampini ierr = VecScatterEnd(mat_ctx->l2g_lambda,x,mat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 849edf7251bSStefano Zampini ierr = MatMultTranspose(mat_ctx->B_delta,mat_ctx->lambda_local,pcis->vec1_B);CHKERRQ(ierr); 850edf7251bSStefano Zampini /* Application of \widetilde{S}^-1 */ 851edf7251bSStefano Zampini ierr = VecSet(pcis->vec1_D,0.0);CHKERRQ(ierr); 852edf7251bSStefano Zampini ierr = PCBDDCApplyInterfacePreconditioner(mat_ctx->pc,PETSC_TRUE);CHKERRQ(ierr); 853edf7251bSStefano Zampini /* Application of B_delta */ 854edf7251bSStefano Zampini ierr = MatMult(mat_ctx->B_delta,pcis->vec1_B,mat_ctx->lambda_local);CHKERRQ(ierr); 855edf7251bSStefano Zampini ierr = VecSet(y,0.0);CHKERRQ(ierr); 856edf7251bSStefano Zampini ierr = VecScatterBegin(mat_ctx->l2g_lambda,mat_ctx->lambda_local,y,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 857edf7251bSStefano Zampini ierr = VecScatterEnd(mat_ctx->l2g_lambda,mat_ctx->lambda_local,y,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 858edf7251bSStefano Zampini PetscFunctionReturn(0); 859edf7251bSStefano Zampini } 860edf7251bSStefano Zampini 861edf7251bSStefano Zampini #undef __FUNCT__ 862674ae819SStefano Zampini #define __FUNCT__ "FETIDPPCApply" 863674ae819SStefano Zampini PetscErrorCode FETIDPPCApply(PC fetipc, Vec x, Vec y) 864674ae819SStefano Zampini { 865674ae819SStefano Zampini FETIDPPC_ctx pc_ctx; 866674ae819SStefano Zampini PC_IS *pcis; 867674ae819SStefano Zampini PetscErrorCode ierr; 868674ae819SStefano Zampini 869674ae819SStefano Zampini PetscFunctionBegin; 870302440fdSBarry Smith ierr = PCShellGetContext(fetipc,(void**)&pc_ctx);CHKERRQ(ierr); 871674ae819SStefano Zampini pcis = (PC_IS*)pc_ctx->pc->data; 872674ae819SStefano Zampini /* Application of B_Ddelta^T */ 873674ae819SStefano Zampini ierr = VecScatterBegin(pc_ctx->l2g_lambda,x,pc_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 874674ae819SStefano Zampini ierr = VecScatterEnd(pc_ctx->l2g_lambda,x,pc_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 875674ae819SStefano Zampini ierr = VecSet(pcis->vec2_B,0.0);CHKERRQ(ierr); 876674ae819SStefano Zampini ierr = MatMultTranspose(pc_ctx->B_Ddelta,pc_ctx->lambda_local,pcis->vec2_B);CHKERRQ(ierr); 877ed6c3d69SStefano Zampini /* Application of local Schur complement */ 878ed6c3d69SStefano Zampini ierr = MatMult(pc_ctx->S_j,pcis->vec2_B,pcis->vec1_B);CHKERRQ(ierr); 879edf7251bSStefano Zampini /* Application of B_Ddelta */ 880edf7251bSStefano Zampini ierr = MatMult(pc_ctx->B_Ddelta,pcis->vec1_B,pc_ctx->lambda_local);CHKERRQ(ierr); 881edf7251bSStefano Zampini ierr = VecSet(y,0.0);CHKERRQ(ierr); 882edf7251bSStefano Zampini ierr = VecScatterBegin(pc_ctx->l2g_lambda,pc_ctx->lambda_local,y,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 883edf7251bSStefano Zampini ierr = VecScatterEnd(pc_ctx->l2g_lambda,pc_ctx->lambda_local,y,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 884edf7251bSStefano Zampini PetscFunctionReturn(0); 885edf7251bSStefano Zampini } 886edf7251bSStefano Zampini 887edf7251bSStefano Zampini #undef __FUNCT__ 888edf7251bSStefano Zampini #define __FUNCT__ "FETIDPPCApplyTranspose" 889edf7251bSStefano Zampini PetscErrorCode FETIDPPCApplyTranspose(PC fetipc, Vec x, Vec y) 890edf7251bSStefano Zampini { 891edf7251bSStefano Zampini FETIDPPC_ctx pc_ctx; 892edf7251bSStefano Zampini PC_IS *pcis; 893edf7251bSStefano Zampini PetscErrorCode ierr; 894edf7251bSStefano Zampini 895edf7251bSStefano Zampini PetscFunctionBegin; 896302440fdSBarry Smith ierr = PCShellGetContext(fetipc,(void**)&pc_ctx);CHKERRQ(ierr); 897edf7251bSStefano Zampini pcis = (PC_IS*)pc_ctx->pc->data; 898edf7251bSStefano Zampini /* Application of B_Ddelta^T */ 899edf7251bSStefano Zampini ierr = VecScatterBegin(pc_ctx->l2g_lambda,x,pc_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 900edf7251bSStefano Zampini ierr = VecScatterEnd(pc_ctx->l2g_lambda,x,pc_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 901edf7251bSStefano Zampini ierr = VecSet(pcis->vec2_B,0.0);CHKERRQ(ierr); 902edf7251bSStefano Zampini ierr = MatMultTranspose(pc_ctx->B_Ddelta,pc_ctx->lambda_local,pcis->vec2_B);CHKERRQ(ierr); 903ed6c3d69SStefano Zampini /* Application of local Schur complement */ 904ed6c3d69SStefano Zampini ierr = MatMultTranspose(pc_ctx->S_j,pcis->vec2_B,pcis->vec1_B);CHKERRQ(ierr); 905674ae819SStefano Zampini /* Application of B_Ddelta */ 906674ae819SStefano Zampini ierr = MatMult(pc_ctx->B_Ddelta,pcis->vec1_B,pc_ctx->lambda_local);CHKERRQ(ierr); 907674ae819SStefano Zampini ierr = VecSet(y,0.0);CHKERRQ(ierr); 908674ae819SStefano Zampini ierr = VecScatterBegin(pc_ctx->l2g_lambda,pc_ctx->lambda_local,y,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 909674ae819SStefano Zampini ierr = VecScatterEnd(pc_ctx->l2g_lambda,pc_ctx->lambda_local,y,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 910674ae819SStefano Zampini PetscFunctionReturn(0); 911674ae819SStefano Zampini } 912