1ab5c6b0cSJed Brown #include <../src/ksp/pc/impls/bddc/bddc.h> 2ab5c6b0cSJed Brown #include <../src/ksp/pc/impls/bddc/bddcprivate.h> 3674ae819SStefano Zampini 434a97f8cSStefano Zampini /* prototypes for deluxe functions */ 534a97f8cSStefano Zampini static PetscErrorCode PCBDDCScalingCreate_Deluxe(PC); 634a97f8cSStefano Zampini static PetscErrorCode PCBDDCScalingDestroy_Deluxe(PC); 734a97f8cSStefano Zampini static PetscErrorCode PCBDDCScalingSetUp_Deluxe(PC); 834a97f8cSStefano Zampini static PetscErrorCode PCBDDCScalingSetUp_Deluxe_Par(PC,PetscInt,PetscInt,PetscInt[],PetscInt[]); 934a97f8cSStefano Zampini static PetscErrorCode PCBDDCScalingSetUp_Deluxe_Seq(PC,PetscInt,PetscInt,PetscInt[],PetscInt[]); 1034a97f8cSStefano Zampini static PetscErrorCode PCBDDCScalingReset_Deluxe_Solvers(PCBDDCDeluxeScaling); 11674ae819SStefano Zampini 12674ae819SStefano Zampini #undef __FUNCT__ 13674ae819SStefano Zampini #define __FUNCT__ "PCBDDCScalingExtension_Basic" 14674ae819SStefano Zampini static PetscErrorCode PCBDDCScalingExtension_Basic(PC pc, Vec local_interface_vector, Vec global_vector) 15674ae819SStefano Zampini { 16674ae819SStefano Zampini PC_IS* pcis = (PC_IS*)pc->data; 17674ae819SStefano Zampini PC_BDDC* pcbddc = (PC_BDDC*)pc->data; 18674ae819SStefano Zampini PetscErrorCode ierr; 19674ae819SStefano Zampini 20674ae819SStefano Zampini PetscFunctionBegin; 21674ae819SStefano Zampini /* Apply partition of unity */ 22674ae819SStefano Zampini ierr = VecPointwiseMult(pcbddc->work_scaling,pcis->D,local_interface_vector);CHKERRQ(ierr); 23674ae819SStefano Zampini ierr = VecSet(global_vector,0.0);CHKERRQ(ierr); 24674ae819SStefano Zampini ierr = VecScatterBegin(pcis->global_to_B,pcbddc->work_scaling,global_vector,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 25674ae819SStefano Zampini ierr = VecScatterEnd(pcis->global_to_B,pcbddc->work_scaling,global_vector,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 26674ae819SStefano Zampini PetscFunctionReturn(0); 27674ae819SStefano Zampini } 28674ae819SStefano Zampini 29674ae819SStefano Zampini #undef __FUNCT__ 30674ae819SStefano Zampini #define __FUNCT__ "PCBDDCScalingExtension_Deluxe" 31674ae819SStefano Zampini static PetscErrorCode PCBDDCScalingExtension_Deluxe(PC pc, Vec x, Vec y) 32674ae819SStefano Zampini { 33674ae819SStefano Zampini PC_IS* pcis=(PC_IS*)pc->data; 34674ae819SStefano Zampini PC_BDDC* pcbddc=(PC_BDDC*)pc->data; 35674ae819SStefano Zampini PCBDDCDeluxeScaling deluxe_ctx = pcbddc->deluxe_ctx; 36b1b3d7a2SStefano Zampini PCBDDCSubSchurs sub_schurs = pcbddc->sub_schurs[1]; 37674ae819SStefano Zampini PetscInt i; 38674ae819SStefano Zampini PetscErrorCode ierr; 39674ae819SStefano Zampini 40674ae819SStefano Zampini /* TODO CHECK STUFF RELATED WITH FAKE WORK */ 41674ae819SStefano Zampini PetscFunctionBegin; 4234a97f8cSStefano Zampini ierr = VecSet(pcbddc->work_scaling,0.0);CHKERRQ(ierr); /* needed by the fake work below */ 4334a97f8cSStefano Zampini if (deluxe_ctx->n_simple) { 44674ae819SStefano Zampini /* scale deluxe vertices using diagonal scaling */ 4534a97f8cSStefano Zampini PetscScalar *array_x,*array_D,*array; 46674ae819SStefano Zampini ierr = VecGetArray(x,&array_x);CHKERRQ(ierr); 47674ae819SStefano Zampini ierr = VecGetArray(pcis->D,&array_D);CHKERRQ(ierr); 48674ae819SStefano Zampini ierr = VecGetArray(pcbddc->work_scaling,&array);CHKERRQ(ierr); 49674ae819SStefano Zampini for (i=0;i<deluxe_ctx->n_simple;i++) { 50674ae819SStefano Zampini array[deluxe_ctx->idx_simple_B[i]] = array_x[deluxe_ctx->idx_simple_B[i]]*array_D[deluxe_ctx->idx_simple_B[i]]; 51674ae819SStefano Zampini } 52674ae819SStefano Zampini ierr = VecRestoreArray(pcbddc->work_scaling,&array);CHKERRQ(ierr); 53674ae819SStefano Zampini ierr = VecRestoreArray(pcis->D,&array_D);CHKERRQ(ierr); 54674ae819SStefano Zampini ierr = VecRestoreArray(x,&array_x);CHKERRQ(ierr); 5534a97f8cSStefano Zampini } 5634a97f8cSStefano Zampini /* sequential part : all problems and Schur applications collapsed into a single matrix vector multiplication and ksp solution */ 57*41c3ba1bSStefano Zampini if (deluxe_ctx->seq_ksp) { 58674ae819SStefano Zampini ierr = VecScatterBegin(deluxe_ctx->seq_scctx,x,deluxe_ctx->seq_work1,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 59674ae819SStefano Zampini ierr = VecScatterEnd(deluxe_ctx->seq_scctx,x,deluxe_ctx->seq_work1,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 60*41c3ba1bSStefano Zampini ierr = MatMult(sub_schurs->S_Ej_all,deluxe_ctx->seq_work1,deluxe_ctx->seq_work2);CHKERRQ(ierr); 61674ae819SStefano Zampini ierr = KSPSolve(deluxe_ctx->seq_ksp,deluxe_ctx->seq_work2,deluxe_ctx->seq_work1);CHKERRQ(ierr); 62674ae819SStefano Zampini /* fake work due to final ADD VALUES and vertices scaling needed? TODO: check it */ 63674ae819SStefano Zampini ierr = VecScatterBegin(deluxe_ctx->seq_scctx,deluxe_ctx->seq_work1,pcbddc->work_scaling,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 64674ae819SStefano Zampini ierr = VecScatterEnd(deluxe_ctx->seq_scctx,deluxe_ctx->seq_work1,pcbddc->work_scaling,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 65674ae819SStefano Zampini } 66674ae819SStefano Zampini /* parallel part */ 67674ae819SStefano Zampini for (i=0;i<deluxe_ctx->par_colors;i++) { 68674ae819SStefano Zampini if (deluxe_ctx->par_ksp[i]) { 6934a97f8cSStefano Zampini PetscMPIInt color_rank; 7034a97f8cSStefano Zampini PetscInt subidx = deluxe_ctx->par_col2sub[i]; 7134a97f8cSStefano Zampini /* restrict on subset */ 7234a97f8cSStefano Zampini ierr = VecScatterBegin(deluxe_ctx->par_scctx_s[i],x,sub_schurs->work1[subidx],INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 7334a97f8cSStefano Zampini ierr = VecScatterEnd(deluxe_ctx->par_scctx_s[i],x,sub_schurs->work1[subidx],INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 7434a97f8cSStefano Zampini /* S_Ej */ 7534a97f8cSStefano Zampini ierr = MatMult(sub_schurs->S_Ej[subidx],sub_schurs->work1[subidx],sub_schurs->work2[subidx]);CHKERRQ(ierr); 7634a97f8cSStefano Zampini /* (\sum_j S_Ej)^-1 */ 7734a97f8cSStefano Zampini ierr = VecSet(deluxe_ctx->par_vec[i],0.0);CHKERRQ(ierr); 7834a97f8cSStefano Zampini ierr = VecScatterBegin(deluxe_ctx->par_scctx_p[i],sub_schurs->work2[subidx],deluxe_ctx->par_vec[i],ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 7934a97f8cSStefano Zampini ierr = VecScatterEnd(deluxe_ctx->par_scctx_p[i],sub_schurs->work2[subidx],deluxe_ctx->par_vec[i],ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 80674ae819SStefano Zampini ierr = KSPSolve(deluxe_ctx->par_ksp[i],deluxe_ctx->par_vec[i],deluxe_ctx->par_vec[i]);CHKERRQ(ierr); 8134a97f8cSStefano Zampini ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)deluxe_ctx->par_ksp[i]),&color_rank);CHKERRQ(ierr); 8234a97f8cSStefano Zampini /* get back solution on subset */ 8334a97f8cSStefano Zampini ierr = VecScatterBegin(deluxe_ctx->par_scctx_p[i],deluxe_ctx->par_vec[i],sub_schurs->work1[subidx],INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 8434a97f8cSStefano Zampini ierr = VecScatterEnd(deluxe_ctx->par_scctx_p[i],deluxe_ctx->par_vec[i],sub_schurs->work1[subidx],INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 8534a97f8cSStefano Zampini if (!color_rank) { /* only the master process in coloured comm copies the computed values */ 8634a97f8cSStefano Zampini ierr = VecScatterBegin(deluxe_ctx->par_scctx_s[i],sub_schurs->work1[subidx],pcbddc->work_scaling,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 8734a97f8cSStefano Zampini ierr = VecScatterEnd(deluxe_ctx->par_scctx_s[i],sub_schurs->work1[subidx],pcbddc->work_scaling,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 88674ae819SStefano Zampini } 89674ae819SStefano Zampini } 90674ae819SStefano Zampini } 91674ae819SStefano Zampini /* put local boundary part in global vector */ 9234a97f8cSStefano Zampini ierr = VecSet(y,0.0);CHKERRQ(ierr); 93674ae819SStefano Zampini ierr = VecScatterBegin(pcis->global_to_B,pcbddc->work_scaling,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 94674ae819SStefano Zampini ierr = VecScatterEnd(pcis->global_to_B,pcbddc->work_scaling,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 95674ae819SStefano Zampini PetscFunctionReturn(0); 96674ae819SStefano Zampini } 97674ae819SStefano Zampini 98674ae819SStefano Zampini #undef __FUNCT__ 99674ae819SStefano Zampini #define __FUNCT__ "PCBDDCScalingExtension" 100674ae819SStefano Zampini PetscErrorCode PCBDDCScalingExtension(PC pc, Vec local_interface_vector, Vec global_vector) 101674ae819SStefano Zampini { 102674ae819SStefano Zampini PC_BDDC *pcbddc=(PC_BDDC*)pc->data; 103674ae819SStefano Zampini PetscErrorCode ierr; 104674ae819SStefano Zampini 105674ae819SStefano Zampini PetscFunctionBegin; 106674ae819SStefano Zampini PetscValidHeaderSpecific(pc,PC_CLASSID,1); 107674ae819SStefano Zampini PetscValidHeaderSpecific(local_interface_vector,VEC_CLASSID,2); 108674ae819SStefano Zampini PetscValidHeaderSpecific(global_vector,VEC_CLASSID,3); 109674ae819SStefano Zampini if (local_interface_vector == pcbddc->work_scaling) { 110674ae819SStefano Zampini SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"Local vector cannot be pcbddc->work_scaling!\n"); 111674ae819SStefano Zampini } 112674ae819SStefano Zampini ierr = PetscTryMethod(pc,"PCBDDCScalingExtension_C",(PC,Vec,Vec),(pc,local_interface_vector,global_vector));CHKERRQ(ierr); 113674ae819SStefano Zampini PetscFunctionReturn(0); 114674ae819SStefano Zampini } 115674ae819SStefano Zampini 116674ae819SStefano Zampini #undef __FUNCT__ 117674ae819SStefano Zampini #define __FUNCT__ "PCBDDCScalingRestriction_Basic" 118674ae819SStefano Zampini static PetscErrorCode PCBDDCScalingRestriction_Basic(PC pc, Vec global_vector, Vec local_interface_vector) 119674ae819SStefano Zampini { 120674ae819SStefano Zampini PetscErrorCode ierr; 121674ae819SStefano Zampini PC_IS* pcis = (PC_IS*)pc->data; 122674ae819SStefano Zampini 123674ae819SStefano Zampini PetscFunctionBegin; 124674ae819SStefano Zampini ierr = VecScatterBegin(pcis->global_to_B,global_vector,local_interface_vector,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 125674ae819SStefano Zampini ierr = VecScatterEnd(pcis->global_to_B,global_vector,local_interface_vector,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 126674ae819SStefano Zampini /* Apply partition of unity */ 127674ae819SStefano Zampini ierr = VecPointwiseMult(local_interface_vector,pcis->D,local_interface_vector);CHKERRQ(ierr); 128674ae819SStefano Zampini PetscFunctionReturn(0); 129674ae819SStefano Zampini } 130674ae819SStefano Zampini 131674ae819SStefano Zampini #undef __FUNCT__ 132674ae819SStefano Zampini #define __FUNCT__ "PCBDDCScalingRestriction_Deluxe" 133674ae819SStefano Zampini static PetscErrorCode PCBDDCScalingRestriction_Deluxe(PC pc, Vec x, Vec y) 134674ae819SStefano Zampini { 135674ae819SStefano Zampini PC_IS* pcis=(PC_IS*)pc->data; 136674ae819SStefano Zampini PC_BDDC* pcbddc=(PC_BDDC*)pc->data; 137674ae819SStefano Zampini PCBDDCDeluxeScaling deluxe_ctx = pcbddc->deluxe_ctx; 138b1b3d7a2SStefano Zampini PCBDDCSubSchurs sub_schurs = pcbddc->sub_schurs[1]; 139674ae819SStefano Zampini PetscInt i; 140674ae819SStefano Zampini PetscErrorCode ierr; 141674ae819SStefano Zampini 142674ae819SStefano Zampini PetscFunctionBegin; 143674ae819SStefano Zampini /* get local boundary part of global vector */ 144674ae819SStefano Zampini ierr = VecScatterBegin(pcis->global_to_B,x,y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 145674ae819SStefano Zampini ierr = VecScatterEnd(pcis->global_to_B,x,y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 14634a97f8cSStefano Zampini if (deluxe_ctx->n_simple) { 14734a97f8cSStefano Zampini /* scale deluxe vertices using diagonal scaling */ 14834a97f8cSStefano Zampini PetscScalar *array_y,*array_D; 149674ae819SStefano Zampini ierr = VecGetArray(y,&array_y);CHKERRQ(ierr); 150674ae819SStefano Zampini ierr = VecGetArray(pcis->D,&array_D);CHKERRQ(ierr); 151674ae819SStefano Zampini for (i=0;i<deluxe_ctx->n_simple;i++) { 152674ae819SStefano Zampini array_y[deluxe_ctx->idx_simple_B[i]] *= array_D[deluxe_ctx->idx_simple_B[i]]; 153674ae819SStefano Zampini } 154674ae819SStefano Zampini ierr = VecRestoreArray(pcis->D,&array_D);CHKERRQ(ierr); 155674ae819SStefano Zampini ierr = VecRestoreArray(y,&array_y);CHKERRQ(ierr); 15634a97f8cSStefano Zampini } 15734a97f8cSStefano Zampini /* sequential part : all problems and Schur applications collapsed into a single matrix vector multiplication and ksp solution */ 158*41c3ba1bSStefano Zampini if (deluxe_ctx->seq_ksp) { 159674ae819SStefano Zampini ierr = VecScatterBegin(deluxe_ctx->seq_scctx,y,deluxe_ctx->seq_work1,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 160674ae819SStefano Zampini ierr = VecScatterEnd(deluxe_ctx->seq_scctx,y,deluxe_ctx->seq_work1,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 161674ae819SStefano Zampini ierr = KSPSolveTranspose(deluxe_ctx->seq_ksp,deluxe_ctx->seq_work1,deluxe_ctx->seq_work2);CHKERRQ(ierr); 162*41c3ba1bSStefano Zampini ierr = MatMultTranspose(sub_schurs->S_Ej_all,deluxe_ctx->seq_work2,deluxe_ctx->seq_work1);CHKERRQ(ierr); 163674ae819SStefano Zampini ierr = VecScatterBegin(deluxe_ctx->seq_scctx,deluxe_ctx->seq_work1,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 164674ae819SStefano Zampini ierr = VecScatterEnd(deluxe_ctx->seq_scctx,deluxe_ctx->seq_work1,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 165674ae819SStefano Zampini } 166674ae819SStefano Zampini /* parallel part */ 167674ae819SStefano Zampini for (i=0;i<deluxe_ctx->par_colors;i++) { 168674ae819SStefano Zampini if (deluxe_ctx->par_ksp[i]) { 16934a97f8cSStefano Zampini PetscInt subidx = deluxe_ctx->par_col2sub[i]; 17034a97f8cSStefano Zampini /* restrict on subset */ 17134a97f8cSStefano Zampini ierr = VecScatterBegin(deluxe_ctx->par_scctx_s[i],y,sub_schurs->work1[subidx],INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 17234a97f8cSStefano Zampini ierr = VecScatterEnd(deluxe_ctx->par_scctx_s[i],y,sub_schurs->work1[subidx],INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 17334a97f8cSStefano Zampini /* (\sum_j S_Ej)^-T */ 17434a97f8cSStefano Zampini ierr = VecSet(deluxe_ctx->par_vec[i],0.0);CHKERRQ(ierr); 17534a97f8cSStefano Zampini ierr = VecScatterBegin(deluxe_ctx->par_scctx_p[i],sub_schurs->work1[subidx],deluxe_ctx->par_vec[i],INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 17634a97f8cSStefano Zampini ierr = VecScatterEnd(deluxe_ctx->par_scctx_p[i],sub_schurs->work1[subidx],deluxe_ctx->par_vec[i],INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 177674ae819SStefano Zampini ierr = KSPSolveTranspose(deluxe_ctx->par_ksp[i],deluxe_ctx->par_vec[i],deluxe_ctx->par_vec[i]);CHKERRQ(ierr); 17834a97f8cSStefano Zampini ierr = VecScatterBegin(deluxe_ctx->par_scctx_p[i],deluxe_ctx->par_vec[i],sub_schurs->work1[subidx],INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 17934a97f8cSStefano Zampini ierr = VecScatterEnd(deluxe_ctx->par_scctx_p[i],deluxe_ctx->par_vec[i],sub_schurs->work1[subidx],INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 18034a97f8cSStefano Zampini /* S_Ej^T */ 18134a97f8cSStefano Zampini ierr = MatMultTranspose(sub_schurs->S_Ej[subidx],sub_schurs->work1[subidx],sub_schurs->work2[subidx]);CHKERRQ(ierr); 18234a97f8cSStefano Zampini /* extend to boundary */ 18334a97f8cSStefano Zampini ierr = VecScatterBegin(deluxe_ctx->par_scctx_s[i],sub_schurs->work2[subidx],y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 18434a97f8cSStefano Zampini ierr = VecScatterEnd(deluxe_ctx->par_scctx_s[i],sub_schurs->work2[subidx],y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 185674ae819SStefano Zampini } 186674ae819SStefano Zampini } 187674ae819SStefano Zampini PetscFunctionReturn(0); 188674ae819SStefano Zampini } 189674ae819SStefano Zampini 190674ae819SStefano Zampini #undef __FUNCT__ 191674ae819SStefano Zampini #define __FUNCT__ "PCBDDCScalingRestriction" 192674ae819SStefano Zampini PetscErrorCode PCBDDCScalingRestriction(PC pc, Vec global_vector, Vec local_interface_vector) 193674ae819SStefano Zampini { 194674ae819SStefano Zampini PC_BDDC *pcbddc=(PC_BDDC*)pc->data; 195674ae819SStefano Zampini PetscErrorCode ierr; 196674ae819SStefano Zampini 197674ae819SStefano Zampini PetscFunctionBegin; 198674ae819SStefano Zampini PetscValidHeaderSpecific(pc,PC_CLASSID,1); 199674ae819SStefano Zampini PetscValidHeaderSpecific(global_vector,VEC_CLASSID,2); 200674ae819SStefano Zampini PetscValidHeaderSpecific(local_interface_vector,VEC_CLASSID,3); 201674ae819SStefano Zampini if (local_interface_vector == pcbddc->work_scaling) { 202674ae819SStefano Zampini SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"Local vector should cannot be pcbddc->work_scaling!\n"); 203674ae819SStefano Zampini } 204674ae819SStefano Zampini ierr = PetscTryMethod(pc,"PCBDDCScalingRestriction_C",(PC,Vec,Vec),(pc,global_vector,local_interface_vector));CHKERRQ(ierr); 205674ae819SStefano Zampini PetscFunctionReturn(0); 206674ae819SStefano Zampini } 207674ae819SStefano Zampini 208674ae819SStefano Zampini #undef __FUNCT__ 209674ae819SStefano Zampini #define __FUNCT__ "PCBDDCScalingSetUp" 210674ae819SStefano Zampini PetscErrorCode PCBDDCScalingSetUp(PC pc) 211674ae819SStefano Zampini { 212674ae819SStefano Zampini PC_IS* pcis=(PC_IS*)pc->data; 213674ae819SStefano Zampini PC_BDDC* pcbddc=(PC_BDDC*)pc->data; 214674ae819SStefano Zampini PetscErrorCode ierr; 215674ae819SStefano Zampini 216674ae819SStefano Zampini PetscFunctionBegin; 217674ae819SStefano Zampini PetscValidHeaderSpecific(pc,PC_CLASSID,1); 218674ae819SStefano Zampini /* create work vector for the operator */ 21934a97f8cSStefano Zampini ierr = VecDestroy(&pcbddc->work_scaling);CHKERRQ(ierr); 220674ae819SStefano Zampini ierr = VecDuplicate(pcis->vec1_B,&pcbddc->work_scaling);CHKERRQ(ierr); 22134a97f8cSStefano Zampini /* always rebuild pcis->D */ 22228d874f6SStefano Zampini if (pcis->use_stiffness_scaling) { 223674ae819SStefano Zampini ierr = MatGetDiagonal(pcbddc->local_mat,pcis->vec1_N);CHKERRQ(ierr); 224674ae819SStefano Zampini ierr = VecScatterBegin(pcis->N_to_B,pcis->vec1_N,pcis->D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 225674ae819SStefano Zampini ierr = VecScatterEnd(pcis->N_to_B,pcis->vec1_N,pcis->D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 226674ae819SStefano Zampini } 227674ae819SStefano Zampini ierr = VecCopy(pcis->D,pcis->vec1_B);CHKERRQ(ierr); 228674ae819SStefano Zampini ierr = VecSet(pcis->vec1_global,0.0);CHKERRQ(ierr); 229674ae819SStefano Zampini ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_B,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 230674ae819SStefano Zampini ierr = VecScatterEnd(pcis->global_to_B,pcis->vec1_B,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 231674ae819SStefano Zampini ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_global,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 232674ae819SStefano Zampini ierr = VecScatterEnd(pcis->global_to_B,pcis->vec1_global,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 233674ae819SStefano Zampini ierr = VecPointwiseDivide(pcis->D,pcis->D,pcis->vec1_B);CHKERRQ(ierr); 234674ae819SStefano Zampini /* now setup */ 235681e7c04SStefano Zampini if (pcbddc->use_deluxe_scaling) { 23634a97f8cSStefano Zampini if (!pcbddc->deluxe_ctx) { 23734a97f8cSStefano Zampini ierr = PCBDDCScalingCreate_Deluxe(pc);CHKERRQ(ierr); 23834a97f8cSStefano Zampini } 23934a97f8cSStefano Zampini ierr = PCBDDCScalingSetUp_Deluxe(pc);CHKERRQ(ierr); 240674ae819SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCScalingRestriction_C",PCBDDCScalingRestriction_Deluxe);CHKERRQ(ierr); 241674ae819SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCScalingExtension_C",PCBDDCScalingExtension_Deluxe);CHKERRQ(ierr); 242674ae819SStefano Zampini } else { 243674ae819SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCScalingRestriction_C",PCBDDCScalingRestriction_Basic);CHKERRQ(ierr); 244674ae819SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCScalingExtension_C",PCBDDCScalingExtension_Basic);CHKERRQ(ierr); 245674ae819SStefano Zampini } 24634a97f8cSStefano Zampini 247674ae819SStefano Zampini /* test */ 248674ae819SStefano Zampini if (pcbddc->dbg_flag) { 24934a97f8cSStefano Zampini Vec vec2_global; 250674ae819SStefano Zampini PetscViewer viewer=pcbddc->dbg_viewer; 25134a97f8cSStefano Zampini PetscReal error; 252674ae819SStefano Zampini 253674ae819SStefano Zampini /* extension -> from local to parallel */ 25434a97f8cSStefano Zampini ierr = VecSet(pcis->vec1_global,0.0);CHKERRQ(ierr); 25534a97f8cSStefano Zampini ierr = VecSetRandom(pcis->vec1_B,NULL);CHKERRQ(ierr); 25634a97f8cSStefano Zampini ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_B,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 25734a97f8cSStefano Zampini ierr = VecScatterEnd(pcis->global_to_B,pcis->vec1_B,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 25834a97f8cSStefano Zampini ierr = VecDuplicate(pcis->vec1_global,&vec2_global);CHKERRQ(ierr); 25934a97f8cSStefano Zampini ierr = VecCopy(pcis->vec1_global,vec2_global);CHKERRQ(ierr); 26034a97f8cSStefano Zampini 261674ae819SStefano Zampini ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_global,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 262674ae819SStefano Zampini ierr = VecScatterEnd(pcis->global_to_B,pcis->vec1_global,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 263674ae819SStefano Zampini ierr = PCBDDCScalingExtension(pc,pcis->vec1_B,pcis->vec1_global);CHKERRQ(ierr); 26434a97f8cSStefano Zampini ierr = VecAXPY(pcis->vec1_global,-1.0,vec2_global);CHKERRQ(ierr); 26534a97f8cSStefano Zampini ierr = VecNorm(pcis->vec1_global,NORM_INFINITY,&error);CHKERRQ(ierr); 266674ae819SStefano Zampini ierr = PetscViewerASCIIPrintf(viewer,"Error scaling extension %1.14e\n",error);CHKERRQ(ierr); 26734a97f8cSStefano Zampini if (error>1.e-8 && pcbddc->dbg_flag>1) { 268674ae819SStefano Zampini ierr = VecView(pcis->vec1_global,viewer);CHKERRQ(ierr); 269674ae819SStefano Zampini } 27034a97f8cSStefano Zampini ierr = VecDestroy(&vec2_global);CHKERRQ(ierr); 27134a97f8cSStefano Zampini 272674ae819SStefano Zampini /* restriction -> from parallel to local */ 273674ae819SStefano Zampini ierr = VecSet(pcis->vec1_global,0.0);CHKERRQ(ierr); 27434a97f8cSStefano Zampini ierr = VecSetRandom(pcis->vec1_B,NULL);CHKERRQ(ierr); 275674ae819SStefano Zampini ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_B,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 276674ae819SStefano Zampini ierr = VecScatterEnd(pcis->global_to_B,pcis->vec1_B,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 27734a97f8cSStefano Zampini 27834a97f8cSStefano Zampini ierr = PCBDDCScalingRestriction(pc,pcis->vec1_global,pcis->vec1_B);CHKERRQ(ierr); 27934a97f8cSStefano Zampini ierr = VecScale(pcis->vec1_B,-1.0);CHKERRQ(ierr); 28034a97f8cSStefano Zampini ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_B,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 28134a97f8cSStefano Zampini ierr = VecScatterEnd(pcis->global_to_B,pcis->vec1_B,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 28234a97f8cSStefano Zampini ierr = VecNorm(pcis->vec1_global,NORM_INFINITY,&error);CHKERRQ(ierr); 28334a97f8cSStefano Zampini ierr = PetscViewerASCIIPrintf(viewer,"Error scaling restriction %1.14e\n",error);CHKERRQ(ierr); 28434a97f8cSStefano Zampini if (error>1.e-8 && pcbddc->dbg_flag>1) { 285674ae819SStefano Zampini ierr = VecView(pcis->vec1_global,viewer);CHKERRQ(ierr); 286674ae819SStefano Zampini } 287674ae819SStefano Zampini } 288674ae819SStefano Zampini PetscFunctionReturn(0); 289674ae819SStefano Zampini } 290674ae819SStefano Zampini 291674ae819SStefano Zampini #undef __FUNCT__ 292674ae819SStefano Zampini #define __FUNCT__ "PCBDDCScalingDestroy" 293674ae819SStefano Zampini PetscErrorCode PCBDDCScalingDestroy(PC pc) 294674ae819SStefano Zampini { 295674ae819SStefano Zampini PC_BDDC* pcbddc=(PC_BDDC*)pc->data; 296674ae819SStefano Zampini PetscErrorCode ierr; 297674ae819SStefano Zampini 298674ae819SStefano Zampini PetscFunctionBegin; 29934a97f8cSStefano Zampini if (pcbddc->deluxe_ctx) { 30034a97f8cSStefano Zampini ierr = PCBDDCScalingDestroy_Deluxe(pc);CHKERRQ(ierr); 301674ae819SStefano Zampini } 302674ae819SStefano Zampini ierr = VecDestroy(&pcbddc->work_scaling);CHKERRQ(ierr); 303674ae819SStefano Zampini /* remove functions */ 304674ae819SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCScalingRestriction_C",NULL);CHKERRQ(ierr); 305674ae819SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCScalingExtension_C",NULL);CHKERRQ(ierr); 306674ae819SStefano Zampini PetscFunctionReturn(0); 307674ae819SStefano Zampini } 308674ae819SStefano Zampini 30934a97f8cSStefano Zampini #undef __FUNCT__ 31034a97f8cSStefano Zampini #define __FUNCT__ "PCBDDCScalingCreate_Deluxe" 31134a97f8cSStefano Zampini static PetscErrorCode PCBDDCScalingCreate_Deluxe(PC pc) 31234a97f8cSStefano Zampini { 31334a97f8cSStefano Zampini PC_BDDC* pcbddc=(PC_BDDC*)pc->data; 31434a97f8cSStefano Zampini PCBDDCDeluxeScaling deluxe_ctx; 31534a97f8cSStefano Zampini PetscErrorCode ierr; 31634a97f8cSStefano Zampini 31734a97f8cSStefano Zampini PetscFunctionBegin; 31834a97f8cSStefano Zampini ierr = PetscNew(&deluxe_ctx);CHKERRQ(ierr); 31934a97f8cSStefano Zampini pcbddc->deluxe_ctx = deluxe_ctx; 32034a97f8cSStefano Zampini PetscFunctionReturn(0); 32134a97f8cSStefano Zampini } 32234a97f8cSStefano Zampini 32334a97f8cSStefano Zampini #undef __FUNCT__ 32434a97f8cSStefano Zampini #define __FUNCT__ "PCBDDCScalingDestroy_Deluxe" 32534a97f8cSStefano Zampini static PetscErrorCode PCBDDCScalingDestroy_Deluxe(PC pc) 32634a97f8cSStefano Zampini { 32734a97f8cSStefano Zampini PC_BDDC* pcbddc=(PC_BDDC*)pc->data; 32834a97f8cSStefano Zampini PetscErrorCode ierr; 32934a97f8cSStefano Zampini 33034a97f8cSStefano Zampini PetscFunctionBegin; 33134a97f8cSStefano Zampini ierr = PCBDDCScalingReset_Deluxe_Solvers(pcbddc->deluxe_ctx);CHKERRQ(ierr); 33234a97f8cSStefano Zampini ierr = PetscFree(pcbddc->deluxe_ctx);CHKERRQ(ierr); 33334a97f8cSStefano Zampini PetscFunctionReturn(0); 33434a97f8cSStefano Zampini } 33534a97f8cSStefano Zampini 33634a97f8cSStefano Zampini #undef __FUNCT__ 33734a97f8cSStefano Zampini #define __FUNCT__ "PCBDDCScalingReset_Deluxe_Solvers" 33834a97f8cSStefano Zampini static PetscErrorCode PCBDDCScalingReset_Deluxe_Solvers(PCBDDCDeluxeScaling deluxe_ctx) 33934a97f8cSStefano Zampini { 34034a97f8cSStefano Zampini PetscErrorCode ierr; 34134a97f8cSStefano Zampini 34234a97f8cSStefano Zampini PetscFunctionBegin; 34334a97f8cSStefano Zampini ierr = PetscFree(deluxe_ctx->idx_simple_B);CHKERRQ(ierr); 34434a97f8cSStefano Zampini deluxe_ctx->n_simple = 0; 345*41c3ba1bSStefano Zampini if (deluxe_ctx->seq_ksp) { 34634a97f8cSStefano Zampini ierr = VecScatterDestroy(&deluxe_ctx->seq_scctx);CHKERRQ(ierr); 34734a97f8cSStefano Zampini ierr = VecDestroy(&deluxe_ctx->seq_work1);CHKERRQ(ierr); 34834a97f8cSStefano Zampini ierr = VecDestroy(&deluxe_ctx->seq_work2);CHKERRQ(ierr); 34934a97f8cSStefano Zampini ierr = KSPDestroy(&deluxe_ctx->seq_ksp);CHKERRQ(ierr); 35034a97f8cSStefano Zampini } 35134a97f8cSStefano Zampini if (deluxe_ctx->par_colors) { 35234a97f8cSStefano Zampini PetscInt i; 35334a97f8cSStefano Zampini for (i=0;i<deluxe_ctx->par_colors;i++) { 35434a97f8cSStefano Zampini ierr = VecScatterDestroy(&deluxe_ctx->par_scctx_s[i]);CHKERRQ(ierr); 35534a97f8cSStefano Zampini ierr = VecScatterDestroy(&deluxe_ctx->par_scctx_p[i]);CHKERRQ(ierr); 35634a97f8cSStefano Zampini ierr = VecDestroy(&deluxe_ctx->par_vec[i]);CHKERRQ(ierr); 35734a97f8cSStefano Zampini ierr = KSPDestroy(&deluxe_ctx->par_ksp[i]);CHKERRQ(ierr); 35834a97f8cSStefano Zampini } 35934a97f8cSStefano Zampini ierr = PetscFree5(deluxe_ctx->par_ksp, 36034a97f8cSStefano Zampini deluxe_ctx->par_scctx_s, 36134a97f8cSStefano Zampini deluxe_ctx->par_scctx_p, 36234a97f8cSStefano Zampini deluxe_ctx->par_vec, 36334a97f8cSStefano Zampini deluxe_ctx->par_col2sub);CHKERRQ(ierr); 36434a97f8cSStefano Zampini } 36534a97f8cSStefano Zampini deluxe_ctx->par_colors = 0; 36634a97f8cSStefano Zampini PetscFunctionReturn(0); 36734a97f8cSStefano Zampini } 36834a97f8cSStefano Zampini 369812aaec9SStefano Zampini #define OLD_CODE 0 37034a97f8cSStefano Zampini #undef __FUNCT__ 37134a97f8cSStefano Zampini #define __FUNCT__ "PCBDDCScalingSetUp_Deluxe" 37234a97f8cSStefano Zampini static PetscErrorCode PCBDDCScalingSetUp_Deluxe(PC pc) 37334a97f8cSStefano Zampini { 37434a97f8cSStefano Zampini PC_IS *pcis=(PC_IS*)pc->data; 37534a97f8cSStefano Zampini PC_BDDC *pcbddc=(PC_BDDC*)pc->data; 37634a97f8cSStefano Zampini PCBDDCDeluxeScaling deluxe_ctx=pcbddc->deluxe_ctx; 377b1b3d7a2SStefano Zampini PCBDDCSubSchurs sub_schurs=pcbddc->sub_schurs[1]; 37834a97f8cSStefano Zampini PCBDDCGraph graph; 3797821e9e7SStefano Zampini PetscBT bitmask; 380b1b3d7a2SStefano Zampini #if OLD_CODE 381b1b3d7a2SStefano Zampini IS *faces,*edges,*all_cc; 38234a97f8cSStefano Zampini PetscInt *index_sequential,*index_parallel; 38334a97f8cSStefano Zampini PetscInt *auxlocal_sequential,*auxlocal_parallel; 38434a97f8cSStefano Zampini PetscInt *auxglobal_sequential,*auxglobal_parallel; 3857821e9e7SStefano Zampini PetscInt *auxmapping,*idxs; 3861cef56d8SStefano Zampini PetscInt i,max_subset_size; 38734a97f8cSStefano Zampini PetscInt n_sequential_problems,n_local_sequential_problems,n_parallel_problems,n_local_parallel_problems; 3887821e9e7SStefano Zampini PetscInt n_faces,n_edges,n_all_cc; 389b1b3d7a2SStefano Zampini #else 390b1b3d7a2SStefano Zampini PetscInt i; 391b1b3d7a2SStefano Zampini const PetscInt* idxs; 392b1b3d7a2SStefano Zampini Mat S_j; 393b1b3d7a2SStefano Zampini PetscBool free_used_adj; 394b1b3d7a2SStefano Zampini PetscInt *used_xadj,*used_adjncy; 395b1b3d7a2SStefano Zampini #endif 39634a97f8cSStefano Zampini PetscErrorCode ierr; 39734a97f8cSStefano Zampini 39834a97f8cSStefano Zampini PetscFunctionBegin; 39934a97f8cSStefano Zampini /* throw away the solvers */ 40034a97f8cSStefano Zampini ierr = PCBDDCScalingReset_Deluxe_Solvers(deluxe_ctx);CHKERRQ(ierr); 40134a97f8cSStefano Zampini 40234a97f8cSStefano Zampini /* attach interface graph for determining subsets */ 40334a97f8cSStefano Zampini if (pcbddc->deluxe_rebuild) { /* in case rebuild has been requested, it uses a graph generated only by the neighbouring information */ 40434a97f8cSStefano Zampini PetscInt *idx_V_N; 40534a97f8cSStefano Zampini IS verticesIS; 40634a97f8cSStefano Zampini ierr = PCBDDCGetPrimalVerticesLocalIdx(pc,&i,&idx_V_N);CHKERRQ(ierr); 40734a97f8cSStefano Zampini ierr = ISCreateGeneral(PETSC_COMM_SELF,i,idx_V_N,PETSC_OWN_POINTER,&verticesIS);CHKERRQ(ierr); 40834a97f8cSStefano Zampini ierr = PCBDDCGraphCreate(&graph);CHKERRQ(ierr); 40934a97f8cSStefano Zampini ierr = PCBDDCGraphInit(graph,pcbddc->mat_graph->l2gmap);CHKERRQ(ierr); 4100a95f60dSStefano Zampini ierr = PCBDDCGraphSetUp(graph,0,NULL,pcbddc->DirichletBoundariesLocal,0,NULL,verticesIS);CHKERRQ(ierr); 41134a97f8cSStefano Zampini ierr = PCBDDCGraphComputeConnectedComponents(graph);CHKERRQ(ierr); 41234a97f8cSStefano Zampini ierr = ISDestroy(&verticesIS);CHKERRQ(ierr); 41334a97f8cSStefano Zampini /* 41434a97f8cSStefano Zampini if (pcbddc->dbg_flag) { 41534a97f8cSStefano Zampini ierr = PCBDDCGraphASCIIView(graph,pcbddc->dbg_flag,pcbddc->dbg_viewer);CHKERRQ(ierr); 41634a97f8cSStefano Zampini } 41734a97f8cSStefano Zampini */ 41834a97f8cSStefano Zampini } else { 41934a97f8cSStefano Zampini graph = pcbddc->mat_graph; 42034a97f8cSStefano Zampini } 42134a97f8cSStefano Zampini 422b1b3d7a2SStefano Zampini #if OLD_CODE 4237821e9e7SStefano Zampini /* get index sets for faces and edges */ 4247821e9e7SStefano Zampini ierr = PCBDDCGraphGetCandidatesIS(graph,&n_faces,&faces,&n_edges,&edges,NULL);CHKERRQ(ierr); 4257821e9e7SStefano Zampini n_all_cc = n_faces+n_edges; 4267821e9e7SStefano Zampini ierr = PetscMalloc1(n_all_cc,&all_cc);CHKERRQ(ierr); 4277821e9e7SStefano Zampini for (i=0;i<n_faces;i++) { 4287821e9e7SStefano Zampini all_cc[i] = faces[i]; 4297821e9e7SStefano Zampini } 4307821e9e7SStefano Zampini for (i=0;i<n_edges;i++) { 4317821e9e7SStefano Zampini all_cc[n_faces+i] = edges[i]; 4327821e9e7SStefano Zampini } 433310ea0faSStefano Zampini ierr = PetscFree(faces);CHKERRQ(ierr); 434310ea0faSStefano Zampini ierr = PetscFree(edges);CHKERRQ(ierr); 4357821e9e7SStefano Zampini 43634a97f8cSStefano Zampini /* map interface's subsets */ 4371cef56d8SStefano Zampini max_subset_size = 0; 4387821e9e7SStefano Zampini for (i=0;i<n_all_cc;i++) { 4397821e9e7SStefano Zampini PetscInt subset_size; 4407821e9e7SStefano Zampini ierr = ISGetLocalSize(all_cc[i],&subset_size);CHKERRQ(ierr); 4417821e9e7SStefano Zampini max_subset_size = PetscMax(max_subset_size,subset_size); 44234a97f8cSStefano Zampini } 4431cef56d8SStefano Zampini ierr = PetscMalloc5(max_subset_size,&auxmapping, 44434a97f8cSStefano Zampini graph->ncc,&auxlocal_sequential, 44534a97f8cSStefano Zampini graph->ncc,&auxlocal_parallel, 44634a97f8cSStefano Zampini graph->ncc,&index_sequential, 44734a97f8cSStefano Zampini graph->ncc,&index_parallel);CHKERRQ(ierr); 44834a97f8cSStefano Zampini 4491cef56d8SStefano Zampini /* if threshold is negative, uses all sequential problems */ 4501cef56d8SStefano Zampini if (pcbddc->deluxe_threshold < 0) pcbddc->deluxe_threshold = max_subset_size; 4511cef56d8SStefano Zampini 4527821e9e7SStefano Zampini /* workspace */ 4537821e9e7SStefano Zampini ierr = PetscBTCreate(pcis->n,&bitmask);CHKERRQ(ierr); 4547821e9e7SStefano Zampini ierr = ISGetIndices(pcis->is_I_local,(const PetscInt**)&idxs);CHKERRQ(ierr); 4557821e9e7SStefano Zampini for (i=0;i<pcis->n-pcis->n_B;i++) { 4567821e9e7SStefano Zampini ierr = PetscBTSet(bitmask,idxs[i]);CHKERRQ(ierr); 4577821e9e7SStefano Zampini } 4587821e9e7SStefano Zampini ierr = ISRestoreIndices(pcis->is_I_local,(const PetscInt**)&idxs);CHKERRQ(ierr); 4597821e9e7SStefano Zampini 4607821e9e7SStefano Zampini /* determine which problem has to be solved in parallel or sequentially */ 46134a97f8cSStefano Zampini n_local_sequential_problems = 0; 46234a97f8cSStefano Zampini n_local_parallel_problems = 0; 4637821e9e7SStefano Zampini for (i=0;i<n_all_cc;i++) { 4647821e9e7SStefano Zampini PetscInt subset_size,j,min_loc = 0; 4657821e9e7SStefano Zampini 4667821e9e7SStefano Zampini ierr = ISGetLocalSize(all_cc[i],&subset_size);CHKERRQ(ierr); 4677821e9e7SStefano Zampini ierr = ISGetIndices(all_cc[i],(const PetscInt**)&idxs);CHKERRQ(ierr); 4687821e9e7SStefano Zampini for (j=0;j<subset_size;j++) { 4697821e9e7SStefano Zampini ierr = PetscBTSet(bitmask,idxs[j]);CHKERRQ(ierr); 4707821e9e7SStefano Zampini } 4717821e9e7SStefano Zampini ierr = ISLocalToGlobalMappingApply(graph->l2gmap,subset_size,idxs,auxmapping);CHKERRQ(ierr); 47234a97f8cSStefano Zampini for (j=1;j<subset_size;j++) { 47334a97f8cSStefano Zampini if (auxmapping[j]<auxmapping[min_loc]) { 47434a97f8cSStefano Zampini min_loc = j; 47534a97f8cSStefano Zampini } 47634a97f8cSStefano Zampini } 47734a97f8cSStefano Zampini if (subset_size > pcbddc->deluxe_threshold) { 47834a97f8cSStefano Zampini index_parallel[n_local_parallel_problems] = i; 4797821e9e7SStefano Zampini auxlocal_parallel[n_local_parallel_problems] = idxs[min_loc]; 48034a97f8cSStefano Zampini n_local_parallel_problems++; 48134a97f8cSStefano Zampini } else { 48234a97f8cSStefano Zampini index_sequential[n_local_sequential_problems] = i; 4837821e9e7SStefano Zampini auxlocal_sequential[n_local_sequential_problems] = idxs[min_loc]; 48434a97f8cSStefano Zampini n_local_sequential_problems++; 48534a97f8cSStefano Zampini } 4867821e9e7SStefano Zampini ierr = ISRestoreIndices(all_cc[i],(const PetscInt**)&idxs);CHKERRQ(ierr); 48734a97f8cSStefano Zampini } 48834a97f8cSStefano Zampini 4897821e9e7SStefano Zampini /* diagonal scaling on interface dofs not contained in cc */ 4907821e9e7SStefano Zampini deluxe_ctx->n_simple = 0; 4917821e9e7SStefano Zampini for (i=0;i<pcis->n;i++) { 4927821e9e7SStefano Zampini if (!PetscBTLookup(bitmask,i)) { 49334a97f8cSStefano Zampini deluxe_ctx->n_simple++; 49434a97f8cSStefano Zampini } 49534a97f8cSStefano Zampini } 49634a97f8cSStefano Zampini ierr = PetscMalloc1(deluxe_ctx->n_simple,&deluxe_ctx->idx_simple_B);CHKERRQ(ierr); 49734a97f8cSStefano Zampini deluxe_ctx->n_simple = 0; 4987821e9e7SStefano Zampini for (i=0;i<pcis->n;i++) { 4997821e9e7SStefano Zampini if (!PetscBTLookup(bitmask,i)) { 50034a97f8cSStefano Zampini deluxe_ctx->idx_simple_B[deluxe_ctx->n_simple++] = i; 50134a97f8cSStefano Zampini } 50234a97f8cSStefano Zampini } 5035e8657edSStefano Zampini ierr = ISGlobalToLocalMappingApply(pcis->BtoNmap,IS_GTOLM_DROP,deluxe_ctx->n_simple,deluxe_ctx->idx_simple_B,&i,deluxe_ctx->idx_simple_B);CHKERRQ(ierr); 50434a97f8cSStefano Zampini if (i != deluxe_ctx->n_simple) { 50534a97f8cSStefano Zampini SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Error when mapping simple scaling dofs! %d != %d",i,deluxe_ctx->n_simple); 50634a97f8cSStefano Zampini } 5077821e9e7SStefano Zampini ierr = PetscBTDestroy(&bitmask);CHKERRQ(ierr); 50834a97f8cSStefano Zampini 50934a97f8cSStefano Zampini /* SetUp local schur complements on subsets TODO better reuse procedure */ 51034a97f8cSStefano Zampini if (!sub_schurs->n_subs) { 51134a97f8cSStefano Zampini Mat S_j; 51234a97f8cSStefano Zampini PetscBool free_used_adj; 51334a97f8cSStefano Zampini PetscInt *used_xadj,*used_adjncy; 51434a97f8cSStefano Zampini 51534a97f8cSStefano Zampini /* decide the adjacency to be used for determining internal problems for local schur on subsets */ 51634a97f8cSStefano Zampini free_used_adj = PETSC_FALSE; 51734a97f8cSStefano Zampini if (pcbddc->deluxe_layers == -1) { 51834a97f8cSStefano Zampini used_xadj = NULL; 51934a97f8cSStefano Zampini used_adjncy = NULL; 52034a97f8cSStefano Zampini } else { 52134a97f8cSStefano Zampini if ((pcbddc->deluxe_use_useradj && pcbddc->mat_graph->xadj) || !pcbddc->deluxe_compute_rowadj) { 52234a97f8cSStefano Zampini used_xadj = pcbddc->mat_graph->xadj; 52334a97f8cSStefano Zampini used_adjncy = pcbddc->mat_graph->adjncy; 52434a97f8cSStefano Zampini } else { 52534a97f8cSStefano Zampini Mat mat_adj; 52634a97f8cSStefano Zampini PetscBool flg_row=PETSC_TRUE; 52734a97f8cSStefano Zampini const PetscInt *xadj,*adjncy; 52834a97f8cSStefano Zampini PetscInt nvtxs; 52934a97f8cSStefano Zampini 53034a97f8cSStefano Zampini ierr = MatConvert(pcbddc->local_mat,MATMPIADJ,MAT_INITIAL_MATRIX,&mat_adj);CHKERRQ(ierr); 53134a97f8cSStefano Zampini ierr = MatGetRowIJ(mat_adj,0,PETSC_TRUE,PETSC_FALSE,&nvtxs,&xadj,&adjncy,&flg_row);CHKERRQ(ierr); 53234a97f8cSStefano Zampini if (!flg_row) { 53334a97f8cSStefano Zampini SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Error in MatGetRowIJ called in %s\n",__FUNCT__); 53434a97f8cSStefano Zampini } 53534a97f8cSStefano Zampini ierr = PetscMalloc2(nvtxs+1,&used_xadj,xadj[nvtxs],&used_adjncy);CHKERRQ(ierr); 53634a97f8cSStefano Zampini ierr = PetscMemcpy(used_xadj,xadj,(nvtxs+1)*sizeof(*xadj));CHKERRQ(ierr); 53734a97f8cSStefano Zampini ierr = PetscMemcpy(used_adjncy,adjncy,(xadj[nvtxs])*sizeof(*adjncy));CHKERRQ(ierr); 53834a97f8cSStefano Zampini ierr = MatRestoreRowIJ(mat_adj,0,PETSC_TRUE,PETSC_FALSE,&nvtxs,&xadj,&adjncy,&flg_row);CHKERRQ(ierr); 53934a97f8cSStefano Zampini if (!flg_row) { 54034a97f8cSStefano Zampini SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Error in MatRestoreRowIJ called in %s\n",__FUNCT__); 54134a97f8cSStefano Zampini } 54234a97f8cSStefano Zampini ierr = MatDestroy(&mat_adj);CHKERRQ(ierr); 54334a97f8cSStefano Zampini free_used_adj = PETSC_TRUE; 54434a97f8cSStefano Zampini } 54534a97f8cSStefano Zampini } 54634a97f8cSStefano Zampini 54734a97f8cSStefano Zampini /* Create Schur complement matrix */ 54834a97f8cSStefano Zampini ierr = MatCreateSchurComplement(pcis->A_II,pcis->A_II,pcis->A_IB,pcis->A_BI,pcis->A_BB,&S_j);CHKERRQ(ierr); 54934a97f8cSStefano Zampini ierr = MatSchurComplementSetKSP(S_j,pcbddc->ksp_D);CHKERRQ(ierr); 5507821e9e7SStefano Zampini 55134a97f8cSStefano Zampini /* setup Schur complements on subsets */ 5527821e9e7SStefano Zampini ierr = PCBDDCSubSchursSetUp(sub_schurs,S_j,pcis->is_I_local,pcis->is_B_local,n_all_cc,all_cc,used_xadj,used_adjncy,pcbddc->deluxe_layers);CHKERRQ(ierr); 55334a97f8cSStefano Zampini ierr = MatDestroy(&S_j);CHKERRQ(ierr); 55434a97f8cSStefano Zampini /* free adjacency */ 55534a97f8cSStefano Zampini if (free_used_adj) { 55634a97f8cSStefano Zampini ierr = PetscFree2(used_xadj,used_adjncy);CHKERRQ(ierr); 55734a97f8cSStefano Zampini } 55834a97f8cSStefano Zampini } 5597821e9e7SStefano Zampini for (i=0;i<n_all_cc;i++) { 5607821e9e7SStefano Zampini ierr = ISDestroy(&all_cc[i]);CHKERRQ(ierr); 5617821e9e7SStefano Zampini } 5627821e9e7SStefano Zampini ierr = PetscFree(all_cc);CHKERRQ(ierr); 56334a97f8cSStefano Zampini 56434a97f8cSStefano Zampini /* Number parallel problems */ 56534a97f8cSStefano Zampini auxglobal_parallel = 0; 56634a97f8cSStefano Zampini ierr = PCBDDCSubsetNumbering(PetscObjectComm((PetscObject)pc),graph->l2gmap,n_local_parallel_problems,auxlocal_parallel,PETSC_NULL,&n_parallel_problems,&auxglobal_parallel);CHKERRQ(ierr); 56734a97f8cSStefano Zampini if (pcbddc->dbg_flag) { 56834a97f8cSStefano Zampini ierr = PetscViewerASCIIPrintf(pcbddc->dbg_viewer,"Deluxe global number of parallel subproblems: %d\n",n_parallel_problems); 56934a97f8cSStefano Zampini } 57034a97f8cSStefano Zampini 57134a97f8cSStefano Zampini /* Compute data structures to solve parallel problems */ 57234a97f8cSStefano Zampini ierr = PCBDDCScalingSetUp_Deluxe_Par(pc,n_local_parallel_problems,n_parallel_problems,auxglobal_parallel,index_parallel);CHKERRQ(ierr); 57334a97f8cSStefano Zampini ierr = PetscFree(auxglobal_parallel);CHKERRQ(ierr); 57434a97f8cSStefano Zampini 57534a97f8cSStefano Zampini 57634a97f8cSStefano Zampini /* Number sequential problems */ 57734a97f8cSStefano Zampini auxglobal_sequential = 0; 57834a97f8cSStefano Zampini ierr = PCBDDCSubsetNumbering(PetscObjectComm((PetscObject)pc),graph->l2gmap,n_local_sequential_problems,auxlocal_sequential,PETSC_NULL,&n_sequential_problems,&auxglobal_sequential);CHKERRQ(ierr); 57934a97f8cSStefano Zampini if (pcbddc->dbg_flag) { 58034a97f8cSStefano Zampini ierr = PetscViewerASCIIPrintf(pcbddc->dbg_viewer,"Deluxe global number of sequential subproblems: %d\n",n_sequential_problems); 58134a97f8cSStefano Zampini } 58234a97f8cSStefano Zampini 58334a97f8cSStefano Zampini /* Compute data structures to solve sequential problems */ 58434a97f8cSStefano Zampini ierr = PCBDDCScalingSetUp_Deluxe_Seq(pc,n_local_sequential_problems,n_sequential_problems,auxglobal_sequential,index_sequential);CHKERRQ(ierr); 58534a97f8cSStefano Zampini ierr = PetscFree(auxglobal_sequential);CHKERRQ(ierr); 58634a97f8cSStefano Zampini 58734a97f8cSStefano Zampini /* free workspace */ 58834a97f8cSStefano Zampini ierr = PetscFree5(auxmapping,auxlocal_sequential,auxlocal_parallel,index_sequential,index_parallel);CHKERRQ(ierr); 589b1b3d7a2SStefano Zampini #else 59034a97f8cSStefano Zampini 591b1b3d7a2SStefano Zampini /* decide the adjacency to be used for determining internal problems for local schur on subsets */ 592b1b3d7a2SStefano Zampini free_used_adj = PETSC_FALSE; 593b1b3d7a2SStefano Zampini if (pcbddc->deluxe_layers == -1) { 594b1b3d7a2SStefano Zampini used_xadj = NULL; 595b1b3d7a2SStefano Zampini used_adjncy = NULL; 596b1b3d7a2SStefano Zampini } else { 597b1b3d7a2SStefano Zampini if ((pcbddc->deluxe_use_useradj && pcbddc->mat_graph->xadj) || !pcbddc->deluxe_compute_rowadj) { 598b1b3d7a2SStefano Zampini used_xadj = pcbddc->mat_graph->xadj; 599b1b3d7a2SStefano Zampini used_adjncy = pcbddc->mat_graph->adjncy; 600b1b3d7a2SStefano Zampini } else { 601b1b3d7a2SStefano Zampini Mat mat_adj; 602b1b3d7a2SStefano Zampini PetscBool flg_row=PETSC_TRUE; 603b1b3d7a2SStefano Zampini const PetscInt *xadj,*adjncy; 604b1b3d7a2SStefano Zampini PetscInt nvtxs; 605b1b3d7a2SStefano Zampini 606b1b3d7a2SStefano Zampini ierr = MatConvert(pcbddc->local_mat,MATMPIADJ,MAT_INITIAL_MATRIX,&mat_adj);CHKERRQ(ierr); 607b1b3d7a2SStefano Zampini ierr = MatGetRowIJ(mat_adj,0,PETSC_TRUE,PETSC_FALSE,&nvtxs,&xadj,&adjncy,&flg_row);CHKERRQ(ierr); 608b1b3d7a2SStefano Zampini if (!flg_row) { 609b1b3d7a2SStefano Zampini SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Error in MatGetRowIJ called in %s\n",__FUNCT__); 610b1b3d7a2SStefano Zampini } 611b1b3d7a2SStefano Zampini ierr = PetscMalloc2(nvtxs+1,&used_xadj,xadj[nvtxs],&used_adjncy);CHKERRQ(ierr); 612b1b3d7a2SStefano Zampini ierr = PetscMemcpy(used_xadj,xadj,(nvtxs+1)*sizeof(*xadj));CHKERRQ(ierr); 613b1b3d7a2SStefano Zampini ierr = PetscMemcpy(used_adjncy,adjncy,(xadj[nvtxs])*sizeof(*adjncy));CHKERRQ(ierr); 614b1b3d7a2SStefano Zampini ierr = MatRestoreRowIJ(mat_adj,0,PETSC_TRUE,PETSC_FALSE,&nvtxs,&xadj,&adjncy,&flg_row);CHKERRQ(ierr); 615b1b3d7a2SStefano Zampini if (!flg_row) { 616b1b3d7a2SStefano Zampini SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Error in MatRestoreRowIJ called in %s\n",__FUNCT__); 617b1b3d7a2SStefano Zampini } 618b1b3d7a2SStefano Zampini ierr = MatDestroy(&mat_adj);CHKERRQ(ierr); 619b1b3d7a2SStefano Zampini free_used_adj = PETSC_TRUE; 620b1b3d7a2SStefano Zampini } 621b1b3d7a2SStefano Zampini } 622b1b3d7a2SStefano Zampini 623b1b3d7a2SStefano Zampini 624b1b3d7a2SStefano Zampini /* Create Schur complement matrix */ 625b1b3d7a2SStefano Zampini ierr = MatCreateSchurComplement(pcis->A_II,pcis->A_II,pcis->A_IB,pcis->A_BI,pcis->A_BB,&S_j);CHKERRQ(ierr); 626b1b3d7a2SStefano Zampini ierr = MatSchurComplementSetKSP(S_j,pcbddc->ksp_D);CHKERRQ(ierr); 627b1b3d7a2SStefano Zampini 6281e9c79c2SStefano Zampini /* sub_schurs init */ /* TODO reuse adaptive one if valid (i.e. pcbddc->local_mat == matis->A and same graph info (HOW?) ) */ 6291e9c79c2SStefano Zampini ierr = PCBDDCSubSchursInit(sub_schurs,pcbddc->local_mat,S_j,pcis->is_I_local,pcis->is_B_local,graph,pcbddc->deluxe_threshold);CHKERRQ(ierr); 630b1b3d7a2SStefano Zampini ierr = MatDestroy(&S_j);CHKERRQ(ierr); 631b1b3d7a2SStefano Zampini ierr = PCBDDCSubSchursSetUpNew(sub_schurs,used_xadj,used_adjncy,pcbddc->deluxe_layers);CHKERRQ(ierr); 632b1b3d7a2SStefano Zampini 633b1b3d7a2SStefano Zampini /* Compute data structures to solve parallel problems */ 634b1b3d7a2SStefano Zampini ierr = PCBDDCScalingSetUp_Deluxe_Par(pc,sub_schurs->n_subs_par,sub_schurs->n_subs_par_g, 635b1b3d7a2SStefano Zampini sub_schurs->auxglobal_parallel, 636b1b3d7a2SStefano Zampini sub_schurs->index_parallel);CHKERRQ(ierr); 637b1b3d7a2SStefano Zampini /* Compute data structures to solve sequential problems */ 638b1b3d7a2SStefano Zampini ierr = PCBDDCScalingSetUp_Deluxe_Seq(pc,sub_schurs->n_subs_seq,sub_schurs->n_subs_seq_g, 639b1b3d7a2SStefano Zampini sub_schurs->auxglobal_sequential, 640b1b3d7a2SStefano Zampini sub_schurs->index_sequential);CHKERRQ(ierr); 641b1b3d7a2SStefano Zampini /* free adjacency */ 642b1b3d7a2SStefano Zampini if (free_used_adj) { 643b1b3d7a2SStefano Zampini ierr = PetscFree2(used_xadj,used_adjncy);CHKERRQ(ierr); 644b1b3d7a2SStefano Zampini } 645b1b3d7a2SStefano Zampini 646b1b3d7a2SStefano Zampini /* diagonal scaling on interface dofs not contained in cc */ 647b1b3d7a2SStefano Zampini ierr = PetscBTCreate(pcis->n,&bitmask);CHKERRQ(ierr); 648b1b3d7a2SStefano Zampini ierr = ISGetIndices(pcis->is_I_local,&idxs);CHKERRQ(ierr); 649b1b3d7a2SStefano Zampini for (i=0;i<pcis->n-pcis->n_B;i++) { 650b1b3d7a2SStefano Zampini ierr = PetscBTSet(bitmask,idxs[i]);CHKERRQ(ierr); 651b1b3d7a2SStefano Zampini } 652b1b3d7a2SStefano Zampini ierr = ISRestoreIndices(pcis->is_I_local,&idxs);CHKERRQ(ierr); 653b1b3d7a2SStefano Zampini 654b1b3d7a2SStefano Zampini for (i=0;i<sub_schurs->n_subs;i++) { 655b1b3d7a2SStefano Zampini PetscInt subset_size,j; 656b1b3d7a2SStefano Zampini 657b1b3d7a2SStefano Zampini ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr); 658b1b3d7a2SStefano Zampini ierr = ISGetIndices(sub_schurs->is_subs[i],&idxs);CHKERRQ(ierr); 659b1b3d7a2SStefano Zampini for (j=0;j<subset_size;j++) { 660b1b3d7a2SStefano Zampini ierr = PetscBTSet(bitmask,idxs[j]);CHKERRQ(ierr); 661b1b3d7a2SStefano Zampini } 662b1b3d7a2SStefano Zampini ierr = ISRestoreIndices(sub_schurs->is_subs[i],&idxs);CHKERRQ(ierr); 663b1b3d7a2SStefano Zampini } 664b1b3d7a2SStefano Zampini 665b1b3d7a2SStefano Zampini deluxe_ctx->n_simple = 0; 666b1b3d7a2SStefano Zampini for (i=0;i<pcis->n;i++) { 667b1b3d7a2SStefano Zampini if (!PetscBTLookup(bitmask,i)) { 668b1b3d7a2SStefano Zampini deluxe_ctx->n_simple++; 669b1b3d7a2SStefano Zampini } 670b1b3d7a2SStefano Zampini } 671b1b3d7a2SStefano Zampini ierr = PetscMalloc1(deluxe_ctx->n_simple,&deluxe_ctx->idx_simple_B);CHKERRQ(ierr); 672b1b3d7a2SStefano Zampini deluxe_ctx->n_simple = 0; 673b1b3d7a2SStefano Zampini for (i=0;i<pcis->n;i++) { 674b1b3d7a2SStefano Zampini if (!PetscBTLookup(bitmask,i)) { 675b1b3d7a2SStefano Zampini deluxe_ctx->idx_simple_B[deluxe_ctx->n_simple++] = i; 676b1b3d7a2SStefano Zampini } 677b1b3d7a2SStefano Zampini } 678b1b3d7a2SStefano Zampini ierr = ISGlobalToLocalMappingApply(pcis->BtoNmap,IS_GTOLM_DROP,deluxe_ctx->n_simple,deluxe_ctx->idx_simple_B,&i,deluxe_ctx->idx_simple_B);CHKERRQ(ierr); 679b1b3d7a2SStefano Zampini if (i != deluxe_ctx->n_simple) { 680b1b3d7a2SStefano Zampini SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Error when mapping simple scaling dofs! %d != %d",i,deluxe_ctx->n_simple); 681b1b3d7a2SStefano Zampini } 682b1b3d7a2SStefano Zampini ierr = PetscBTDestroy(&bitmask);CHKERRQ(ierr); 683b1b3d7a2SStefano Zampini 684b1b3d7a2SStefano Zampini #endif 68534a97f8cSStefano Zampini /* free graph struct */ 68634a97f8cSStefano Zampini if (pcbddc->deluxe_rebuild) { 68734a97f8cSStefano Zampini ierr = PCBDDCGraphDestroy(&graph);CHKERRQ(ierr); 68834a97f8cSStefano Zampini } 68934a97f8cSStefano Zampini PetscFunctionReturn(0); 69034a97f8cSStefano Zampini } 69134a97f8cSStefano Zampini 69234a97f8cSStefano Zampini #undef __FUNCT__ 69334a97f8cSStefano Zampini #define __FUNCT__ "PCBDDCScalingSetUp_Deluxe_Par" 69434a97f8cSStefano Zampini static PetscErrorCode PCBDDCScalingSetUp_Deluxe_Par(PC pc, PetscInt n_local_parallel_problems,PetscInt n_parallel_problems,PetscInt global_parallel[],PetscInt index_parallel[]) 69534a97f8cSStefano Zampini { 69634a97f8cSStefano Zampini PC_BDDC *pcbddc=(PC_BDDC*)pc->data; 69734a97f8cSStefano Zampini PCBDDCDeluxeScaling deluxe_ctx=pcbddc->deluxe_ctx; 69834a97f8cSStefano Zampini /* coloring */ 69934a97f8cSStefano Zampini Mat parallel_problems; 70034a97f8cSStefano Zampini MatColoring coloring_obj; 70134a97f8cSStefano Zampini ISColoring coloring_parallel_problems; 70234a97f8cSStefano Zampini IS *par_is_colors,*is_colors; 70334a97f8cSStefano Zampini /* working stuff */ 70434a97f8cSStefano Zampini PetscInt i,j; 70534a97f8cSStefano Zampini PetscErrorCode ierr; 70634a97f8cSStefano Zampini 70734a97f8cSStefano Zampini PetscFunctionBegin; 70834a97f8cSStefano Zampini if (!n_parallel_problems) { 70934a97f8cSStefano Zampini PetscFunctionReturn(0); 71034a97f8cSStefano Zampini } 71134a97f8cSStefano Zampini /* Color parallel subproblems */ 71234a97f8cSStefano Zampini ierr = MatCreate(PetscObjectComm((PetscObject)pc),¶llel_problems);CHKERRQ(ierr); 71334a97f8cSStefano Zampini ierr = MatSetSizes(parallel_problems,PETSC_DECIDE,PETSC_DECIDE,n_parallel_problems,n_parallel_problems);CHKERRQ(ierr); 71434a97f8cSStefano Zampini ierr = MatSetType(parallel_problems,MATAIJ);CHKERRQ(ierr); 71534a97f8cSStefano Zampini ierr = MatSetUp(parallel_problems);CHKERRQ(ierr); 71634a97f8cSStefano Zampini ierr = MatSetOption(parallel_problems,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_FALSE);CHKERRQ(ierr); 71734a97f8cSStefano Zampini ierr = MatSetOption(parallel_problems,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_FALSE);CHKERRQ(ierr); 71834a97f8cSStefano Zampini for (i=0;i<n_local_parallel_problems;i++) { 71934a97f8cSStefano Zampini PetscInt row = global_parallel[i]; 72034a97f8cSStefano Zampini for (j=0;j<n_local_parallel_problems;j++) { 72134a97f8cSStefano Zampini PetscInt col = global_parallel[j]; 72234a97f8cSStefano Zampini if (row != col) { 72334a97f8cSStefano Zampini ierr = MatSetValue(parallel_problems,row,col,1.0,INSERT_VALUES);CHKERRQ(ierr); 72434a97f8cSStefano Zampini } 72534a97f8cSStefano Zampini } 72634a97f8cSStefano Zampini } 72734a97f8cSStefano Zampini ierr = MatAssemblyBegin(parallel_problems,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 72834a97f8cSStefano Zampini ierr = MatAssemblyEnd(parallel_problems,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 72934a97f8cSStefano Zampini if (pcbddc->dbg_flag > 1) { 73034a97f8cSStefano Zampini ierr = PetscViewerFlush(pcbddc->dbg_viewer);CHKERRQ(ierr); 73134a97f8cSStefano Zampini ierr = PetscViewerASCIIPrintf(pcbddc->dbg_viewer,"Adj matrix for deluxe parallel problems\n");CHKERRQ(ierr); 73234a97f8cSStefano Zampini ierr = MatView(parallel_problems,pcbddc->dbg_viewer);CHKERRQ(ierr); 73334a97f8cSStefano Zampini } 73434a97f8cSStefano Zampini ierr = MatColoringCreate(parallel_problems,&coloring_obj);CHKERRQ(ierr); 73534a97f8cSStefano Zampini ierr = MatColoringSetDistance(coloring_obj,1);CHKERRQ(ierr); 73634a97f8cSStefano Zampini ierr = MatColoringSetType(coloring_obj,MATCOLORINGJP);CHKERRQ(ierr); 73734a97f8cSStefano Zampini ierr = MatColoringApply(coloring_obj,&coloring_parallel_problems);CHKERRQ(ierr); 73834a97f8cSStefano Zampini ierr = ISColoringGetIS(coloring_parallel_problems,&deluxe_ctx->par_colors,&par_is_colors);CHKERRQ(ierr); 73934a97f8cSStefano Zampini if (pcbddc->dbg_flag) { 74034a97f8cSStefano Zampini ierr = PetscViewerFlush(pcbddc->dbg_viewer);CHKERRQ(ierr); 74134a97f8cSStefano Zampini ierr = PetscViewerASCIIPrintf(pcbddc->dbg_viewer,"Number of colors %d for parallel part of deluxe\n",deluxe_ctx->par_colors);CHKERRQ(ierr); 74234a97f8cSStefano Zampini } 74334a97f8cSStefano Zampini 74434a97f8cSStefano Zampini /* all procs should know the color distribution */ 74534a97f8cSStefano Zampini ierr = PetscMalloc1(deluxe_ctx->par_colors,&is_colors);CHKERRQ(ierr); 74634a97f8cSStefano Zampini for (i=0;i<deluxe_ctx->par_colors;i++) { 74734a97f8cSStefano Zampini if (pcbddc->dbg_flag) { 74834a97f8cSStefano Zampini ierr = PetscViewerASCIIPrintf(pcbddc->dbg_viewer,"Global problem indexes for color %d\n",i);CHKERRQ(ierr); 74934a97f8cSStefano Zampini ierr = ISView(par_is_colors[i],pcbddc->dbg_viewer);CHKERRQ(ierr); 75034a97f8cSStefano Zampini ierr = PetscViewerFlush(pcbddc->dbg_viewer);CHKERRQ(ierr); 75134a97f8cSStefano Zampini } 75234a97f8cSStefano Zampini ierr = ISAllGather(par_is_colors[i],&is_colors[i]);CHKERRQ(ierr); 75334a97f8cSStefano Zampini } 75434a97f8cSStefano Zampini 75534a97f8cSStefano Zampini /* free unneeded objects */ 75634a97f8cSStefano Zampini ierr = ISColoringRestoreIS(coloring_parallel_problems,&par_is_colors);CHKERRQ(ierr); 75734a97f8cSStefano Zampini ierr = ISColoringDestroy(&coloring_parallel_problems);CHKERRQ(ierr); 75834a97f8cSStefano Zampini ierr = MatColoringDestroy(&coloring_obj);CHKERRQ(ierr); 75934a97f8cSStefano Zampini ierr = MatDestroy(¶llel_problems);CHKERRQ(ierr); 76034a97f8cSStefano Zampini 76134a97f8cSStefano Zampini /* allocate deluxe arrays for parallel problems */ 76234a97f8cSStefano Zampini ierr = PetscMalloc5(deluxe_ctx->par_colors,&deluxe_ctx->par_ksp, 76334a97f8cSStefano Zampini deluxe_ctx->par_colors,&deluxe_ctx->par_scctx_s, 76434a97f8cSStefano Zampini deluxe_ctx->par_colors,&deluxe_ctx->par_scctx_p, 76534a97f8cSStefano Zampini deluxe_ctx->par_colors,&deluxe_ctx->par_vec, 76634a97f8cSStefano Zampini deluxe_ctx->par_colors,&deluxe_ctx->par_col2sub);CHKERRQ(ierr); 76734a97f8cSStefano Zampini 76834a97f8cSStefano Zampini /* cycle on colors */ 76934a97f8cSStefano Zampini for (i=0;i<deluxe_ctx->par_colors;i++) { 77034a97f8cSStefano Zampini PetscSubcomm par_subcomm; 77134a97f8cSStefano Zampini const PetscInt* idxs_subproblems; 77234a97f8cSStefano Zampini PetscInt color_size; 77334a97f8cSStefano Zampini PetscMPIInt rank,active_color; 77434a97f8cSStefano Zampini 77534a97f8cSStefano Zampini /* get local index of i-th parallel colored problem */ 77634a97f8cSStefano Zampini ierr = ISGetLocalSize(is_colors[i],&color_size);CHKERRQ(ierr); 77734a97f8cSStefano Zampini ierr = ISGetIndices(is_colors[i],&idxs_subproblems);CHKERRQ(ierr); 77834a97f8cSStefano Zampini /* split comm for computing parallel problems for this color */ 77934a97f8cSStefano Zampini /* Processes not partecipating at this stage will have color = color_size */ 78034a97f8cSStefano Zampini /* because PetscCommDuplicate does not handle MPI_COMM_NULL */ 78134a97f8cSStefano Zampini active_color = color_size; 78234a97f8cSStefano Zampini deluxe_ctx->par_col2sub[i] = -1; 78334a97f8cSStefano Zampini for (j=0;j<n_local_parallel_problems;j++) { 78434a97f8cSStefano Zampini PetscInt local_idx; 78534a97f8cSStefano Zampini ierr = PetscFindInt(global_parallel[j],color_size,idxs_subproblems,&local_idx);CHKERRQ(ierr); 78634a97f8cSStefano Zampini if (local_idx > -1) { 78734a97f8cSStefano Zampini ierr = PetscMPIIntCast(local_idx,&active_color);CHKERRQ(ierr); 78834a97f8cSStefano Zampini deluxe_ctx->par_col2sub[i] = index_parallel[j]; 78934a97f8cSStefano Zampini break; 79034a97f8cSStefano Zampini } 79134a97f8cSStefano Zampini } 79234a97f8cSStefano Zampini ierr = ISRestoreIndices(is_colors[i],&idxs_subproblems);CHKERRQ(ierr); 79334a97f8cSStefano Zampini ierr = PetscSubcommCreate(PetscObjectComm((PetscObject)pc),&par_subcomm);CHKERRQ(ierr); 79434a97f8cSStefano Zampini ierr = PetscSubcommSetNumber(par_subcomm,color_size+1);CHKERRQ(ierr); 79534a97f8cSStefano Zampini ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)pc),&rank);CHKERRQ(ierr); 79634a97f8cSStefano Zampini ierr = PetscSubcommSetTypeGeneral(par_subcomm,active_color,rank);CHKERRQ(ierr); 79734a97f8cSStefano Zampini /* print debug info */ 79834a97f8cSStefano Zampini if (pcbddc->dbg_flag) { 79934a97f8cSStefano Zampini PetscMPIInt crank,csize; 80034a97f8cSStefano Zampini ierr = MPI_Comm_rank(par_subcomm->comm,&crank);CHKERRQ(ierr); 80134a97f8cSStefano Zampini ierr = MPI_Comm_size(par_subcomm->comm,&csize);CHKERRQ(ierr); 80234a97f8cSStefano Zampini ierr = PetscViewerASCIIPrintf(pcbddc->dbg_viewer,"Color %d: size %d, details follows.\n",i,color_size);CHKERRQ(ierr); 80334a97f8cSStefano Zampini ierr = PetscViewerFlush(pcbddc->dbg_viewer);CHKERRQ(ierr); 80434a97f8cSStefano Zampini ierr = PetscViewerASCIISynchronizedAllow(pcbddc->dbg_viewer,PETSC_TRUE);CHKERRQ(ierr); 80534a97f8cSStefano Zampini ierr = PetscViewerASCIISynchronizedPrintf(pcbddc->dbg_viewer," Subdomain %d: color in subcomm %d (rank %d out of %d) (lidx %d)\n",PetscGlobalRank,par_subcomm->color,crank,csize,deluxe_ctx->par_col2sub[i]);CHKERRQ(ierr); 80634a97f8cSStefano Zampini ierr = PetscViewerFlush(pcbddc->dbg_viewer);CHKERRQ(ierr); 80734a97f8cSStefano Zampini } 80834a97f8cSStefano Zampini 80934a97f8cSStefano Zampini if (deluxe_ctx->par_col2sub[i] >= 0) { 8105e8657edSStefano Zampini PC pctemp; 8115e8657edSStefano Zampini PC_IS *pcis=(PC_IS*)pc->data; 81234a97f8cSStefano Zampini Mat color_mat,color_mat_is,temp_mat; 81334a97f8cSStefano Zampini ISLocalToGlobalMapping WtoNmap,l2gmap_subset; 81434a97f8cSStefano Zampini IS is_local_numbering,isB_local,isW_local,isW; 815b1b3d7a2SStefano Zampini PCBDDCSubSchurs sub_schurs = pcbddc->sub_schurs[1]; 81634a97f8cSStefano Zampini PetscInt subidx,n_local_dofs,n_global_dofs; 81734a97f8cSStefano Zampini PetscInt *global_numbering,*local_numbering; 81834a97f8cSStefano Zampini char ksp_prefix[256]; 81934a97f8cSStefano Zampini size_t len; 82034a97f8cSStefano Zampini 82134a97f8cSStefano Zampini /* Local index for schur complement on subset */ 82234a97f8cSStefano Zampini subidx = deluxe_ctx->par_col2sub[i]; 82334a97f8cSStefano Zampini 82434a97f8cSStefano Zampini /* Parallel numbering for dofs in colored subset */ 82534a97f8cSStefano Zampini ierr = ISSum(sub_schurs->is_AEj_I[subidx],sub_schurs->is_AEj_B[subidx],&is_local_numbering);CHKERRQ(ierr); 82634a97f8cSStefano Zampini ierr = ISGetLocalSize(is_local_numbering,&n_local_dofs);CHKERRQ(ierr); 82734a97f8cSStefano Zampini ierr = ISGetIndices(is_local_numbering,(const PetscInt **)&local_numbering);CHKERRQ(ierr); 82834a97f8cSStefano Zampini ierr = PCBDDCSubsetNumbering(par_subcomm->comm,pcbddc->mat_graph->l2gmap,n_local_dofs,local_numbering,PETSC_NULL,&n_global_dofs,&global_numbering);CHKERRQ(ierr); 82934a97f8cSStefano Zampini ierr = ISRestoreIndices(is_local_numbering,(const PetscInt **)&local_numbering);CHKERRQ(ierr); 83034a97f8cSStefano Zampini 83134a97f8cSStefano Zampini /* L2Gmap from relevant dofs to local dofs */ 83234a97f8cSStefano Zampini ierr = ISLocalToGlobalMappingCreateIS(is_local_numbering,&WtoNmap);CHKERRQ(ierr); 83334a97f8cSStefano Zampini 83434a97f8cSStefano Zampini /* L2Gmap from local to global dofs */ 83534a97f8cSStefano Zampini ierr = ISLocalToGlobalMappingCreate(par_subcomm->comm,1,n_local_dofs,global_numbering,PETSC_COPY_VALUES,&l2gmap_subset);CHKERRQ(ierr); 83634a97f8cSStefano Zampini 83734a97f8cSStefano Zampini /* compute parallel matrix (extended dirichlet problem on subset) */ 83834a97f8cSStefano Zampini ierr = MatCreateIS(par_subcomm->comm,1,PETSC_DECIDE,PETSC_DECIDE,n_global_dofs,n_global_dofs,l2gmap_subset,&color_mat_is);CHKERRQ(ierr); 83934a97f8cSStefano Zampini ierr = MatGetSubMatrix(pcbddc->local_mat,is_local_numbering,is_local_numbering,MAT_INITIAL_MATRIX,&temp_mat);CHKERRQ(ierr); 84034a97f8cSStefano Zampini ierr = MatISSetLocalMat(color_mat_is,temp_mat);CHKERRQ(ierr); 84134a97f8cSStefano Zampini ierr = MatDestroy(&temp_mat);CHKERRQ(ierr); 84234a97f8cSStefano Zampini ierr = MatISGetMPIXAIJ(color_mat_is,MAT_INITIAL_MATRIX,&color_mat);CHKERRQ(ierr); 84334a97f8cSStefano Zampini ierr = MatDestroy(&color_mat_is);CHKERRQ(ierr); 84434a97f8cSStefano Zampini 84534a97f8cSStefano Zampini /* work vector for (parallel) extended dirichlet problem */ 8468a26ef87SStefano Zampini ierr = MatCreateVecs(color_mat,&deluxe_ctx->par_vec[i],NULL);CHKERRQ(ierr); 84734a97f8cSStefano Zampini 84834a97f8cSStefano Zampini /* compute scatters */ 84934a97f8cSStefano Zampini /* deluxe_ctx->par_scctx_p[i] extension from local subset to extended dirichlet problem 85034a97f8cSStefano Zampini deluxe_ctx->par_scctx_s[i] restriction from local boundary to subset -> simple copy of selected values */ 8515e8657edSStefano Zampini ierr = ISGlobalToLocalMappingApplyIS(pcis->BtoNmap,IS_GTOLM_DROP,sub_schurs->is_AEj_B[subidx],&isB_local);CHKERRQ(ierr); 85234a97f8cSStefano Zampini ierr = VecScatterCreate(pcbddc->work_scaling,isB_local,sub_schurs->work1[subidx],NULL,&deluxe_ctx->par_scctx_s[i]);CHKERRQ(ierr); 85334a97f8cSStefano Zampini ierr = ISGlobalToLocalMappingApplyIS(WtoNmap,IS_GTOLM_DROP,sub_schurs->is_AEj_B[subidx],&isW_local);CHKERRQ(ierr); 85434a97f8cSStefano Zampini ierr = ISLocalToGlobalMappingApplyIS(l2gmap_subset,isW_local,&isW);CHKERRQ(ierr); 85534a97f8cSStefano Zampini ierr = VecScatterCreate(sub_schurs->work1[subidx],NULL,deluxe_ctx->par_vec[i],isW,&deluxe_ctx->par_scctx_p[i]);CHKERRQ(ierr); 85634a97f8cSStefano Zampini 85734a97f8cSStefano Zampini /* free objects no longer neeeded */ 85834a97f8cSStefano Zampini ierr = ISDestroy(&isW);CHKERRQ(ierr); 85934a97f8cSStefano Zampini ierr = ISDestroy(&isW_local);CHKERRQ(ierr); 86034a97f8cSStefano Zampini ierr = ISDestroy(&isB_local);CHKERRQ(ierr); 86134a97f8cSStefano Zampini ierr = ISLocalToGlobalMappingDestroy(&WtoNmap);CHKERRQ(ierr); 86234a97f8cSStefano Zampini ierr = ISLocalToGlobalMappingDestroy(&l2gmap_subset);CHKERRQ(ierr); 86334a97f8cSStefano Zampini ierr = ISDestroy(&is_local_numbering);CHKERRQ(ierr); 86434a97f8cSStefano Zampini ierr = PetscFree(global_numbering);CHKERRQ(ierr); 86534a97f8cSStefano Zampini 86634a97f8cSStefano Zampini /* KSP for extended dirichlet problem */ 86734a97f8cSStefano Zampini ierr = KSPCreate(par_subcomm->comm,&deluxe_ctx->par_ksp[i]);CHKERRQ(ierr); 86834a97f8cSStefano Zampini ierr = KSPSetOperators(deluxe_ctx->par_ksp[i],color_mat,color_mat);CHKERRQ(ierr); 86934a97f8cSStefano Zampini ierr = KSPSetTolerances(deluxe_ctx->par_ksp[i],1.e-12,1.e-12,1.e10,10000);CHKERRQ(ierr); 87034a97f8cSStefano Zampini ierr = KSPSetType(deluxe_ctx->par_ksp[i],KSPPREONLY);CHKERRQ(ierr); 8715e8657edSStefano Zampini ierr = KSPGetPC(deluxe_ctx->par_ksp[i],&pctemp);CHKERRQ(ierr); 8725e8657edSStefano Zampini ierr = PCSetType(pctemp,PCREDUNDANT);CHKERRQ(ierr); 87334a97f8cSStefano Zampini ierr = PetscStrlen(((PetscObject)(pcbddc->ksp_D))->prefix,&len);CHKERRQ(ierr); 8748856534dSStefano Zampini len -= 10; /* remove "dirichlet_" */ 8758856534dSStefano Zampini ierr = PetscStrncpy(ksp_prefix,((PetscObject)(pcbddc->ksp_D))->prefix,len+1);CHKERRQ(ierr); /* PetscStrncpy puts a terminating char at the end */ 87634a97f8cSStefano Zampini ierr = PetscStrcat(ksp_prefix,"deluxe_par_");CHKERRQ(ierr); 87734a97f8cSStefano Zampini ierr = KSPSetOptionsPrefix(deluxe_ctx->par_ksp[i],ksp_prefix);CHKERRQ(ierr); 87834a97f8cSStefano Zampini ierr = KSPSetFromOptions(deluxe_ctx->par_ksp[i]);CHKERRQ(ierr); 87934a97f8cSStefano Zampini ierr = KSPSetUp(deluxe_ctx->par_ksp[i]);CHKERRQ(ierr); 88034a97f8cSStefano Zampini ierr = MatDestroy(&color_mat);CHKERRQ(ierr); 88134a97f8cSStefano Zampini } else { /* not partecipating in color */ 88234a97f8cSStefano Zampini deluxe_ctx->par_ksp[i] = 0; 88334a97f8cSStefano Zampini deluxe_ctx->par_vec[i] = 0; 88434a97f8cSStefano Zampini deluxe_ctx->par_scctx_p[i] = 0; 88534a97f8cSStefano Zampini deluxe_ctx->par_scctx_s[i] = 0; 88634a97f8cSStefano Zampini } 88734a97f8cSStefano Zampini ierr = PetscSubcommDestroy(&par_subcomm);CHKERRQ(ierr); 88834a97f8cSStefano Zampini } 88934a97f8cSStefano Zampini for (i=0;i<deluxe_ctx->par_colors;i++) { 89034a97f8cSStefano Zampini ierr = ISDestroy(&is_colors[i]);CHKERRQ(ierr); 89134a97f8cSStefano Zampini } 89234a97f8cSStefano Zampini ierr = PetscFree(is_colors);CHKERRQ(ierr); 89334a97f8cSStefano Zampini 89434a97f8cSStefano Zampini if (pcbddc->dbg_flag) { 89534a97f8cSStefano Zampini Vec test_vec; 89634a97f8cSStefano Zampini PetscReal error; 897b1b3d7a2SStefano Zampini PCBDDCSubSchurs sub_schurs = pcbddc->sub_schurs[1]; 89834a97f8cSStefano Zampini /* test partition of unity of coloured schur complements */ 89934a97f8cSStefano Zampini for (i=0;i<deluxe_ctx->par_colors;i++) { 90034a97f8cSStefano Zampini PetscInt subidx = deluxe_ctx->par_col2sub[i]; 90134a97f8cSStefano Zampini PetscBool error_found = PETSC_FALSE; 90234a97f8cSStefano Zampini ierr = PetscViewerASCIISynchronizedAllow(pcbddc->dbg_viewer,PETSC_TRUE);CHKERRQ(ierr); 90334a97f8cSStefano Zampini 90434a97f8cSStefano Zampini if (deluxe_ctx->par_ksp[i]) { 90534a97f8cSStefano Zampini /* create random test vec being zero on internal nodes of the extende dirichlet problem */ 90634a97f8cSStefano Zampini ierr = VecDuplicate(deluxe_ctx->par_vec[i],&test_vec);CHKERRQ(ierr); 90734a97f8cSStefano Zampini ierr = VecSetRandom(sub_schurs->work1[subidx],PETSC_NULL);CHKERRQ(ierr); 90834a97f8cSStefano Zampini ierr = VecSet(test_vec,0.0);CHKERRQ(ierr); 90934a97f8cSStefano Zampini ierr = VecScatterBegin(deluxe_ctx->par_scctx_p[i],sub_schurs->work1[subidx],test_vec,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 91034a97f8cSStefano Zampini ierr = VecScatterEnd(deluxe_ctx->par_scctx_p[i],sub_schurs->work1[subidx],test_vec,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 91134a97f8cSStefano Zampini /* w_j */ 91234a97f8cSStefano Zampini ierr = VecScatterBegin(deluxe_ctx->par_scctx_p[i],test_vec,sub_schurs->work1[subidx],INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 91334a97f8cSStefano Zampini ierr = VecScatterEnd(deluxe_ctx->par_scctx_p[i],test_vec,sub_schurs->work1[subidx],INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 91434a97f8cSStefano Zampini /* S_j*w_j */ 91534a97f8cSStefano Zampini ierr = MatMult(sub_schurs->S_Ej[subidx],sub_schurs->work1[subidx],sub_schurs->work2[subidx]);CHKERRQ(ierr); 91634a97f8cSStefano Zampini /* \sum_j S_j*w_j */ 91734a97f8cSStefano Zampini ierr = VecSet(deluxe_ctx->par_vec[i],0.0);CHKERRQ(ierr); 91834a97f8cSStefano Zampini ierr = VecScatterBegin(deluxe_ctx->par_scctx_p[i],sub_schurs->work2[subidx],deluxe_ctx->par_vec[i],ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 91934a97f8cSStefano Zampini ierr = VecScatterEnd(deluxe_ctx->par_scctx_p[i],sub_schurs->work2[subidx],deluxe_ctx->par_vec[i],ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 92034a97f8cSStefano Zampini /* (\sum_j S_j)^(-1)(\sum_j S_j*w_j) */ 92134a97f8cSStefano Zampini ierr = KSPSolve(deluxe_ctx->par_ksp[i],deluxe_ctx->par_vec[i],deluxe_ctx->par_vec[i]);CHKERRQ(ierr); 92234a97f8cSStefano Zampini ierr = VecScatterBegin(deluxe_ctx->par_scctx_p[i],deluxe_ctx->par_vec[i],sub_schurs->work1[subidx],INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 92334a97f8cSStefano Zampini ierr = VecScatterEnd(deluxe_ctx->par_scctx_p[i],deluxe_ctx->par_vec[i],sub_schurs->work1[subidx],INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 92434a97f8cSStefano Zampini ierr = VecSet(deluxe_ctx->par_vec[i],0.0);CHKERRQ(ierr); 92534a97f8cSStefano Zampini ierr = VecScatterBegin(deluxe_ctx->par_scctx_p[i],sub_schurs->work1[subidx],deluxe_ctx->par_vec[i],INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 92634a97f8cSStefano Zampini ierr = VecScatterEnd(deluxe_ctx->par_scctx_p[i],sub_schurs->work1[subidx],deluxe_ctx->par_vec[i],INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 92734a97f8cSStefano Zampini /* test partition of unity */ 92834a97f8cSStefano Zampini ierr = VecAXPY(test_vec,-1.0,deluxe_ctx->par_vec[i]);CHKERRQ(ierr); 92934a97f8cSStefano Zampini ierr = VecNorm(test_vec,NORM_INFINITY,&error);CHKERRQ(ierr); 930c63f45b2SStefano Zampini if (PetscAbsReal(error) > 1.e-2) { 93134a97f8cSStefano Zampini /* ierr = VecView(test_vec,0);CHKERRQ(ierr); */ 93234a97f8cSStefano Zampini error_found = PETSC_TRUE; 93334a97f8cSStefano Zampini } 93434a97f8cSStefano Zampini ierr = VecDestroy(&test_vec);CHKERRQ(ierr); 93534a97f8cSStefano Zampini } 93634a97f8cSStefano Zampini if (error_found) { 93734a97f8cSStefano Zampini ierr = PetscViewerASCIISynchronizedPrintf(pcbddc->dbg_viewer,"Error testing local schur for color %d and subdomain %d\n",i,PetscGlobalRank);CHKERRQ(ierr); 93834a97f8cSStefano Zampini } 93934a97f8cSStefano Zampini ierr = PetscViewerFlush(pcbddc->dbg_viewer);CHKERRQ(ierr); 94034a97f8cSStefano Zampini } 94134a97f8cSStefano Zampini } 94234a97f8cSStefano Zampini PetscFunctionReturn(0); 94334a97f8cSStefano Zampini } 94434a97f8cSStefano Zampini 94534a97f8cSStefano Zampini 94634a97f8cSStefano Zampini #undef __FUNCT__ 94734a97f8cSStefano Zampini #define __FUNCT__ "PCBDDCScalingSetUp_Deluxe_Seq" 948212decf3SStefano Zampini static PetscErrorCode PCBDDCScalingSetUp_Deluxe_Seq(PC pc,PetscInt n_local_sequential_problems,PetscInt n_sequential_problems,PetscInt global_sequential[],PetscInt local_sequential[]) 94934a97f8cSStefano Zampini { 95034a97f8cSStefano Zampini PC_BDDC *pcbddc=(PC_BDDC*)pc->data; 95134a97f8cSStefano Zampini PCBDDCDeluxeScaling deluxe_ctx=pcbddc->deluxe_ctx; 952b1b3d7a2SStefano Zampini PCBDDCSubSchurs sub_schurs = pcbddc->sub_schurs[1]; 953812aaec9SStefano Zampini ISLocalToGlobalMapping l2gmap_subsets; 95434a97f8cSStefano Zampini Mat global_schur_subsets,*submat_global_schur_subsets,work_mat; 95534a97f8cSStefano Zampini IS is_to,is_from; 956212decf3SStefano Zampini PetscScalar *fill_vals; 957812aaec9SStefano Zampini PetscInt *nnz,*all_local_idx_G,*all_local_idx_B,*all_local_idx_N,*all_permutation_G,*dummy_idx; 958212decf3SStefano Zampini PetscInt i,j,local_problem_index; 959812aaec9SStefano Zampini PetscInt subset_size,max_subset_size; 96034a97f8cSStefano Zampini PetscInt local_size,global_size; 961ccdc0be0SStefano Zampini PC pc_temp; 962ccdc0be0SStefano Zampini MatSolverPackage solver=NULL; 96334a97f8cSStefano Zampini char ksp_prefix[256]; 96434a97f8cSStefano Zampini size_t len; 96534a97f8cSStefano Zampini PetscErrorCode ierr; 96634a97f8cSStefano Zampini 96734a97f8cSStefano Zampini PetscFunctionBegin; 96834a97f8cSStefano Zampini if (!n_sequential_problems) { 96934a97f8cSStefano Zampini PetscFunctionReturn(0); 97034a97f8cSStefano Zampini } 971812aaec9SStefano Zampini 97234a97f8cSStefano Zampini /* Get info on subset sizes and sum of all subsets sizes */ 97334a97f8cSStefano Zampini max_subset_size = 0; 97434a97f8cSStefano Zampini local_size = 0; 97534a97f8cSStefano Zampini for (i=0;i<n_local_sequential_problems;i++) { 97634a97f8cSStefano Zampini local_problem_index = local_sequential[i]; 97734a97f8cSStefano Zampini ierr = ISGetLocalSize(sub_schurs->is_AEj_B[local_problem_index],&subset_size);CHKERRQ(ierr); 97834a97f8cSStefano Zampini max_subset_size = PetscMax(subset_size,max_subset_size); 97934a97f8cSStefano Zampini local_size += subset_size; 98034a97f8cSStefano Zampini } 98134a97f8cSStefano Zampini 98234a97f8cSStefano Zampini /* Work arrays for local indices */ 98334a97f8cSStefano Zampini ierr = PetscMalloc1(local_size,&all_local_idx_B);CHKERRQ(ierr); 98434a97f8cSStefano Zampini ierr = PetscMalloc1(local_size,&all_local_idx_N);CHKERRQ(ierr); 985812aaec9SStefano Zampini ierr = PetscMalloc1(local_size,&nnz);CHKERRQ(ierr); 98634a97f8cSStefano Zampini 98734a97f8cSStefano Zampini /* Get local indices in local whole numbering and local boundary numbering */ 98834a97f8cSStefano Zampini local_size = 0; 98934a97f8cSStefano Zampini for (i=0;i<n_local_sequential_problems;i++) { 9905e8657edSStefano Zampini PC_IS *pcis=(PC_IS*)pc->data; 99134a97f8cSStefano Zampini PetscInt *idxs; 99234a97f8cSStefano Zampini /* get info on local problem */ 99334a97f8cSStefano Zampini local_problem_index = local_sequential[i]; 99434a97f8cSStefano Zampini ierr = ISGetLocalSize(sub_schurs->is_AEj_B[local_problem_index],&subset_size);CHKERRQ(ierr); 99534a97f8cSStefano Zampini ierr = ISGetIndices(sub_schurs->is_AEj_B[local_problem_index],(const PetscInt**)&idxs);CHKERRQ(ierr); 99634a97f8cSStefano Zampini /* subset indices in local numbering */ 99734a97f8cSStefano Zampini ierr = PetscMemcpy(all_local_idx_N+local_size,idxs,subset_size*sizeof(PetscInt));CHKERRQ(ierr); 99834a97f8cSStefano Zampini /* subset indices in local boundary numbering */ 9995e8657edSStefano Zampini ierr = ISGlobalToLocalMappingApply(pcis->BtoNmap,IS_GTOLM_DROP,subset_size,idxs,&j,&all_local_idx_B[local_size]);CHKERRQ(ierr); 100034a97f8cSStefano Zampini ierr = ISRestoreIndices(sub_schurs->is_AEj_B[local_problem_index],(const PetscInt**)&idxs);CHKERRQ(ierr); 100134a97f8cSStefano Zampini if (j != subset_size) { 100234a97f8cSStefano Zampini SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Error in BDDC deluxe serial %d (BtoNmap)! %d != %d\n",local_problem_index,subset_size,j); 100334a97f8cSStefano Zampini } 1004812aaec9SStefano Zampini for (j=0;j<subset_size;j++) nnz[local_size+j] = subset_size; 100534a97f8cSStefano Zampini local_size += subset_size; 100634a97f8cSStefano Zampini } 100734a97f8cSStefano Zampini 100834a97f8cSStefano Zampini /* Number dofs on all subsets (parallel) and sort numbering */ 100934a97f8cSStefano Zampini ierr = PCBDDCSubsetNumbering(PetscObjectComm((PetscObject)pc),pcbddc->mat_graph->l2gmap,local_size,all_local_idx_N,PETSC_NULL,&global_size,&all_local_idx_G);CHKERRQ(ierr); 101034a97f8cSStefano Zampini ierr = PetscMalloc1(local_size,&all_permutation_G);CHKERRQ(ierr); 101134a97f8cSStefano Zampini for (i=0;i<local_size;i++) { 101234a97f8cSStefano Zampini all_permutation_G[i]=i; 101334a97f8cSStefano Zampini } 101434a97f8cSStefano Zampini ierr = PetscSortIntWithPermutation(local_size,all_local_idx_G,all_permutation_G);CHKERRQ(ierr); 101534a97f8cSStefano Zampini 101634a97f8cSStefano Zampini /* Local matrix of all local Schur on subsets */ 1017*41c3ba1bSStefano Zampini ierr = MatCreate(PETSC_COMM_SELF,&sub_schurs->S_Ej_all);CHKERRQ(ierr); 1018*41c3ba1bSStefano Zampini ierr = MatSetSizes(sub_schurs->S_Ej_all,PETSC_DECIDE,PETSC_DECIDE,local_size,local_size);CHKERRQ(ierr); 1019*41c3ba1bSStefano Zampini ierr = MatSetType(sub_schurs->S_Ej_all,MATAIJ);CHKERRQ(ierr); 1020*41c3ba1bSStefano Zampini ierr = MatSetOption(sub_schurs->S_Ej_all,MAT_ROW_ORIENTED,PETSC_FALSE);CHKERRQ(ierr); 1021*41c3ba1bSStefano Zampini ierr = MatSeqAIJSetPreallocation(sub_schurs->S_Ej_all,0,nnz);CHKERRQ(ierr); 1022812aaec9SStefano Zampini ierr = PetscFree(nnz);CHKERRQ(ierr); 102334a97f8cSStefano Zampini 102434a97f8cSStefano Zampini /* Work arrays */ 102534a97f8cSStefano Zampini ierr = PetscMalloc2(max_subset_size,&dummy_idx,max_subset_size*max_subset_size,&fill_vals);CHKERRQ(ierr); 102634a97f8cSStefano Zampini 102734a97f8cSStefano Zampini /* Loop on local problems to compute Schur complements explicitly */ 102834a97f8cSStefano Zampini local_size = 0; 102934a97f8cSStefano Zampini for (i=0;i<n_local_sequential_problems;i++) { 103034a97f8cSStefano Zampini /* get info on local problem */ 103134a97f8cSStefano Zampini local_problem_index = local_sequential[i]; 103234a97f8cSStefano Zampini ierr = ISGetLocalSize(sub_schurs->is_AEj_B[local_problem_index],&subset_size);CHKERRQ(ierr); 103334a97f8cSStefano Zampini /* local Schur */ 103434a97f8cSStefano Zampini for (j=0;j<subset_size;j++) { 103534a97f8cSStefano Zampini ierr = VecSet(sub_schurs->work1[local_problem_index],0.0);CHKERRQ(ierr); 103634a97f8cSStefano Zampini ierr = VecSetValue(sub_schurs->work1[local_problem_index],j,1.0,INSERT_VALUES);CHKERRQ(ierr); 1037212decf3SStefano Zampini ierr = VecPlaceArray(sub_schurs->work2[local_problem_index],&fill_vals[j*subset_size]);CHKERRQ(ierr); 103834a97f8cSStefano Zampini ierr = MatMult(sub_schurs->S_Ej[local_problem_index],sub_schurs->work1[local_problem_index],sub_schurs->work2[local_problem_index]);CHKERRQ(ierr); 1039212decf3SStefano Zampini ierr = VecResetArray(sub_schurs->work2[local_problem_index]);CHKERRQ(ierr); 104034a97f8cSStefano Zampini } 104134a97f8cSStefano Zampini for (j=0;j<subset_size;j++) { 104234a97f8cSStefano Zampini dummy_idx[j]=local_size+j; 104334a97f8cSStefano Zampini } 1044*41c3ba1bSStefano Zampini ierr = MatSetValues(sub_schurs->S_Ej_all,subset_size,dummy_idx,subset_size,dummy_idx,fill_vals,INSERT_VALUES);CHKERRQ(ierr); 104534a97f8cSStefano Zampini local_size += subset_size; 104634a97f8cSStefano Zampini } 1047812aaec9SStefano Zampini ierr = PetscFree2(dummy_idx,fill_vals);CHKERRQ(ierr); 1048*41c3ba1bSStefano Zampini ierr = MatAssemblyBegin(sub_schurs->S_Ej_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1049*41c3ba1bSStefano Zampini ierr = MatAssemblyEnd(sub_schurs->S_Ej_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1050812aaec9SStefano Zampini 1051812aaec9SStefano Zampini /* Global matrix of all assembled Schur on subsets */ 1052812aaec9SStefano Zampini ierr = ISLocalToGlobalMappingCreate(PetscObjectComm((PetscObject)pc),1,local_size,all_local_idx_G,PETSC_COPY_VALUES,&l2gmap_subsets);CHKERRQ(ierr); 1053812aaec9SStefano Zampini ierr = MatCreateIS(PetscObjectComm((PetscObject)pc),1,PETSC_DECIDE,PETSC_DECIDE,global_size,global_size,l2gmap_subsets,&work_mat);CHKERRQ(ierr); 1054*41c3ba1bSStefano Zampini ierr = ISLocalToGlobalMappingDestroy(&l2gmap_subsets);CHKERRQ(ierr); 1055*41c3ba1bSStefano Zampini ierr = MatISSetLocalMat(work_mat,sub_schurs->S_Ej_all);CHKERRQ(ierr); 1056812aaec9SStefano Zampini ierr = MatISGetMPIXAIJ(work_mat,MAT_INITIAL_MATRIX,&global_schur_subsets);CHKERRQ(ierr); 1057812aaec9SStefano Zampini ierr = MatDestroy(&work_mat);CHKERRQ(ierr); 105834a97f8cSStefano Zampini /* Create work vectors for sequential part of deluxe */ 1059*41c3ba1bSStefano Zampini ierr = MatCreateVecs(sub_schurs->S_Ej_all,&deluxe_ctx->seq_work1,&deluxe_ctx->seq_work2);CHKERRQ(ierr); 106034a97f8cSStefano Zampini 106134a97f8cSStefano Zampini /* Compute deluxe sequential scatter */ 106234a97f8cSStefano Zampini ierr = ISCreateGeneral(PETSC_COMM_SELF,local_size,all_local_idx_B,PETSC_OWN_POINTER,&is_from);CHKERRQ(ierr); 106334a97f8cSStefano Zampini ierr = VecScatterCreate(pcbddc->work_scaling,is_from,deluxe_ctx->seq_work1,NULL,&deluxe_ctx->seq_scctx);CHKERRQ(ierr); 106434a97f8cSStefano Zampini ierr = ISDestroy(&is_from);CHKERRQ(ierr); 106534a97f8cSStefano Zampini 106634a97f8cSStefano Zampini /* Get local part of (\sum_j S_Ej) */ 106734a97f8cSStefano Zampini for (i=0;i<local_size;i++) { 106834a97f8cSStefano Zampini all_local_idx_N[i] = all_local_idx_G[all_permutation_G[i]]; 106934a97f8cSStefano Zampini } 107034a97f8cSStefano Zampini ierr = ISCreateGeneral(PetscObjectComm((PetscObject)pc),local_size,all_local_idx_N,PETSC_OWN_POINTER,&is_to);CHKERRQ(ierr); 107134a97f8cSStefano Zampini ierr = MatGetSubMatrices(global_schur_subsets,1,&is_to,&is_to,MAT_INITIAL_MATRIX,&submat_global_schur_subsets);CHKERRQ(ierr); 107234a97f8cSStefano Zampini ierr = MatDestroy(&global_schur_subsets);CHKERRQ(ierr); 107334a97f8cSStefano Zampini ierr = ISDestroy(&is_to);CHKERRQ(ierr); 107434a97f8cSStefano Zampini for (i=0;i<local_size;i++) { 107534a97f8cSStefano Zampini all_local_idx_G[all_permutation_G[i]] = i; 107634a97f8cSStefano Zampini } 107734a97f8cSStefano Zampini ierr = ISCreateGeneral(PETSC_COMM_SELF,local_size,all_local_idx_G,PETSC_OWN_POINTER,&is_from);CHKERRQ(ierr); 107834a97f8cSStefano Zampini ierr = ISSetPermutation(is_from);CHKERRQ(ierr); 1079*41c3ba1bSStefano Zampini ierr = MatPermute(submat_global_schur_subsets[0],is_from,is_from,&sub_schurs->sum_S_Ej_all);CHKERRQ(ierr); 108034a97f8cSStefano Zampini ierr = MatDestroyMatrices(1,&submat_global_schur_subsets);CHKERRQ(ierr); 108134a97f8cSStefano Zampini ierr = ISDestroy(&is_from);CHKERRQ(ierr); 108234a97f8cSStefano Zampini ierr = PetscFree(all_permutation_G);CHKERRQ(ierr); 108334a97f8cSStefano Zampini 108434a97f8cSStefano Zampini /* Create KSP object for sequential part of deluxe scaling */ 108534a97f8cSStefano Zampini ierr = KSPCreate(PETSC_COMM_SELF,&deluxe_ctx->seq_ksp);CHKERRQ(ierr); 1086*41c3ba1bSStefano Zampini ierr = KSPSetOperators(deluxe_ctx->seq_ksp,sub_schurs->sum_S_Ej_all,sub_schurs->sum_S_Ej_all);CHKERRQ(ierr); 108734a97f8cSStefano Zampini ierr = KSPSetType(deluxe_ctx->seq_ksp,KSPPREONLY);CHKERRQ(ierr); 1088ccdc0be0SStefano Zampini ierr = KSPGetPC(deluxe_ctx->seq_ksp,&pc_temp);CHKERRQ(ierr); 1089ccdc0be0SStefano Zampini ierr = PCSetType(pc_temp,PCLU);CHKERRQ(ierr); 1090ccdc0be0SStefano Zampini ierr = KSPGetPC(pcbddc->ksp_D,&pc_temp);CHKERRQ(ierr); 1091ccdc0be0SStefano Zampini ierr = PCFactorGetMatSolverPackage(pc_temp,(const MatSolverPackage*)&solver);CHKERRQ(ierr); 10921e9c79c2SStefano Zampini if (solver && local_size) { /* if local_size is null, some external packages will report errors */ 1093a4e17c67SStefano Zampini PC new_pc; 1094a4e17c67SStefano Zampini PCType type; 1095a4e17c67SStefano Zampini ierr = PCGetType(pc_temp,&type);CHKERRQ(ierr); 1096a4e17c67SStefano Zampini ierr = KSPGetPC(deluxe_ctx->seq_ksp,&new_pc);CHKERRQ(ierr); 1097a4e17c67SStefano Zampini ierr = PCSetType(new_pc,type);CHKERRQ(ierr); 1098a4e17c67SStefano Zampini ierr = PCFactorSetMatSolverPackage(new_pc,solver);CHKERRQ(ierr); 1099ccdc0be0SStefano Zampini } 110034a97f8cSStefano Zampini ierr = PetscStrlen(((PetscObject)(pcbddc->ksp_D))->prefix,&len);CHKERRQ(ierr); 110134a97f8cSStefano Zampini len -= 10; /* remove "dirichlet_" */ 11028856534dSStefano Zampini ierr = PetscStrncpy(ksp_prefix,((PetscObject)(pcbddc->ksp_D))->prefix,len+1);CHKERRQ(ierr); 110334a97f8cSStefano Zampini ierr = PetscStrcat(ksp_prefix,"deluxe_seq_");CHKERRQ(ierr); 110434a97f8cSStefano Zampini ierr = KSPSetOptionsPrefix(deluxe_ctx->seq_ksp,ksp_prefix);CHKERRQ(ierr); 11051e9c79c2SStefano Zampini if (local_size) { 110634a97f8cSStefano Zampini ierr = KSPSetFromOptions(deluxe_ctx->seq_ksp);CHKERRQ(ierr); 11071e9c79c2SStefano Zampini } 110834a97f8cSStefano Zampini ierr = KSPSetUp(deluxe_ctx->seq_ksp);CHKERRQ(ierr); 110934a97f8cSStefano Zampini PetscFunctionReturn(0); 111034a97f8cSStefano Zampini } 1111