1ab5c6b0cSJed Brown #include <../src/ksp/pc/impls/bddc/bddc.h> 2ab5c6b0cSJed Brown #include <../src/ksp/pc/impls/bddc/bddcprivate.h> 3*72b8c272SStefano Zampini #include <petscblaslapack.h> 4674ae819SStefano Zampini 534a97f8cSStefano Zampini /* prototypes for deluxe functions */ 634a97f8cSStefano Zampini static PetscErrorCode PCBDDCScalingCreate_Deluxe(PC); 734a97f8cSStefano Zampini static PetscErrorCode PCBDDCScalingDestroy_Deluxe(PC); 834a97f8cSStefano Zampini static PetscErrorCode PCBDDCScalingSetUp_Deluxe(PC); 95a95e1ceSStefano Zampini static PetscErrorCode PCBDDCScalingSetUp_Deluxe_Private(PC); 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; 36674ae819SStefano Zampini PetscErrorCode ierr; 37674ae819SStefano Zampini 38674ae819SStefano Zampini PetscFunctionBegin; 395a95e1ceSStefano Zampini ierr = VecSet(pcbddc->work_scaling,0.0);CHKERRQ(ierr); 405a95e1ceSStefano Zampini ierr = VecSet(y,0.0);CHKERRQ(ierr); 415a95e1ceSStefano Zampini if (deluxe_ctx->n_simple) { /* scale deluxe vertices using diagonal scaling */ 425a95e1ceSStefano Zampini PetscInt i; 432b095fd8SStefano Zampini const PetscScalar *array_x,*array_D; 442b095fd8SStefano Zampini PetscScalar *array; 452b095fd8SStefano Zampini ierr = VecGetArrayRead(x,&array_x);CHKERRQ(ierr); 462b095fd8SStefano Zampini ierr = VecGetArrayRead(pcis->D,&array_D);CHKERRQ(ierr); 47674ae819SStefano Zampini ierr = VecGetArray(pcbddc->work_scaling,&array);CHKERRQ(ierr); 48674ae819SStefano Zampini for (i=0;i<deluxe_ctx->n_simple;i++) { 49674ae819SStefano Zampini array[deluxe_ctx->idx_simple_B[i]] = array_x[deluxe_ctx->idx_simple_B[i]]*array_D[deluxe_ctx->idx_simple_B[i]]; 50674ae819SStefano Zampini } 51674ae819SStefano Zampini ierr = VecRestoreArray(pcbddc->work_scaling,&array);CHKERRQ(ierr); 522b095fd8SStefano Zampini ierr = VecRestoreArrayRead(pcis->D,&array_D);CHKERRQ(ierr); 532b095fd8SStefano Zampini ierr = VecRestoreArrayRead(x,&array_x);CHKERRQ(ierr); 5434a97f8cSStefano Zampini } 55ac632422SStefano Zampini /* sequential part : all problems and Schur applications collapsed into a single matrix vector multiplication or a matvec and a solve */ 565a95e1ceSStefano Zampini if (deluxe_ctx->seq_mat) { 57*72b8c272SStefano Zampini PetscInt i; 58*72b8c272SStefano Zampini for (i=0;i<deluxe_ctx->seq_n;i++) { 5916e386b8SStefano Zampini if (deluxe_ctx->change) { 60*72b8c272SStefano Zampini ierr = VecScatterBegin(deluxe_ctx->seq_scctx[i],x,deluxe_ctx->seq_work2[i],INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 61*72b8c272SStefano Zampini ierr = VecScatterEnd(deluxe_ctx->seq_scctx[i],x,deluxe_ctx->seq_work2[i],INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 62*72b8c272SStefano Zampini if (deluxe_ctx->change_with_qr) { 63*72b8c272SStefano Zampini Mat change; 64*72b8c272SStefano Zampini 65*72b8c272SStefano Zampini ierr = KSPGetOperators(deluxe_ctx->change[i],&change,NULL);CHKERRQ(ierr); 66*72b8c272SStefano Zampini ierr = MatMultTranspose(change,deluxe_ctx->seq_work2[i],deluxe_ctx->seq_work1[i]);CHKERRQ(ierr); 67*72b8c272SStefano Zampini } else { 68*72b8c272SStefano Zampini ierr = KSPSolve(deluxe_ctx->change[i],deluxe_ctx->seq_work2[i],deluxe_ctx->seq_work1[i]);CHKERRQ(ierr); 6916e386b8SStefano Zampini } 70*72b8c272SStefano Zampini } else { 71*72b8c272SStefano Zampini ierr = VecScatterBegin(deluxe_ctx->seq_scctx[i],x,deluxe_ctx->seq_work1[i],INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 72*72b8c272SStefano Zampini ierr = VecScatterEnd(deluxe_ctx->seq_scctx[i],x,deluxe_ctx->seq_work1[i],INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 73*72b8c272SStefano Zampini } 74*72b8c272SStefano Zampini ierr = MatMultTranspose(deluxe_ctx->seq_mat[i],deluxe_ctx->seq_work1[i],deluxe_ctx->seq_work2[i]);CHKERRQ(ierr); 75*72b8c272SStefano Zampini if (deluxe_ctx->seq_mat_inv_sum[i]) { 76*72b8c272SStefano Zampini PetscScalar *x; 77*72b8c272SStefano Zampini 78*72b8c272SStefano Zampini ierr = VecGetArray(deluxe_ctx->seq_work2[i],&x);CHKERRQ(ierr); 79*72b8c272SStefano Zampini ierr = VecPlaceArray(deluxe_ctx->seq_work1[i],x);CHKERRQ(ierr); 80*72b8c272SStefano Zampini ierr = VecRestoreArray(deluxe_ctx->seq_work2[i],&x);CHKERRQ(ierr); 81*72b8c272SStefano Zampini ierr = MatSolveTranspose(deluxe_ctx->seq_mat_inv_sum[i],deluxe_ctx->seq_work1[i],deluxe_ctx->seq_work2[i]);CHKERRQ(ierr); 82*72b8c272SStefano Zampini ierr = VecResetArray(deluxe_ctx->seq_work1[i]);CHKERRQ(ierr); 83ac632422SStefano Zampini } 8416e386b8SStefano Zampini if (deluxe_ctx->change) { 8516e386b8SStefano Zampini Mat change; 8616e386b8SStefano Zampini 87*72b8c272SStefano Zampini ierr = KSPGetOperators(deluxe_ctx->change[i],&change,NULL);CHKERRQ(ierr); 88*72b8c272SStefano Zampini ierr = MatMult(change,deluxe_ctx->seq_work2[i],deluxe_ctx->seq_work1[i]);CHKERRQ(ierr); 89*72b8c272SStefano Zampini ierr = VecScatterBegin(deluxe_ctx->seq_scctx[i],deluxe_ctx->seq_work1[i],pcbddc->work_scaling,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 90*72b8c272SStefano Zampini ierr = VecScatterEnd(deluxe_ctx->seq_scctx[i],deluxe_ctx->seq_work1[i],pcbddc->work_scaling,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 9116e386b8SStefano Zampini } else { 92*72b8c272SStefano Zampini ierr = VecScatterBegin(deluxe_ctx->seq_scctx[i],deluxe_ctx->seq_work2[i],pcbddc->work_scaling,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 93*72b8c272SStefano Zampini ierr = VecScatterEnd(deluxe_ctx->seq_scctx[i],deluxe_ctx->seq_work2[i],pcbddc->work_scaling,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 94*72b8c272SStefano Zampini } 95674ae819SStefano Zampini } 9616e386b8SStefano Zampini } 97674ae819SStefano Zampini /* put local boundary part in global vector */ 98674ae819SStefano Zampini ierr = VecScatterBegin(pcis->global_to_B,pcbddc->work_scaling,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 99674ae819SStefano Zampini ierr = VecScatterEnd(pcis->global_to_B,pcbddc->work_scaling,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 100674ae819SStefano Zampini PetscFunctionReturn(0); 101674ae819SStefano Zampini } 102674ae819SStefano Zampini 103674ae819SStefano Zampini #undef __FUNCT__ 104674ae819SStefano Zampini #define __FUNCT__ "PCBDDCScalingExtension" 105674ae819SStefano Zampini PetscErrorCode PCBDDCScalingExtension(PC pc, Vec local_interface_vector, Vec global_vector) 106674ae819SStefano Zampini { 107674ae819SStefano Zampini PC_BDDC *pcbddc=(PC_BDDC*)pc->data; 108674ae819SStefano Zampini PetscErrorCode ierr; 109674ae819SStefano Zampini 110674ae819SStefano Zampini PetscFunctionBegin; 111674ae819SStefano Zampini PetscValidHeaderSpecific(pc,PC_CLASSID,1); 112674ae819SStefano Zampini PetscValidHeaderSpecific(local_interface_vector,VEC_CLASSID,2); 113674ae819SStefano Zampini PetscValidHeaderSpecific(global_vector,VEC_CLASSID,3); 114674ae819SStefano Zampini if (local_interface_vector == pcbddc->work_scaling) { 115a07ea27aSStefano Zampini SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Local vector cannot be pcbddc->work_scaling!\n"); 116674ae819SStefano Zampini } 117674ae819SStefano Zampini ierr = PetscTryMethod(pc,"PCBDDCScalingExtension_C",(PC,Vec,Vec),(pc,local_interface_vector,global_vector));CHKERRQ(ierr); 118674ae819SStefano Zampini PetscFunctionReturn(0); 119674ae819SStefano Zampini } 120674ae819SStefano Zampini 121674ae819SStefano Zampini #undef __FUNCT__ 122674ae819SStefano Zampini #define __FUNCT__ "PCBDDCScalingRestriction_Basic" 123674ae819SStefano Zampini static PetscErrorCode PCBDDCScalingRestriction_Basic(PC pc, Vec global_vector, Vec local_interface_vector) 124674ae819SStefano Zampini { 125674ae819SStefano Zampini PetscErrorCode ierr; 126674ae819SStefano Zampini PC_IS* pcis = (PC_IS*)pc->data; 127674ae819SStefano Zampini 128674ae819SStefano Zampini PetscFunctionBegin; 129674ae819SStefano Zampini ierr = VecScatterBegin(pcis->global_to_B,global_vector,local_interface_vector,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 130674ae819SStefano Zampini ierr = VecScatterEnd(pcis->global_to_B,global_vector,local_interface_vector,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 131674ae819SStefano Zampini /* Apply partition of unity */ 132674ae819SStefano Zampini ierr = VecPointwiseMult(local_interface_vector,pcis->D,local_interface_vector);CHKERRQ(ierr); 133674ae819SStefano Zampini PetscFunctionReturn(0); 134674ae819SStefano Zampini } 135674ae819SStefano Zampini 136674ae819SStefano Zampini #undef __FUNCT__ 137674ae819SStefano Zampini #define __FUNCT__ "PCBDDCScalingRestriction_Deluxe" 138674ae819SStefano Zampini static PetscErrorCode PCBDDCScalingRestriction_Deluxe(PC pc, Vec x, Vec y) 139674ae819SStefano Zampini { 140674ae819SStefano Zampini PC_IS* pcis=(PC_IS*)pc->data; 141674ae819SStefano Zampini PC_BDDC* pcbddc=(PC_BDDC*)pc->data; 142674ae819SStefano Zampini PCBDDCDeluxeScaling deluxe_ctx = pcbddc->deluxe_ctx; 143674ae819SStefano Zampini PetscErrorCode ierr; 144674ae819SStefano Zampini 145674ae819SStefano Zampini PetscFunctionBegin; 146674ae819SStefano Zampini /* get local boundary part of global vector */ 147674ae819SStefano Zampini ierr = VecScatterBegin(pcis->global_to_B,x,y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 148674ae819SStefano Zampini ierr = VecScatterEnd(pcis->global_to_B,x,y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1495a95e1ceSStefano Zampini if (deluxe_ctx->n_simple) { /* scale deluxe vertices using diagonal scaling */ 1505a95e1ceSStefano Zampini PetscInt i; 1512b095fd8SStefano Zampini PetscScalar *array_y; 1522b095fd8SStefano Zampini const PetscScalar *array_D; 153674ae819SStefano Zampini ierr = VecGetArray(y,&array_y);CHKERRQ(ierr); 1542b095fd8SStefano Zampini ierr = VecGetArrayRead(pcis->D,&array_D);CHKERRQ(ierr); 155674ae819SStefano Zampini for (i=0;i<deluxe_ctx->n_simple;i++) { 156674ae819SStefano Zampini array_y[deluxe_ctx->idx_simple_B[i]] *= array_D[deluxe_ctx->idx_simple_B[i]]; 157674ae819SStefano Zampini } 1582b095fd8SStefano Zampini ierr = VecRestoreArrayRead(pcis->D,&array_D);CHKERRQ(ierr); 159674ae819SStefano Zampini ierr = VecRestoreArray(y,&array_y);CHKERRQ(ierr); 16034a97f8cSStefano Zampini } 16134a97f8cSStefano Zampini /* sequential part : all problems and Schur applications collapsed into a single matrix vector multiplication and ksp solution */ 1625a95e1ceSStefano Zampini if (deluxe_ctx->seq_mat) { 163*72b8c272SStefano Zampini PetscInt i; 164*72b8c272SStefano Zampini for (i=0;i<deluxe_ctx->seq_n;i++) { 16516e386b8SStefano Zampini if (deluxe_ctx->change) { 16616e386b8SStefano Zampini Mat change; 16716e386b8SStefano Zampini 168*72b8c272SStefano Zampini ierr = VecScatterBegin(deluxe_ctx->seq_scctx[i],y,deluxe_ctx->seq_work2[i],INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 169*72b8c272SStefano Zampini ierr = VecScatterEnd(deluxe_ctx->seq_scctx[i],y,deluxe_ctx->seq_work2[i],INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 170*72b8c272SStefano Zampini ierr = KSPGetOperators(deluxe_ctx->change[i],&change,NULL);CHKERRQ(ierr); 171*72b8c272SStefano Zampini ierr = MatMultTranspose(change,deluxe_ctx->seq_work2[i],deluxe_ctx->seq_work1[i]);CHKERRQ(ierr); 17216e386b8SStefano Zampini } else { 173*72b8c272SStefano Zampini ierr = VecScatterBegin(deluxe_ctx->seq_scctx[i],y,deluxe_ctx->seq_work1[i],INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 174*72b8c272SStefano Zampini ierr = VecScatterEnd(deluxe_ctx->seq_scctx[i],y,deluxe_ctx->seq_work1[i],INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 17516e386b8SStefano Zampini } 176*72b8c272SStefano Zampini if (deluxe_ctx->seq_mat_inv_sum[i]) { 177*72b8c272SStefano Zampini PetscScalar *x; 178*72b8c272SStefano Zampini 179*72b8c272SStefano Zampini ierr = VecGetArray(deluxe_ctx->seq_work1[i],&x);CHKERRQ(ierr); 180*72b8c272SStefano Zampini ierr = VecPlaceArray(deluxe_ctx->seq_work2[i],x);CHKERRQ(ierr); 181*72b8c272SStefano Zampini ierr = VecRestoreArray(deluxe_ctx->seq_work1[i],&x);CHKERRQ(ierr); 182*72b8c272SStefano Zampini ierr = MatSolve(deluxe_ctx->seq_mat_inv_sum[i],deluxe_ctx->seq_work1[i],deluxe_ctx->seq_work2[i]);CHKERRQ(ierr); 183*72b8c272SStefano Zampini ierr = VecResetArray(deluxe_ctx->seq_work2[i]);CHKERRQ(ierr); 184ac632422SStefano Zampini } 185*72b8c272SStefano Zampini ierr = MatMult(deluxe_ctx->seq_mat[i],deluxe_ctx->seq_work1[i],deluxe_ctx->seq_work2[i]);CHKERRQ(ierr); 18616e386b8SStefano Zampini if (deluxe_ctx->change) { 187*72b8c272SStefano Zampini if (deluxe_ctx->change_with_qr) { 188*72b8c272SStefano Zampini Mat change; 189*72b8c272SStefano Zampini 190*72b8c272SStefano Zampini ierr = KSPGetOperators(deluxe_ctx->change[i],&change,NULL);CHKERRQ(ierr); 191*72b8c272SStefano Zampini ierr = MatMult(change,deluxe_ctx->seq_work2[i],deluxe_ctx->seq_work1[i]);CHKERRQ(ierr); 192*72b8c272SStefano Zampini } else { 193*72b8c272SStefano Zampini ierr = KSPSolveTranspose(deluxe_ctx->change[i],deluxe_ctx->seq_work2[i],deluxe_ctx->seq_work1[i]);CHKERRQ(ierr); 19416e386b8SStefano Zampini } 195*72b8c272SStefano Zampini ierr = VecScatterBegin(deluxe_ctx->seq_scctx[i],deluxe_ctx->seq_work1[i],y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 196*72b8c272SStefano Zampini ierr = VecScatterEnd(deluxe_ctx->seq_scctx[i],deluxe_ctx->seq_work1[i],y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 197*72b8c272SStefano Zampini } else { 198*72b8c272SStefano Zampini ierr = VecScatterBegin(deluxe_ctx->seq_scctx[i],deluxe_ctx->seq_work2[i],y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 199*72b8c272SStefano Zampini ierr = VecScatterEnd(deluxe_ctx->seq_scctx[i],deluxe_ctx->seq_work2[i],y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 200*72b8c272SStefano Zampini } 201*72b8c272SStefano Zampini } 202674ae819SStefano Zampini } 203674ae819SStefano Zampini PetscFunctionReturn(0); 204674ae819SStefano Zampini } 205674ae819SStefano Zampini 206674ae819SStefano Zampini #undef __FUNCT__ 207674ae819SStefano Zampini #define __FUNCT__ "PCBDDCScalingRestriction" 208674ae819SStefano Zampini PetscErrorCode PCBDDCScalingRestriction(PC pc, Vec global_vector, Vec local_interface_vector) 209674ae819SStefano Zampini { 210674ae819SStefano Zampini PC_BDDC *pcbddc=(PC_BDDC*)pc->data; 211674ae819SStefano Zampini PetscErrorCode ierr; 212674ae819SStefano Zampini 213674ae819SStefano Zampini PetscFunctionBegin; 214674ae819SStefano Zampini PetscValidHeaderSpecific(pc,PC_CLASSID,1); 215674ae819SStefano Zampini PetscValidHeaderSpecific(global_vector,VEC_CLASSID,2); 216674ae819SStefano Zampini PetscValidHeaderSpecific(local_interface_vector,VEC_CLASSID,3); 217674ae819SStefano Zampini if (local_interface_vector == pcbddc->work_scaling) { 218a07ea27aSStefano Zampini SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Local vector cannot be pcbddc->work_scaling!\n"); 219674ae819SStefano Zampini } 220674ae819SStefano Zampini ierr = PetscTryMethod(pc,"PCBDDCScalingRestriction_C",(PC,Vec,Vec),(pc,global_vector,local_interface_vector));CHKERRQ(ierr); 221674ae819SStefano Zampini PetscFunctionReturn(0); 222674ae819SStefano Zampini } 223674ae819SStefano Zampini 224674ae819SStefano Zampini #undef __FUNCT__ 225674ae819SStefano Zampini #define __FUNCT__ "PCBDDCScalingSetUp" 226674ae819SStefano Zampini PetscErrorCode PCBDDCScalingSetUp(PC pc) 227674ae819SStefano Zampini { 228674ae819SStefano Zampini PC_IS* pcis=(PC_IS*)pc->data; 229674ae819SStefano Zampini PC_BDDC* pcbddc=(PC_BDDC*)pc->data; 230674ae819SStefano Zampini PetscErrorCode ierr; 231674ae819SStefano Zampini 232674ae819SStefano Zampini PetscFunctionBegin; 233674ae819SStefano Zampini PetscValidHeaderSpecific(pc,PC_CLASSID,1); 234674ae819SStefano Zampini /* create work vector for the operator */ 23534a97f8cSStefano Zampini ierr = VecDestroy(&pcbddc->work_scaling);CHKERRQ(ierr); 236674ae819SStefano Zampini ierr = VecDuplicate(pcis->vec1_B,&pcbddc->work_scaling);CHKERRQ(ierr); 23734a97f8cSStefano Zampini /* always rebuild pcis->D */ 23828d874f6SStefano Zampini if (pcis->use_stiffness_scaling) { 239674ae819SStefano Zampini ierr = MatGetDiagonal(pcbddc->local_mat,pcis->vec1_N);CHKERRQ(ierr); 240674ae819SStefano Zampini ierr = VecScatterBegin(pcis->N_to_B,pcis->vec1_N,pcis->D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 241674ae819SStefano Zampini ierr = VecScatterEnd(pcis->N_to_B,pcis->vec1_N,pcis->D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 242674ae819SStefano Zampini } 243674ae819SStefano Zampini ierr = VecCopy(pcis->D,pcis->vec1_B);CHKERRQ(ierr); 244674ae819SStefano Zampini ierr = VecSet(pcis->vec1_global,0.0);CHKERRQ(ierr); 245674ae819SStefano Zampini ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_B,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 246674ae819SStefano Zampini ierr = VecScatterEnd(pcis->global_to_B,pcis->vec1_B,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 247674ae819SStefano Zampini ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_global,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 248674ae819SStefano Zampini ierr = VecScatterEnd(pcis->global_to_B,pcis->vec1_global,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 249674ae819SStefano Zampini ierr = VecPointwiseDivide(pcis->D,pcis->D,pcis->vec1_B);CHKERRQ(ierr); 250674ae819SStefano Zampini /* now setup */ 251681e7c04SStefano Zampini if (pcbddc->use_deluxe_scaling) { 25234a97f8cSStefano Zampini if (!pcbddc->deluxe_ctx) { 25334a97f8cSStefano Zampini ierr = PCBDDCScalingCreate_Deluxe(pc);CHKERRQ(ierr); 25434a97f8cSStefano Zampini } 25534a97f8cSStefano Zampini ierr = PCBDDCScalingSetUp_Deluxe(pc);CHKERRQ(ierr); 256674ae819SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCScalingRestriction_C",PCBDDCScalingRestriction_Deluxe);CHKERRQ(ierr); 257674ae819SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCScalingExtension_C",PCBDDCScalingExtension_Deluxe);CHKERRQ(ierr); 258674ae819SStefano Zampini } else { 259674ae819SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCScalingRestriction_C",PCBDDCScalingRestriction_Basic);CHKERRQ(ierr); 260674ae819SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCScalingExtension_C",PCBDDCScalingExtension_Basic);CHKERRQ(ierr); 261674ae819SStefano Zampini } 26234a97f8cSStefano Zampini 263674ae819SStefano Zampini /* test */ 264674ae819SStefano Zampini if (pcbddc->dbg_flag) { 265*72b8c272SStefano Zampini Mat B0_B = NULL; 266*72b8c272SStefano Zampini Vec B0_Bv = NULL, B0_Bv2 = NULL; 26734a97f8cSStefano Zampini Vec vec2_global; 268674ae819SStefano Zampini PetscViewer viewer = pcbddc->dbg_viewer; 26934a97f8cSStefano Zampini PetscReal error; 270674ae819SStefano Zampini 271674ae819SStefano Zampini /* extension -> from local to parallel */ 27234a97f8cSStefano Zampini ierr = VecSet(pcis->vec1_global,0.0);CHKERRQ(ierr); 27334a97f8cSStefano Zampini ierr = VecSetRandom(pcis->vec1_B,NULL);CHKERRQ(ierr); 27434a97f8cSStefano Zampini ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_B,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 27534a97f8cSStefano Zampini ierr = VecScatterEnd(pcis->global_to_B,pcis->vec1_B,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 27634a97f8cSStefano Zampini ierr = VecDuplicate(pcis->vec1_global,&vec2_global);CHKERRQ(ierr); 27734a97f8cSStefano Zampini ierr = VecCopy(pcis->vec1_global,vec2_global);CHKERRQ(ierr); 278674ae819SStefano Zampini ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_global,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 279674ae819SStefano Zampini ierr = VecScatterEnd(pcis->global_to_B,pcis->vec1_global,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 280*72b8c272SStefano Zampini if (pcbddc->benign_n) { 281*72b8c272SStefano Zampini IS is_dummy; 282*72b8c272SStefano Zampini 283*72b8c272SStefano Zampini ierr = ISCreateStride(PETSC_COMM_SELF,pcbddc->benign_n,0,1,&is_dummy);CHKERRQ(ierr); 284*72b8c272SStefano Zampini ierr = MatGetSubMatrix(pcbddc->benign_B0,is_dummy,pcis->is_B_local,MAT_INITIAL_MATRIX,&B0_B);CHKERRQ(ierr); 285*72b8c272SStefano Zampini ierr = ISDestroy(&is_dummy);CHKERRQ(ierr); 286*72b8c272SStefano Zampini ierr = MatCreateVecs(B0_B,NULL,&B0_Bv);CHKERRQ(ierr); 287*72b8c272SStefano Zampini ierr = VecDuplicate(B0_Bv,&B0_Bv2);CHKERRQ(ierr); 288*72b8c272SStefano Zampini ierr = MatMult(B0_B,pcis->vec1_B,B0_Bv);CHKERRQ(ierr); 289*72b8c272SStefano Zampini } 290674ae819SStefano Zampini ierr = PCBDDCScalingExtension(pc,pcis->vec1_B,pcis->vec1_global);CHKERRQ(ierr); 291*72b8c272SStefano Zampini if (pcbddc->benign_saddle_point) { 292*72b8c272SStefano Zampini PetscReal errorl = 0.; 293*72b8c272SStefano Zampini ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_global,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 294*72b8c272SStefano Zampini ierr = VecScatterEnd(pcis->global_to_B,pcis->vec1_global,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 295*72b8c272SStefano Zampini if (pcbddc->benign_n) { 296*72b8c272SStefano Zampini ierr = MatMult(B0_B,pcis->vec1_B,B0_Bv2);CHKERRQ(ierr); 297*72b8c272SStefano Zampini ierr = VecAXPY(B0_Bv,-1.0,B0_Bv2);CHKERRQ(ierr); 298*72b8c272SStefano Zampini ierr = VecNorm(B0_Bv,NORM_INFINITY,&errorl);CHKERRQ(ierr); 299*72b8c272SStefano Zampini } 300*72b8c272SStefano Zampini ierr = MPI_Allreduce(&errorl,&error,1,MPIU_REAL,MPI_SUM,PetscObjectComm((PetscObject)pc));CHKERRQ(ierr); 301*72b8c272SStefano Zampini ierr = PetscViewerASCIIPrintf(viewer,"Error benign extension %1.14e\n",error);CHKERRQ(ierr); 302*72b8c272SStefano Zampini } 30334a97f8cSStefano Zampini ierr = VecAXPY(pcis->vec1_global,-1.0,vec2_global);CHKERRQ(ierr); 30434a97f8cSStefano Zampini ierr = VecNorm(pcis->vec1_global,NORM_INFINITY,&error);CHKERRQ(ierr); 305674ae819SStefano Zampini ierr = PetscViewerASCIIPrintf(viewer,"Error scaling extension %1.14e\n",error);CHKERRQ(ierr); 30634a97f8cSStefano Zampini ierr = VecDestroy(&vec2_global);CHKERRQ(ierr); 30734a97f8cSStefano Zampini 308674ae819SStefano Zampini /* restriction -> from parallel to local */ 309674ae819SStefano Zampini ierr = VecSet(pcis->vec1_global,0.0);CHKERRQ(ierr); 31034a97f8cSStefano Zampini ierr = VecSetRandom(pcis->vec1_B,NULL);CHKERRQ(ierr); 311674ae819SStefano Zampini ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_B,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 312674ae819SStefano Zampini ierr = VecScatterEnd(pcis->global_to_B,pcis->vec1_B,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 31334a97f8cSStefano Zampini ierr = PCBDDCScalingRestriction(pc,pcis->vec1_global,pcis->vec1_B);CHKERRQ(ierr); 31434a97f8cSStefano Zampini ierr = VecScale(pcis->vec1_B,-1.0);CHKERRQ(ierr); 31534a97f8cSStefano Zampini ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_B,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 31634a97f8cSStefano Zampini ierr = VecScatterEnd(pcis->global_to_B,pcis->vec1_B,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 31734a97f8cSStefano Zampini ierr = VecNorm(pcis->vec1_global,NORM_INFINITY,&error);CHKERRQ(ierr); 31834a97f8cSStefano Zampini ierr = PetscViewerASCIIPrintf(viewer,"Error scaling restriction %1.14e\n",error);CHKERRQ(ierr); 319*72b8c272SStefano Zampini ierr = MatDestroy(&B0_B);CHKERRQ(ierr); 320*72b8c272SStefano Zampini ierr = VecDestroy(&B0_Bv);CHKERRQ(ierr); 321*72b8c272SStefano Zampini ierr = VecDestroy(&B0_Bv2);CHKERRQ(ierr); 322674ae819SStefano Zampini } 323674ae819SStefano Zampini PetscFunctionReturn(0); 324674ae819SStefano Zampini } 325674ae819SStefano Zampini 326674ae819SStefano Zampini #undef __FUNCT__ 327674ae819SStefano Zampini #define __FUNCT__ "PCBDDCScalingDestroy" 328674ae819SStefano Zampini PetscErrorCode PCBDDCScalingDestroy(PC pc) 329674ae819SStefano Zampini { 330674ae819SStefano Zampini PC_BDDC* pcbddc=(PC_BDDC*)pc->data; 331674ae819SStefano Zampini PetscErrorCode ierr; 332674ae819SStefano Zampini 333674ae819SStefano Zampini PetscFunctionBegin; 33434a97f8cSStefano Zampini if (pcbddc->deluxe_ctx) { 33534a97f8cSStefano Zampini ierr = PCBDDCScalingDestroy_Deluxe(pc);CHKERRQ(ierr); 336674ae819SStefano Zampini } 337674ae819SStefano Zampini ierr = VecDestroy(&pcbddc->work_scaling);CHKERRQ(ierr); 338674ae819SStefano Zampini /* remove functions */ 339674ae819SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCScalingRestriction_C",NULL);CHKERRQ(ierr); 340674ae819SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCScalingExtension_C",NULL);CHKERRQ(ierr); 341674ae819SStefano Zampini PetscFunctionReturn(0); 342674ae819SStefano Zampini } 343674ae819SStefano Zampini 34434a97f8cSStefano Zampini #undef __FUNCT__ 34534a97f8cSStefano Zampini #define __FUNCT__ "PCBDDCScalingCreate_Deluxe" 34634a97f8cSStefano Zampini static PetscErrorCode PCBDDCScalingCreate_Deluxe(PC pc) 34734a97f8cSStefano Zampini { 34834a97f8cSStefano Zampini PC_BDDC* pcbddc=(PC_BDDC*)pc->data; 34934a97f8cSStefano Zampini PCBDDCDeluxeScaling deluxe_ctx; 35034a97f8cSStefano Zampini PetscErrorCode ierr; 35134a97f8cSStefano Zampini 35234a97f8cSStefano Zampini PetscFunctionBegin; 35334a97f8cSStefano Zampini ierr = PetscNew(&deluxe_ctx);CHKERRQ(ierr); 35434a97f8cSStefano Zampini pcbddc->deluxe_ctx = deluxe_ctx; 35534a97f8cSStefano Zampini PetscFunctionReturn(0); 35634a97f8cSStefano Zampini } 35734a97f8cSStefano Zampini 35834a97f8cSStefano Zampini #undef __FUNCT__ 35934a97f8cSStefano Zampini #define __FUNCT__ "PCBDDCScalingDestroy_Deluxe" 36034a97f8cSStefano Zampini static PetscErrorCode PCBDDCScalingDestroy_Deluxe(PC pc) 36134a97f8cSStefano Zampini { 36234a97f8cSStefano Zampini PC_BDDC* pcbddc=(PC_BDDC*)pc->data; 36334a97f8cSStefano Zampini PetscErrorCode ierr; 36434a97f8cSStefano Zampini 36534a97f8cSStefano Zampini PetscFunctionBegin; 36634a97f8cSStefano Zampini ierr = PCBDDCScalingReset_Deluxe_Solvers(pcbddc->deluxe_ctx);CHKERRQ(ierr); 36734a97f8cSStefano Zampini ierr = PetscFree(pcbddc->deluxe_ctx);CHKERRQ(ierr); 36834a97f8cSStefano Zampini PetscFunctionReturn(0); 36934a97f8cSStefano Zampini } 37034a97f8cSStefano Zampini 37134a97f8cSStefano Zampini #undef __FUNCT__ 37234a97f8cSStefano Zampini #define __FUNCT__ "PCBDDCScalingReset_Deluxe_Solvers" 37334a97f8cSStefano Zampini static PetscErrorCode PCBDDCScalingReset_Deluxe_Solvers(PCBDDCDeluxeScaling deluxe_ctx) 37434a97f8cSStefano Zampini { 375*72b8c272SStefano Zampini PetscInt i; 37634a97f8cSStefano Zampini PetscErrorCode ierr; 37734a97f8cSStefano Zampini 37834a97f8cSStefano Zampini PetscFunctionBegin; 37934a97f8cSStefano Zampini ierr = PetscFree(deluxe_ctx->idx_simple_B);CHKERRQ(ierr); 38034a97f8cSStefano Zampini deluxe_ctx->n_simple = 0; 381*72b8c272SStefano Zampini for (i=0;i<deluxe_ctx->seq_n;i++) { 382*72b8c272SStefano Zampini ierr = VecScatterDestroy(&deluxe_ctx->seq_scctx[i]);CHKERRQ(ierr); 383*72b8c272SStefano Zampini ierr = VecDestroy(&deluxe_ctx->seq_work1[i]);CHKERRQ(ierr); 384*72b8c272SStefano Zampini ierr = VecDestroy(&deluxe_ctx->seq_work2[i]);CHKERRQ(ierr); 385*72b8c272SStefano Zampini ierr = MatDestroy(&deluxe_ctx->seq_mat[i]);CHKERRQ(ierr); 386*72b8c272SStefano Zampini ierr = MatDestroy(&deluxe_ctx->seq_mat_inv_sum[i]);CHKERRQ(ierr); 387*72b8c272SStefano Zampini } 388*72b8c272SStefano Zampini ierr = PetscFree5(deluxe_ctx->seq_scctx,deluxe_ctx->seq_work1,deluxe_ctx->seq_work2,deluxe_ctx->seq_mat,deluxe_ctx->seq_mat_inv_sum);CHKERRQ(ierr); 389*72b8c272SStefano Zampini ierr = PetscFree(deluxe_ctx->workspace);CHKERRQ(ierr); 390*72b8c272SStefano Zampini deluxe_ctx->seq_n = 0; 39134a97f8cSStefano Zampini PetscFunctionReturn(0); 39234a97f8cSStefano Zampini } 39334a97f8cSStefano Zampini 39434a97f8cSStefano Zampini #undef __FUNCT__ 39534a97f8cSStefano Zampini #define __FUNCT__ "PCBDDCScalingSetUp_Deluxe" 39634a97f8cSStefano Zampini static PetscErrorCode PCBDDCScalingSetUp_Deluxe(PC pc) 39734a97f8cSStefano Zampini { 39834a97f8cSStefano Zampini PC_IS *pcis=(PC_IS*)pc->data; 39934a97f8cSStefano Zampini PC_BDDC *pcbddc=(PC_BDDC*)pc->data; 40034a97f8cSStefano Zampini PCBDDCDeluxeScaling deluxe_ctx=pcbddc->deluxe_ctx; 401b96c3477SStefano Zampini PCBDDCSubSchurs sub_schurs=pcbddc->sub_schurs; 40234a97f8cSStefano Zampini PetscErrorCode ierr; 40334a97f8cSStefano Zampini 40434a97f8cSStefano Zampini PetscFunctionBegin; 405*72b8c272SStefano Zampini /* reset data structures if the topology has changed */ 406*72b8c272SStefano Zampini if (pcbddc->recompute_topography) { 40734a97f8cSStefano Zampini ierr = PCBDDCScalingReset_Deluxe_Solvers(deluxe_ctx);CHKERRQ(ierr); 408*72b8c272SStefano Zampini } 40934a97f8cSStefano Zampini 410b1b3d7a2SStefano Zampini /* Compute data structures to solve sequential problems */ 4115a95e1ceSStefano Zampini ierr = PCBDDCScalingSetUp_Deluxe_Private(pc);CHKERRQ(ierr); 4125db18549SStefano Zampini 413b1b3d7a2SStefano Zampini /* diagonal scaling on interface dofs not contained in cc */ 414d62866d3SStefano Zampini if (sub_schurs->is_vertices || sub_schurs->is_dir) { 4151b968477SStefano Zampini PetscInt n_com,n_dir; 4161b968477SStefano Zampini n_com = 0; 417d62866d3SStefano Zampini if (sub_schurs->is_vertices) { 418d62866d3SStefano Zampini ierr = ISGetLocalSize(sub_schurs->is_vertices,&n_com);CHKERRQ(ierr); 4191b968477SStefano Zampini } 4201b968477SStefano Zampini n_dir = 0; 421d62866d3SStefano Zampini if (sub_schurs->is_dir) { 422d62866d3SStefano Zampini ierr = ISGetLocalSize(sub_schurs->is_dir,&n_dir);CHKERRQ(ierr); 4231b968477SStefano Zampini } 424*72b8c272SStefano Zampini if (!deluxe_ctx->n_simple) { 4251b968477SStefano Zampini deluxe_ctx->n_simple = n_dir + n_com; 4261b968477SStefano Zampini ierr = PetscMalloc1(deluxe_ctx->n_simple,&deluxe_ctx->idx_simple_B);CHKERRQ(ierr); 427d62866d3SStefano Zampini if (sub_schurs->is_vertices) { 4289bb4a8caSStefano Zampini PetscInt nmap; 4299bb4a8caSStefano Zampini const PetscInt *idxs; 4301b968477SStefano Zampini 431d62866d3SStefano Zampini ierr = ISGetIndices(sub_schurs->is_vertices,&idxs);CHKERRQ(ierr); 4321b968477SStefano Zampini ierr = ISGlobalToLocalMappingApply(pcis->BtoNmap,IS_GTOLM_DROP,n_com,idxs,&nmap,deluxe_ctx->idx_simple_B);CHKERRQ(ierr); 4331b968477SStefano Zampini if (nmap != n_com) { 434d62866d3SStefano Zampini SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Error when mapping simply scaled dofs (is_vertices)! %d != %d",nmap,n_com); 4359bb4a8caSStefano Zampini } 436d62866d3SStefano Zampini ierr = ISRestoreIndices(sub_schurs->is_vertices,&idxs);CHKERRQ(ierr); 4371b968477SStefano Zampini } 438d62866d3SStefano Zampini if (sub_schurs->is_dir) { 4391b968477SStefano Zampini PetscInt nmap; 4401b968477SStefano Zampini const PetscInt *idxs; 4411b968477SStefano Zampini 442d62866d3SStefano Zampini ierr = ISGetIndices(sub_schurs->is_dir,&idxs);CHKERRQ(ierr); 4431b968477SStefano Zampini ierr = ISGlobalToLocalMappingApply(pcis->BtoNmap,IS_GTOLM_DROP,n_dir,idxs,&nmap,deluxe_ctx->idx_simple_B+n_com);CHKERRQ(ierr); 444d62866d3SStefano Zampini if (nmap != n_dir) { 445d62866d3SStefano Zampini SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Error when mapping simply scaled dofs (sub_schurs->is_dir)! %d != %d",nmap,n_dir); 446d62866d3SStefano Zampini } 447d62866d3SStefano Zampini ierr = ISRestoreIndices(sub_schurs->is_dir,&idxs);CHKERRQ(ierr); 4481b968477SStefano Zampini } 4491b968477SStefano Zampini ierr = PetscSortInt(deluxe_ctx->n_simple,deluxe_ctx->idx_simple_B);CHKERRQ(ierr); 4509bb4a8caSStefano Zampini } else { 451*72b8c272SStefano Zampini if (deluxe_ctx->n_simple != n_dir + n_com) { 452*72b8c272SStefano Zampini SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Number of simply scaled dofs %d is different from the previous one computed %d",n_dir + n_com,deluxe_ctx->n_simple); 453*72b8c272SStefano Zampini } 454*72b8c272SStefano Zampini } 455*72b8c272SStefano Zampini } else { 456b1b3d7a2SStefano Zampini deluxe_ctx->n_simple = 0; 4579bb4a8caSStefano Zampini deluxe_ctx->idx_simple_B = 0; 458b1b3d7a2SStefano Zampini } 45934a97f8cSStefano Zampini PetscFunctionReturn(0); 46034a97f8cSStefano Zampini } 46134a97f8cSStefano Zampini 46234a97f8cSStefano Zampini #undef __FUNCT__ 4635a95e1ceSStefano Zampini #define __FUNCT__ "PCBDDCScalingSetUp_Deluxe_Private" 4645a95e1ceSStefano Zampini static PetscErrorCode PCBDDCScalingSetUp_Deluxe_Private(PC pc) 4655db18549SStefano Zampini { 4665db18549SStefano Zampini PC_BDDC *pcbddc=(PC_BDDC*)pc->data; 4675db18549SStefano Zampini PCBDDCDeluxeScaling deluxe_ctx=pcbddc->deluxe_ctx; 468b96c3477SStefano Zampini PCBDDCSubSchurs sub_schurs = pcbddc->sub_schurs; 469*72b8c272SStefano Zampini PetscScalar *matdata,*matdata2; 470*72b8c272SStefano Zampini PetscInt i,max_subset_size,cum,cum2; 471*72b8c272SStefano Zampini const PetscInt *idxs; 472*72b8c272SStefano Zampini PetscBool newsetup = PETSC_FALSE; 4735db18549SStefano Zampini PetscErrorCode ierr; 4745db18549SStefano Zampini 4755db18549SStefano Zampini PetscFunctionBegin; 4765a95e1ceSStefano Zampini if (!sub_schurs->n_subs) { 4779221af80SStefano Zampini PetscFunctionReturn(0); 4789221af80SStefano Zampini } 4799221af80SStefano Zampini 480*72b8c272SStefano Zampini /* Allocate arrays for subproblems */ 481*72b8c272SStefano Zampini if (!deluxe_ctx->seq_n) { 482*72b8c272SStefano Zampini deluxe_ctx->seq_n = sub_schurs->n_subs; 483*72b8c272SStefano Zampini ierr = PetscCalloc5(deluxe_ctx->seq_n,&deluxe_ctx->seq_scctx,deluxe_ctx->seq_n,&deluxe_ctx->seq_work1,deluxe_ctx->seq_n,&deluxe_ctx->seq_work2,deluxe_ctx->seq_n,&deluxe_ctx->seq_mat,deluxe_ctx->seq_n,&deluxe_ctx->seq_mat_inv_sum);CHKERRQ(ierr); 484*72b8c272SStefano Zampini newsetup = PETSC_TRUE; 485*72b8c272SStefano Zampini } else if (deluxe_ctx->seq_n != sub_schurs->n_subs) { 486*72b8c272SStefano Zampini SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Number of deluxe subproblems %d is different from the sub_schurs %d",deluxe_ctx->seq_n,sub_schurs->n_subs); 487*72b8c272SStefano Zampini } 488*72b8c272SStefano Zampini deluxe_ctx->change_with_qr = sub_schurs->change_with_qr; 4895db18549SStefano Zampini 490*72b8c272SStefano Zampini /* Create objects for deluxe */ 491*72b8c272SStefano Zampini max_subset_size = 0; 492*72b8c272SStefano Zampini for (i=0;i<sub_schurs->n_subs;i++) { 493*72b8c272SStefano Zampini PetscInt subset_size; 494*72b8c272SStefano Zampini ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr); 495*72b8c272SStefano Zampini max_subset_size = PetscMax(subset_size,max_subset_size); 496*72b8c272SStefano Zampini } 497*72b8c272SStefano Zampini if (newsetup) { 498*72b8c272SStefano Zampini ierr = PetscMalloc1(2*max_subset_size,&deluxe_ctx->workspace);CHKERRQ(ierr); 499*72b8c272SStefano Zampini } 500*72b8c272SStefano Zampini cum = cum2 = 0; 501*72b8c272SStefano Zampini ierr = ISGetIndices(sub_schurs->is_Ej_all,&idxs);CHKERRQ(ierr); 502*72b8c272SStefano Zampini ierr = MatSeqAIJGetArray(sub_schurs->S_Ej_all,&matdata);CHKERRQ(ierr); 503*72b8c272SStefano Zampini matdata2 = NULL; 504*72b8c272SStefano Zampini if (sub_schurs->sum_S_Ej_all) { 505*72b8c272SStefano Zampini ierr = MatSeqAIJGetArray(sub_schurs->sum_S_Ej_all,&matdata2);CHKERRQ(ierr); 506*72b8c272SStefano Zampini } 507*72b8c272SStefano Zampini for (i=0;i<deluxe_ctx->seq_n;i++) { 508*72b8c272SStefano Zampini PetscInt subset_size; 509925dfd53SStefano Zampini 510*72b8c272SStefano Zampini ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr); 511*72b8c272SStefano Zampini if (newsetup) { 512*72b8c272SStefano Zampini IS sub; 513*72b8c272SStefano Zampini /* work vectors */ 514*72b8c272SStefano Zampini ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,1,subset_size,deluxe_ctx->workspace,&deluxe_ctx->seq_work1[i]);CHKERRQ(ierr); 515*72b8c272SStefano Zampini ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,1,subset_size,deluxe_ctx->workspace+subset_size,&deluxe_ctx->seq_work2[i]);CHKERRQ(ierr); 516*72b8c272SStefano Zampini 517*72b8c272SStefano Zampini /* scatters */ 518*72b8c272SStefano Zampini ierr = ISCreateGeneral(PETSC_COMM_SELF,subset_size,idxs+cum,PETSC_COPY_VALUES,&sub);CHKERRQ(ierr); 519*72b8c272SStefano Zampini ierr = VecScatterCreate(pcbddc->work_scaling,sub,deluxe_ctx->seq_work1[i],NULL,&deluxe_ctx->seq_scctx[i]);CHKERRQ(ierr); 520*72b8c272SStefano Zampini ierr = ISDestroy(&sub);CHKERRQ(ierr); 521*72b8c272SStefano Zampini 522*72b8c272SStefano Zampini /* S_E_j */ 523*72b8c272SStefano Zampini ierr = MatCreateSeqDense(PETSC_COMM_SELF,subset_size,subset_size,matdata+cum2,&deluxe_ctx->seq_mat[i]);CHKERRQ(ierr); 524*72b8c272SStefano Zampini } 525*72b8c272SStefano Zampini /* \sum_k S^k_E_j */ 526*72b8c272SStefano Zampini if (matdata2) { 527*72b8c272SStefano Zampini ierr = MatDestroy(&deluxe_ctx->seq_mat_inv_sum[i]);CHKERRQ(ierr); 528*72b8c272SStefano Zampini ierr = MatCreateSeqDense(PETSC_COMM_SELF,subset_size,subset_size,matdata2+cum2,&deluxe_ctx->seq_mat_inv_sum[i]);CHKERRQ(ierr); 529*72b8c272SStefano Zampini if (sub_schurs->is_hermitian && sub_schurs->is_posdef) { 530*72b8c272SStefano Zampini ierr = MatCholeskyFactor(deluxe_ctx->seq_mat_inv_sum[i],NULL,NULL);CHKERRQ(ierr); 531d62866d3SStefano Zampini } else { 532*72b8c272SStefano Zampini ierr = MatLUFactor(deluxe_ctx->seq_mat_inv_sum[i],NULL,NULL,NULL);CHKERRQ(ierr); 533d62866d3SStefano Zampini } 534*72b8c272SStefano Zampini } 535*72b8c272SStefano Zampini cum += subset_size; 536*72b8c272SStefano Zampini cum2 += subset_size*subset_size; 537*72b8c272SStefano Zampini } 538*72b8c272SStefano Zampini ierr = ISRestoreIndices(sub_schurs->is_Ej_all,&idxs);CHKERRQ(ierr); 539*72b8c272SStefano Zampini ierr = MatSeqAIJRestoreArray(sub_schurs->S_Ej_all,&matdata);CHKERRQ(ierr); 540*72b8c272SStefano Zampini if (sub_schurs->sum_S_Ej_all) { 541*72b8c272SStefano Zampini ierr = MatSeqAIJRestoreArray(sub_schurs->sum_S_Ej_all,&matdata2);CHKERRQ(ierr); 542925dfd53SStefano Zampini } 5435db18549SStefano Zampini 544*72b8c272SStefano Zampini /* the change of basis is just a reference to sub_schurs->change (if any) */ 54516e386b8SStefano Zampini deluxe_ctx->change = sub_schurs->change; 546*72b8c272SStefano Zampini if (deluxe_ctx->change && !deluxe_ctx->change_with_qr) { 547*72b8c272SStefano Zampini for (i=0;i<deluxe_ctx->seq_n;i++) { 548*72b8c272SStefano Zampini if (newsetup) { 549*72b8c272SStefano Zampini PC pc; 550*72b8c272SStefano Zampini 551*72b8c272SStefano Zampini ierr = KSPGetPC(deluxe_ctx->change[i],&pc);CHKERRQ(ierr); 552*72b8c272SStefano Zampini ierr = PCSetType(pc,PCLU);CHKERRQ(ierr); 553*72b8c272SStefano Zampini ierr = KSPSetFromOptions(deluxe_ctx->change[i]);CHKERRQ(ierr); 554*72b8c272SStefano Zampini } 555*72b8c272SStefano Zampini ierr = KSPSetUp(deluxe_ctx->change[i]);CHKERRQ(ierr); 556*72b8c272SStefano Zampini } 55716e386b8SStefano Zampini } 5585db18549SStefano Zampini PetscFunctionReturn(0); 5595db18549SStefano Zampini } 560