xref: /petsc/src/ksp/pc/impls/bddc/bddcscalingbasic.c (revision 5e5bbd0a99ea6ef6b187abb9b4a6405ccdf6db1e)
1*5e5bbd0aSStefano Zampini #include <petsc/private/pcbddcimpl.h>
2*5e5bbd0aSStefano Zampini #include <petsc/private/pcbddcprivateimpl.h>
372b8c272SStefano Zampini #include <petscblaslapack.h>
4839e9adbSstefano_zampini #include <../src/mat/impls/dense/seq/dense.h>
5674ae819SStefano Zampini 
634a97f8cSStefano Zampini /* prototypes for deluxe functions */
734a97f8cSStefano Zampini static PetscErrorCode PCBDDCScalingCreate_Deluxe(PC);
834a97f8cSStefano Zampini static PetscErrorCode PCBDDCScalingDestroy_Deluxe(PC);
934a97f8cSStefano Zampini static PetscErrorCode PCBDDCScalingSetUp_Deluxe(PC);
105a95e1ceSStefano Zampini static PetscErrorCode PCBDDCScalingSetUp_Deluxe_Private(PC);
1134a97f8cSStefano Zampini static PetscErrorCode PCBDDCScalingReset_Deluxe_Solvers(PCBDDCDeluxeScaling);
12674ae819SStefano Zampini 
13839e9adbSstefano_zampini static PetscErrorCode PCBDDCMatTransposeMatSolve_SeqDense(Mat A,Mat B,Mat X)
14839e9adbSstefano_zampini {
15839e9adbSstefano_zampini   Mat_SeqDense      *mat = (Mat_SeqDense*)A->data;
161683a169SBarry Smith   const PetscScalar *b;
171683a169SBarry Smith   PetscScalar       *x;
18839e9adbSstefano_zampini   PetscInt          n;
19839e9adbSstefano_zampini   PetscBLASInt      nrhs,info,m;
20839e9adbSstefano_zampini   PetscBool         flg;
21839e9adbSstefano_zampini 
22839e9adbSstefano_zampini   PetscFunctionBegin;
239566063dSJacob Faibussowitsch   PetscCall(PetscBLASIntCast(A->rmap->n,&m));
249566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompareAny((PetscObject)B,&flg,MATSEQDENSE,MATMPIDENSE,NULL));
2528b400f6SJacob Faibussowitsch   PetscCheck(flg,PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONG,"Matrix B must be MATDENSE matrix");
269566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompareAny((PetscObject)X,&flg,MATSEQDENSE,MATMPIDENSE,NULL));
2728b400f6SJacob Faibussowitsch   PetscCheck(flg,PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONG,"Matrix X must be MATDENSE matrix");
28839e9adbSstefano_zampini 
299566063dSJacob Faibussowitsch   PetscCall(MatGetSize(B,NULL,&n));
309566063dSJacob Faibussowitsch   PetscCall(PetscBLASIntCast(n,&nrhs));
319566063dSJacob Faibussowitsch   PetscCall(MatDenseGetArrayRead(B,&b));
329566063dSJacob Faibussowitsch   PetscCall(MatDenseGetArray(X,&x));
339566063dSJacob Faibussowitsch   PetscCall(PetscArraycpy(x,b,m*nrhs));
349566063dSJacob Faibussowitsch   PetscCall(MatDenseRestoreArrayRead(B,&b));
35839e9adbSstefano_zampini 
36839e9adbSstefano_zampini   if (A->factortype == MAT_FACTOR_LU) {
37839e9adbSstefano_zampini     PetscStackCallBLAS("LAPACKgetrs",LAPACKgetrs_("T",&m,&nrhs,mat->v,&mat->lda,mat->pivots,x,&m,&info));
3828b400f6SJacob Faibussowitsch     PetscCheck(!info,PETSC_COMM_SELF,PETSC_ERR_LIB,"GETRS - Bad solve");
393829e7e6SStefano Zampini   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Only LU factor supported");
40839e9adbSstefano_zampini 
419566063dSJacob Faibussowitsch   PetscCall(MatDenseRestoreArray(X,&x));
429566063dSJacob Faibussowitsch   PetscCall(PetscLogFlops(nrhs*(2.0*m*m - m)));
43839e9adbSstefano_zampini   PetscFunctionReturn(0);
44839e9adbSstefano_zampini }
45839e9adbSstefano_zampini 
46674ae819SStefano Zampini static PetscErrorCode PCBDDCScalingExtension_Basic(PC pc, Vec local_interface_vector, Vec global_vector)
47674ae819SStefano Zampini {
48674ae819SStefano Zampini   PC_IS*         pcis = (PC_IS*)pc->data;
49674ae819SStefano Zampini   PC_BDDC*       pcbddc = (PC_BDDC*)pc->data;
50674ae819SStefano Zampini 
51674ae819SStefano Zampini   PetscFunctionBegin;
52674ae819SStefano Zampini   /* Apply partition of unity */
539566063dSJacob Faibussowitsch   PetscCall(VecPointwiseMult(pcbddc->work_scaling,pcis->D,local_interface_vector));
549566063dSJacob Faibussowitsch   PetscCall(VecSet(global_vector,0.0));
559566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(pcis->global_to_B,pcbddc->work_scaling,global_vector,ADD_VALUES,SCATTER_REVERSE));
569566063dSJacob Faibussowitsch   PetscCall(VecScatterEnd(pcis->global_to_B,pcbddc->work_scaling,global_vector,ADD_VALUES,SCATTER_REVERSE));
57674ae819SStefano Zampini   PetscFunctionReturn(0);
58674ae819SStefano Zampini }
59674ae819SStefano Zampini 
60674ae819SStefano Zampini static PetscErrorCode PCBDDCScalingExtension_Deluxe(PC pc, Vec x, Vec y)
61674ae819SStefano Zampini {
62674ae819SStefano Zampini   PC_IS*              pcis=(PC_IS*)pc->data;
63674ae819SStefano Zampini   PC_BDDC*            pcbddc=(PC_BDDC*)pc->data;
64674ae819SStefano Zampini   PCBDDCDeluxeScaling deluxe_ctx = pcbddc->deluxe_ctx;
65674ae819SStefano Zampini 
66674ae819SStefano Zampini   PetscFunctionBegin;
679566063dSJacob Faibussowitsch   PetscCall(VecSet(pcbddc->work_scaling,0.0));
689566063dSJacob Faibussowitsch   PetscCall(VecSet(y,0.0));
695a95e1ceSStefano Zampini   if (deluxe_ctx->n_simple) { /* scale deluxe vertices using diagonal scaling */
705a95e1ceSStefano Zampini     PetscInt          i;
712b095fd8SStefano Zampini     const PetscScalar *array_x,*array_D;
722b095fd8SStefano Zampini     PetscScalar       *array;
739566063dSJacob Faibussowitsch     PetscCall(VecGetArrayRead(x,&array_x));
749566063dSJacob Faibussowitsch     PetscCall(VecGetArrayRead(pcis->D,&array_D));
759566063dSJacob Faibussowitsch     PetscCall(VecGetArray(pcbddc->work_scaling,&array));
76674ae819SStefano Zampini     for (i=0;i<deluxe_ctx->n_simple;i++) {
77674ae819SStefano Zampini       array[deluxe_ctx->idx_simple_B[i]] = array_x[deluxe_ctx->idx_simple_B[i]]*array_D[deluxe_ctx->idx_simple_B[i]];
78674ae819SStefano Zampini     }
799566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(pcbddc->work_scaling,&array));
809566063dSJacob Faibussowitsch     PetscCall(VecRestoreArrayRead(pcis->D,&array_D));
819566063dSJacob Faibussowitsch     PetscCall(VecRestoreArrayRead(x,&array_x));
8234a97f8cSStefano Zampini   }
83ac632422SStefano Zampini   /* sequential part : all problems and Schur applications collapsed into a single matrix vector multiplication or a matvec and a solve */
845a95e1ceSStefano Zampini   if (deluxe_ctx->seq_mat) {
8572b8c272SStefano Zampini     PetscInt i;
8672b8c272SStefano Zampini     for (i=0;i<deluxe_ctx->seq_n;i++) {
8716e386b8SStefano Zampini       if (deluxe_ctx->change) {
889566063dSJacob Faibussowitsch         PetscCall(VecScatterBegin(deluxe_ctx->seq_scctx[i],x,deluxe_ctx->seq_work2[i],INSERT_VALUES,SCATTER_FORWARD));
899566063dSJacob Faibussowitsch         PetscCall(VecScatterEnd(deluxe_ctx->seq_scctx[i],x,deluxe_ctx->seq_work2[i],INSERT_VALUES,SCATTER_FORWARD));
9072b8c272SStefano Zampini         if (deluxe_ctx->change_with_qr) {
9172b8c272SStefano Zampini           Mat change;
9272b8c272SStefano Zampini 
939566063dSJacob Faibussowitsch           PetscCall(KSPGetOperators(deluxe_ctx->change[i],&change,NULL));
949566063dSJacob Faibussowitsch           PetscCall(MatMultTranspose(change,deluxe_ctx->seq_work2[i],deluxe_ctx->seq_work1[i]));
9572b8c272SStefano Zampini         } else {
969566063dSJacob Faibussowitsch           PetscCall(KSPSolve(deluxe_ctx->change[i],deluxe_ctx->seq_work2[i],deluxe_ctx->seq_work1[i]));
9716e386b8SStefano Zampini         }
9872b8c272SStefano Zampini       } else {
999566063dSJacob Faibussowitsch         PetscCall(VecScatterBegin(deluxe_ctx->seq_scctx[i],x,deluxe_ctx->seq_work1[i],INSERT_VALUES,SCATTER_FORWARD));
1009566063dSJacob Faibussowitsch         PetscCall(VecScatterEnd(deluxe_ctx->seq_scctx[i],x,deluxe_ctx->seq_work1[i],INSERT_VALUES,SCATTER_FORWARD));
10172b8c272SStefano Zampini       }
1029566063dSJacob Faibussowitsch       PetscCall(MatMultTranspose(deluxe_ctx->seq_mat[i],deluxe_ctx->seq_work1[i],deluxe_ctx->seq_work2[i]));
10372b8c272SStefano Zampini       if (deluxe_ctx->seq_mat_inv_sum[i]) {
10472b8c272SStefano Zampini         PetscScalar *x;
10572b8c272SStefano Zampini 
1069566063dSJacob Faibussowitsch         PetscCall(VecGetArray(deluxe_ctx->seq_work2[i],&x));
1079566063dSJacob Faibussowitsch         PetscCall(VecPlaceArray(deluxe_ctx->seq_work1[i],x));
1089566063dSJacob Faibussowitsch         PetscCall(VecRestoreArray(deluxe_ctx->seq_work2[i],&x));
1099566063dSJacob Faibussowitsch         PetscCall(MatSolveTranspose(deluxe_ctx->seq_mat_inv_sum[i],deluxe_ctx->seq_work1[i],deluxe_ctx->seq_work2[i]));
1109566063dSJacob Faibussowitsch         PetscCall(VecResetArray(deluxe_ctx->seq_work1[i]));
111ac632422SStefano Zampini       }
11216e386b8SStefano Zampini       if (deluxe_ctx->change) {
11316e386b8SStefano Zampini         Mat change;
11416e386b8SStefano Zampini 
1159566063dSJacob Faibussowitsch         PetscCall(KSPGetOperators(deluxe_ctx->change[i],&change,NULL));
1169566063dSJacob Faibussowitsch         PetscCall(MatMult(change,deluxe_ctx->seq_work2[i],deluxe_ctx->seq_work1[i]));
1179566063dSJacob Faibussowitsch         PetscCall(VecScatterBegin(deluxe_ctx->seq_scctx[i],deluxe_ctx->seq_work1[i],pcbddc->work_scaling,INSERT_VALUES,SCATTER_REVERSE));
1189566063dSJacob Faibussowitsch         PetscCall(VecScatterEnd(deluxe_ctx->seq_scctx[i],deluxe_ctx->seq_work1[i],pcbddc->work_scaling,INSERT_VALUES,SCATTER_REVERSE));
11916e386b8SStefano Zampini       } else {
1209566063dSJacob Faibussowitsch         PetscCall(VecScatterBegin(deluxe_ctx->seq_scctx[i],deluxe_ctx->seq_work2[i],pcbddc->work_scaling,INSERT_VALUES,SCATTER_REVERSE));
1219566063dSJacob Faibussowitsch         PetscCall(VecScatterEnd(deluxe_ctx->seq_scctx[i],deluxe_ctx->seq_work2[i],pcbddc->work_scaling,INSERT_VALUES,SCATTER_REVERSE));
12272b8c272SStefano Zampini       }
123674ae819SStefano Zampini     }
12416e386b8SStefano Zampini   }
125674ae819SStefano Zampini   /* put local boundary part in global vector */
1269566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(pcis->global_to_B,pcbddc->work_scaling,y,ADD_VALUES,SCATTER_REVERSE));
1279566063dSJacob Faibussowitsch   PetscCall(VecScatterEnd(pcis->global_to_B,pcbddc->work_scaling,y,ADD_VALUES,SCATTER_REVERSE));
128674ae819SStefano Zampini   PetscFunctionReturn(0);
129674ae819SStefano Zampini }
130674ae819SStefano Zampini 
131674ae819SStefano Zampini PetscErrorCode PCBDDCScalingExtension(PC pc, Vec local_interface_vector, Vec global_vector)
132674ae819SStefano Zampini {
133674ae819SStefano Zampini   PC_BDDC        *pcbddc=(PC_BDDC*)pc->data;
134674ae819SStefano Zampini 
135674ae819SStefano Zampini   PetscFunctionBegin;
136674ae819SStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
137674ae819SStefano Zampini   PetscValidHeaderSpecific(local_interface_vector,VEC_CLASSID,2);
138674ae819SStefano Zampini   PetscValidHeaderSpecific(global_vector,VEC_CLASSID,3);
13908401ef6SPierre Jolivet   PetscCheck(local_interface_vector != pcbddc->work_scaling,PETSC_COMM_SELF,PETSC_ERR_SUP,"Local vector cannot be pcbddc->work_scaling!");
140cac4c232SBarry Smith   PetscUseMethod(pc,"PCBDDCScalingExtension_C",(PC,Vec,Vec),(pc,local_interface_vector,global_vector));
141674ae819SStefano Zampini   PetscFunctionReturn(0);
142674ae819SStefano Zampini }
143674ae819SStefano Zampini 
144674ae819SStefano Zampini static PetscErrorCode PCBDDCScalingRestriction_Basic(PC pc, Vec global_vector, Vec local_interface_vector)
145674ae819SStefano Zampini {
146674ae819SStefano Zampini   PC_IS          *pcis = (PC_IS*)pc->data;
147674ae819SStefano Zampini 
148674ae819SStefano Zampini   PetscFunctionBegin;
1499566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(pcis->global_to_B,global_vector,local_interface_vector,INSERT_VALUES,SCATTER_FORWARD));
1509566063dSJacob Faibussowitsch   PetscCall(VecScatterEnd(pcis->global_to_B,global_vector,local_interface_vector,INSERT_VALUES,SCATTER_FORWARD));
151674ae819SStefano Zampini   /* Apply partition of unity */
1529566063dSJacob Faibussowitsch   PetscCall(VecPointwiseMult(local_interface_vector,pcis->D,local_interface_vector));
153674ae819SStefano Zampini   PetscFunctionReturn(0);
154674ae819SStefano Zampini }
155674ae819SStefano Zampini 
156674ae819SStefano Zampini static PetscErrorCode PCBDDCScalingRestriction_Deluxe(PC pc, Vec x, Vec y)
157674ae819SStefano Zampini {
158674ae819SStefano Zampini   PC_IS*              pcis=(PC_IS*)pc->data;
159674ae819SStefano Zampini   PC_BDDC*            pcbddc=(PC_BDDC*)pc->data;
160674ae819SStefano Zampini   PCBDDCDeluxeScaling deluxe_ctx = pcbddc->deluxe_ctx;
161674ae819SStefano Zampini 
162674ae819SStefano Zampini   PetscFunctionBegin;
163674ae819SStefano Zampini   /* get local boundary part of global vector */
1649566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(pcis->global_to_B,x,y,INSERT_VALUES,SCATTER_FORWARD));
1659566063dSJacob Faibussowitsch   PetscCall(VecScatterEnd(pcis->global_to_B,x,y,INSERT_VALUES,SCATTER_FORWARD));
1665a95e1ceSStefano Zampini   if (deluxe_ctx->n_simple) { /* scale deluxe vertices using diagonal scaling */
1675a95e1ceSStefano Zampini     PetscInt          i;
1682b095fd8SStefano Zampini     PetscScalar       *array_y;
1692b095fd8SStefano Zampini     const PetscScalar *array_D;
1709566063dSJacob Faibussowitsch     PetscCall(VecGetArray(y,&array_y));
1719566063dSJacob Faibussowitsch     PetscCall(VecGetArrayRead(pcis->D,&array_D));
172674ae819SStefano Zampini     for (i=0;i<deluxe_ctx->n_simple;i++) {
173674ae819SStefano Zampini       array_y[deluxe_ctx->idx_simple_B[i]] *= array_D[deluxe_ctx->idx_simple_B[i]];
174674ae819SStefano Zampini     }
1759566063dSJacob Faibussowitsch     PetscCall(VecRestoreArrayRead(pcis->D,&array_D));
1769566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(y,&array_y));
17734a97f8cSStefano Zampini   }
178729c86d3Sstefano_zampini   /* sequential part : all problems and Schur applications collapsed into a single matrix vector multiplication or a matvec and a solve */
1795a95e1ceSStefano Zampini   if (deluxe_ctx->seq_mat) {
18072b8c272SStefano Zampini     PetscInt i;
18172b8c272SStefano Zampini     for (i=0;i<deluxe_ctx->seq_n;i++) {
18216e386b8SStefano Zampini       if (deluxe_ctx->change) {
18316e386b8SStefano Zampini         Mat change;
18416e386b8SStefano Zampini 
1859566063dSJacob Faibussowitsch         PetscCall(VecScatterBegin(deluxe_ctx->seq_scctx[i],y,deluxe_ctx->seq_work2[i],INSERT_VALUES,SCATTER_FORWARD));
1869566063dSJacob Faibussowitsch         PetscCall(VecScatterEnd(deluxe_ctx->seq_scctx[i],y,deluxe_ctx->seq_work2[i],INSERT_VALUES,SCATTER_FORWARD));
1879566063dSJacob Faibussowitsch         PetscCall(KSPGetOperators(deluxe_ctx->change[i],&change,NULL));
1889566063dSJacob Faibussowitsch         PetscCall(MatMultTranspose(change,deluxe_ctx->seq_work2[i],deluxe_ctx->seq_work1[i]));
18916e386b8SStefano Zampini       } else {
1909566063dSJacob Faibussowitsch         PetscCall(VecScatterBegin(deluxe_ctx->seq_scctx[i],y,deluxe_ctx->seq_work1[i],INSERT_VALUES,SCATTER_FORWARD));
1919566063dSJacob Faibussowitsch         PetscCall(VecScatterEnd(deluxe_ctx->seq_scctx[i],y,deluxe_ctx->seq_work1[i],INSERT_VALUES,SCATTER_FORWARD));
19216e386b8SStefano Zampini       }
19372b8c272SStefano Zampini       if (deluxe_ctx->seq_mat_inv_sum[i]) {
19472b8c272SStefano Zampini         PetscScalar *x;
19572b8c272SStefano Zampini 
1969566063dSJacob Faibussowitsch         PetscCall(VecGetArray(deluxe_ctx->seq_work1[i],&x));
1979566063dSJacob Faibussowitsch         PetscCall(VecPlaceArray(deluxe_ctx->seq_work2[i],x));
1989566063dSJacob Faibussowitsch         PetscCall(VecRestoreArray(deluxe_ctx->seq_work1[i],&x));
1999566063dSJacob Faibussowitsch         PetscCall(MatSolve(deluxe_ctx->seq_mat_inv_sum[i],deluxe_ctx->seq_work1[i],deluxe_ctx->seq_work2[i]));
2009566063dSJacob Faibussowitsch         PetscCall(VecResetArray(deluxe_ctx->seq_work2[i]));
201ac632422SStefano Zampini       }
2029566063dSJacob Faibussowitsch       PetscCall(MatMult(deluxe_ctx->seq_mat[i],deluxe_ctx->seq_work1[i],deluxe_ctx->seq_work2[i]));
20316e386b8SStefano Zampini       if (deluxe_ctx->change) {
20472b8c272SStefano Zampini         if (deluxe_ctx->change_with_qr) {
20572b8c272SStefano Zampini           Mat change;
20672b8c272SStefano Zampini 
2079566063dSJacob Faibussowitsch           PetscCall(KSPGetOperators(deluxe_ctx->change[i],&change,NULL));
2089566063dSJacob Faibussowitsch           PetscCall(MatMult(change,deluxe_ctx->seq_work2[i],deluxe_ctx->seq_work1[i]));
20972b8c272SStefano Zampini         } else {
2109566063dSJacob Faibussowitsch           PetscCall(KSPSolveTranspose(deluxe_ctx->change[i],deluxe_ctx->seq_work2[i],deluxe_ctx->seq_work1[i]));
21116e386b8SStefano Zampini         }
2129566063dSJacob Faibussowitsch         PetscCall(VecScatterBegin(deluxe_ctx->seq_scctx[i],deluxe_ctx->seq_work1[i],y,INSERT_VALUES,SCATTER_REVERSE));
2139566063dSJacob Faibussowitsch         PetscCall(VecScatterEnd(deluxe_ctx->seq_scctx[i],deluxe_ctx->seq_work1[i],y,INSERT_VALUES,SCATTER_REVERSE));
21472b8c272SStefano Zampini       } else {
2159566063dSJacob Faibussowitsch         PetscCall(VecScatterBegin(deluxe_ctx->seq_scctx[i],deluxe_ctx->seq_work2[i],y,INSERT_VALUES,SCATTER_REVERSE));
2169566063dSJacob Faibussowitsch         PetscCall(VecScatterEnd(deluxe_ctx->seq_scctx[i],deluxe_ctx->seq_work2[i],y,INSERT_VALUES,SCATTER_REVERSE));
21772b8c272SStefano Zampini       }
21872b8c272SStefano Zampini     }
219674ae819SStefano Zampini   }
220674ae819SStefano Zampini   PetscFunctionReturn(0);
221674ae819SStefano Zampini }
222674ae819SStefano Zampini 
223674ae819SStefano Zampini PetscErrorCode PCBDDCScalingRestriction(PC pc, Vec global_vector, Vec local_interface_vector)
224674ae819SStefano Zampini {
225674ae819SStefano Zampini   PC_BDDC        *pcbddc=(PC_BDDC*)pc->data;
226674ae819SStefano Zampini 
227674ae819SStefano Zampini   PetscFunctionBegin;
228674ae819SStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
229674ae819SStefano Zampini   PetscValidHeaderSpecific(global_vector,VEC_CLASSID,2);
230674ae819SStefano Zampini   PetscValidHeaderSpecific(local_interface_vector,VEC_CLASSID,3);
23108401ef6SPierre Jolivet   PetscCheck(local_interface_vector != pcbddc->work_scaling,PETSC_COMM_SELF,PETSC_ERR_SUP,"Local vector cannot be pcbddc->work_scaling!");
232cac4c232SBarry Smith   PetscUseMethod(pc,"PCBDDCScalingRestriction_C",(PC,Vec,Vec),(pc,global_vector,local_interface_vector));
233674ae819SStefano Zampini   PetscFunctionReturn(0);
234674ae819SStefano Zampini }
235674ae819SStefano Zampini 
236674ae819SStefano Zampini PetscErrorCode PCBDDCScalingSetUp(PC pc)
237674ae819SStefano Zampini {
238674ae819SStefano Zampini   PC_IS*         pcis=(PC_IS*)pc->data;
239674ae819SStefano Zampini   PC_BDDC*       pcbddc=(PC_BDDC*)pc->data;
240674ae819SStefano Zampini 
241674ae819SStefano Zampini   PetscFunctionBegin;
242674ae819SStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
2439566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(PC_BDDC_Scaling[pcbddc->current_level],pc,0,0,0));
244674ae819SStefano Zampini   /* create work vector for the operator */
2459566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&pcbddc->work_scaling));
2469566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(pcis->vec1_B,&pcbddc->work_scaling));
24734a97f8cSStefano Zampini   /* always rebuild pcis->D */
24828d874f6SStefano Zampini   if (pcis->use_stiffness_scaling) {
24915f235b8SStefano Zampini     PetscScalar *a;
25015f235b8SStefano Zampini     PetscInt    i,n;
25115f235b8SStefano Zampini 
2529566063dSJacob Faibussowitsch     PetscCall(MatGetDiagonal(pcbddc->local_mat,pcis->vec1_N));
2539566063dSJacob Faibussowitsch     PetscCall(VecScatterBegin(pcis->N_to_B,pcis->vec1_N,pcis->D,INSERT_VALUES,SCATTER_FORWARD));
2549566063dSJacob Faibussowitsch     PetscCall(VecScatterEnd(pcis->N_to_B,pcis->vec1_N,pcis->D,INSERT_VALUES,SCATTER_FORWARD));
2559566063dSJacob Faibussowitsch     PetscCall(VecAbs(pcis->D));
2569566063dSJacob Faibussowitsch     PetscCall(VecGetLocalSize(pcis->D,&n));
2579566063dSJacob Faibussowitsch     PetscCall(VecGetArray(pcis->D,&a));
25815f235b8SStefano Zampini     for (i=0;i<n;i++) if (PetscAbsScalar(a[i])<PETSC_SMALL) a[i] = 1.0;
2599566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(pcis->D,&a));
260674ae819SStefano Zampini   }
2619566063dSJacob Faibussowitsch   PetscCall(VecSet(pcis->vec1_global,0.0));
2629566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(pcis->global_to_B,pcis->D,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE));
2639566063dSJacob Faibussowitsch   PetscCall(VecScatterEnd(pcis->global_to_B,pcis->D,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE));
2649566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(pcis->global_to_B,pcis->vec1_global,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD));
2659566063dSJacob Faibussowitsch   PetscCall(VecScatterEnd(pcis->global_to_B,pcis->vec1_global,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD));
2669566063dSJacob Faibussowitsch   PetscCall(VecPointwiseDivide(pcis->D,pcis->D,pcis->vec1_B));
267674ae819SStefano Zampini   /* now setup */
268681e7c04SStefano Zampini   if (pcbddc->use_deluxe_scaling) {
26934a97f8cSStefano Zampini     if (!pcbddc->deluxe_ctx) {
2709566063dSJacob Faibussowitsch       PetscCall(PCBDDCScalingCreate_Deluxe(pc));
27134a97f8cSStefano Zampini     }
2729566063dSJacob Faibussowitsch     PetscCall(PCBDDCScalingSetUp_Deluxe(pc));
2739566063dSJacob Faibussowitsch     PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCBDDCScalingRestriction_C",PCBDDCScalingRestriction_Deluxe));
2749566063dSJacob Faibussowitsch     PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCBDDCScalingExtension_C",PCBDDCScalingExtension_Deluxe));
275674ae819SStefano Zampini   } else {
2769566063dSJacob Faibussowitsch     PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCBDDCScalingRestriction_C",PCBDDCScalingRestriction_Basic));
2779566063dSJacob Faibussowitsch     PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCBDDCScalingExtension_C",PCBDDCScalingExtension_Basic));
278674ae819SStefano Zampini   }
2799566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(PC_BDDC_Scaling[pcbddc->current_level],pc,0,0,0));
28034a97f8cSStefano Zampini 
281674ae819SStefano Zampini   /* test */
282674ae819SStefano Zampini   if (pcbddc->dbg_flag) {
28372b8c272SStefano Zampini     Mat         B0_B = NULL;
28472b8c272SStefano Zampini     Vec         B0_Bv = NULL, B0_Bv2 = NULL;
28534a97f8cSStefano Zampini     Vec         vec2_global;
286674ae819SStefano Zampini     PetscViewer viewer = pcbddc->dbg_viewer;
28734a97f8cSStefano Zampini     PetscReal   error;
288674ae819SStefano Zampini 
289674ae819SStefano Zampini     /* extension -> from local to parallel */
2909566063dSJacob Faibussowitsch     PetscCall(VecSet(pcis->vec1_global,0.0));
2919566063dSJacob Faibussowitsch     PetscCall(VecSetRandom(pcis->vec1_B,NULL));
2929566063dSJacob Faibussowitsch     PetscCall(VecScatterBegin(pcis->global_to_B,pcis->vec1_B,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE));
2939566063dSJacob Faibussowitsch     PetscCall(VecScatterEnd(pcis->global_to_B,pcis->vec1_B,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE));
2949566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(pcis->vec1_global,&vec2_global));
2959566063dSJacob Faibussowitsch     PetscCall(VecCopy(pcis->vec1_global,vec2_global));
2969566063dSJacob Faibussowitsch     PetscCall(VecScatterBegin(pcis->global_to_B,pcis->vec1_global,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD));
2979566063dSJacob Faibussowitsch     PetscCall(VecScatterEnd(pcis->global_to_B,pcis->vec1_global,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD));
29872b8c272SStefano Zampini     if (pcbddc->benign_n) {
29972b8c272SStefano Zampini       IS is_dummy;
30072b8c272SStefano Zampini 
3019566063dSJacob Faibussowitsch       PetscCall(ISCreateStride(PETSC_COMM_SELF,pcbddc->benign_n,0,1,&is_dummy));
3029566063dSJacob Faibussowitsch       PetscCall(MatCreateSubMatrix(pcbddc->benign_B0,is_dummy,pcis->is_B_local,MAT_INITIAL_MATRIX,&B0_B));
3039566063dSJacob Faibussowitsch       PetscCall(ISDestroy(&is_dummy));
3049566063dSJacob Faibussowitsch       PetscCall(MatCreateVecs(B0_B,NULL,&B0_Bv));
3059566063dSJacob Faibussowitsch       PetscCall(VecDuplicate(B0_Bv,&B0_Bv2));
3069566063dSJacob Faibussowitsch       PetscCall(MatMult(B0_B,pcis->vec1_B,B0_Bv));
30772b8c272SStefano Zampini     }
3089566063dSJacob Faibussowitsch     PetscCall(PCBDDCScalingExtension(pc,pcis->vec1_B,pcis->vec1_global));
30972b8c272SStefano Zampini     if (pcbddc->benign_saddle_point) {
31072b8c272SStefano Zampini       PetscReal errorl = 0.;
3119566063dSJacob Faibussowitsch       PetscCall(VecScatterBegin(pcis->global_to_B,pcis->vec1_global,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD));
3129566063dSJacob Faibussowitsch       PetscCall(VecScatterEnd(pcis->global_to_B,pcis->vec1_global,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD));
31372b8c272SStefano Zampini       if (pcbddc->benign_n) {
3149566063dSJacob Faibussowitsch         PetscCall(MatMult(B0_B,pcis->vec1_B,B0_Bv2));
3159566063dSJacob Faibussowitsch         PetscCall(VecAXPY(B0_Bv,-1.0,B0_Bv2));
3169566063dSJacob Faibussowitsch         PetscCall(VecNorm(B0_Bv,NORM_INFINITY,&errorl));
31772b8c272SStefano Zampini       }
3189566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Allreduce(&errorl,&error,1,MPIU_REAL,MPI_SUM,PetscObjectComm((PetscObject)pc)));
31963a3b9bcSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer,"Error benign extension %1.14e\n",(double)error));
32072b8c272SStefano Zampini     }
3219566063dSJacob Faibussowitsch     PetscCall(VecAXPY(pcis->vec1_global,-1.0,vec2_global));
3229566063dSJacob Faibussowitsch     PetscCall(VecNorm(pcis->vec1_global,NORM_INFINITY,&error));
32363a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer,"Error scaling extension %1.14e\n",(double)error));
3249566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&vec2_global));
32534a97f8cSStefano Zampini 
326674ae819SStefano Zampini     /* restriction -> from parallel to local */
3279566063dSJacob Faibussowitsch     PetscCall(VecSet(pcis->vec1_global,0.0));
3289566063dSJacob Faibussowitsch     PetscCall(VecSetRandom(pcis->vec1_B,NULL));
3299566063dSJacob Faibussowitsch     PetscCall(VecScatterBegin(pcis->global_to_B,pcis->vec1_B,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE));
3309566063dSJacob Faibussowitsch     PetscCall(VecScatterEnd(pcis->global_to_B,pcis->vec1_B,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE));
3319566063dSJacob Faibussowitsch     PetscCall(PCBDDCScalingRestriction(pc,pcis->vec1_global,pcis->vec1_B));
3329566063dSJacob Faibussowitsch     PetscCall(VecScale(pcis->vec1_B,-1.0));
3339566063dSJacob Faibussowitsch     PetscCall(VecScatterBegin(pcis->global_to_B,pcis->vec1_B,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE));
3349566063dSJacob Faibussowitsch     PetscCall(VecScatterEnd(pcis->global_to_B,pcis->vec1_B,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE));
3359566063dSJacob Faibussowitsch     PetscCall(VecNorm(pcis->vec1_global,NORM_INFINITY,&error));
33663a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer,"Error scaling restriction %1.14e\n",(double)error));
3379566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&B0_B));
3389566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&B0_Bv));
3399566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&B0_Bv2));
340674ae819SStefano Zampini   }
341674ae819SStefano Zampini   PetscFunctionReturn(0);
342674ae819SStefano Zampini }
343674ae819SStefano Zampini 
344674ae819SStefano Zampini PetscErrorCode PCBDDCScalingDestroy(PC pc)
345674ae819SStefano Zampini {
346674ae819SStefano Zampini   PC_BDDC*       pcbddc=(PC_BDDC*)pc->data;
347674ae819SStefano Zampini 
348674ae819SStefano Zampini   PetscFunctionBegin;
34934a97f8cSStefano Zampini   if (pcbddc->deluxe_ctx) {
3509566063dSJacob Faibussowitsch     PetscCall(PCBDDCScalingDestroy_Deluxe(pc));
351674ae819SStefano Zampini   }
3529566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&pcbddc->work_scaling));
3539566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCBDDCScalingRestriction_C",NULL));
3549566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCBDDCScalingExtension_C",NULL));
355674ae819SStefano Zampini   PetscFunctionReturn(0);
356674ae819SStefano Zampini }
357674ae819SStefano Zampini 
35834a97f8cSStefano Zampini static PetscErrorCode PCBDDCScalingCreate_Deluxe(PC pc)
35934a97f8cSStefano Zampini {
36034a97f8cSStefano Zampini   PC_BDDC*            pcbddc=(PC_BDDC*)pc->data;
36134a97f8cSStefano Zampini   PCBDDCDeluxeScaling deluxe_ctx;
36234a97f8cSStefano Zampini 
36334a97f8cSStefano Zampini   PetscFunctionBegin;
3649566063dSJacob Faibussowitsch   PetscCall(PetscNew(&deluxe_ctx));
36534a97f8cSStefano Zampini   pcbddc->deluxe_ctx = deluxe_ctx;
36634a97f8cSStefano Zampini   PetscFunctionReturn(0);
36734a97f8cSStefano Zampini }
36834a97f8cSStefano Zampini 
36934a97f8cSStefano Zampini static PetscErrorCode PCBDDCScalingDestroy_Deluxe(PC pc)
37034a97f8cSStefano Zampini {
37134a97f8cSStefano Zampini   PC_BDDC*            pcbddc=(PC_BDDC*)pc->data;
37234a97f8cSStefano Zampini 
37334a97f8cSStefano Zampini   PetscFunctionBegin;
3749566063dSJacob Faibussowitsch   PetscCall(PCBDDCScalingReset_Deluxe_Solvers(pcbddc->deluxe_ctx));
3759566063dSJacob Faibussowitsch   PetscCall(PetscFree(pcbddc->deluxe_ctx));
37634a97f8cSStefano Zampini   PetscFunctionReturn(0);
37734a97f8cSStefano Zampini }
37834a97f8cSStefano Zampini 
37934a97f8cSStefano Zampini static PetscErrorCode PCBDDCScalingReset_Deluxe_Solvers(PCBDDCDeluxeScaling deluxe_ctx)
38034a97f8cSStefano Zampini {
38172b8c272SStefano Zampini   PetscInt       i;
38234a97f8cSStefano Zampini 
38334a97f8cSStefano Zampini   PetscFunctionBegin;
3849566063dSJacob Faibussowitsch   PetscCall(PetscFree(deluxe_ctx->idx_simple_B));
38534a97f8cSStefano Zampini   deluxe_ctx->n_simple = 0;
38672b8c272SStefano Zampini   for (i=0;i<deluxe_ctx->seq_n;i++) {
3879566063dSJacob Faibussowitsch     PetscCall(VecScatterDestroy(&deluxe_ctx->seq_scctx[i]));
3889566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&deluxe_ctx->seq_work1[i]));
3899566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&deluxe_ctx->seq_work2[i]));
3909566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&deluxe_ctx->seq_mat[i]));
3919566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&deluxe_ctx->seq_mat_inv_sum[i]));
39272b8c272SStefano Zampini   }
3939566063dSJacob Faibussowitsch   PetscCall(PetscFree5(deluxe_ctx->seq_scctx,deluxe_ctx->seq_work1,deluxe_ctx->seq_work2,deluxe_ctx->seq_mat,deluxe_ctx->seq_mat_inv_sum));
3949566063dSJacob Faibussowitsch   PetscCall(PetscFree(deluxe_ctx->workspace));
39572b8c272SStefano Zampini   deluxe_ctx->seq_n = 0;
39634a97f8cSStefano Zampini   PetscFunctionReturn(0);
39734a97f8cSStefano Zampini }
39834a97f8cSStefano Zampini 
39934a97f8cSStefano Zampini static PetscErrorCode PCBDDCScalingSetUp_Deluxe(PC pc)
40034a97f8cSStefano Zampini {
40134a97f8cSStefano Zampini   PC_IS               *pcis=(PC_IS*)pc->data;
40234a97f8cSStefano Zampini   PC_BDDC             *pcbddc=(PC_BDDC*)pc->data;
40334a97f8cSStefano Zampini   PCBDDCDeluxeScaling deluxe_ctx=pcbddc->deluxe_ctx;
404b96c3477SStefano Zampini   PCBDDCSubSchurs     sub_schurs=pcbddc->sub_schurs;
40534a97f8cSStefano Zampini 
40634a97f8cSStefano Zampini   PetscFunctionBegin;
40772b8c272SStefano Zampini   /* reset data structures if the topology has changed */
40872b8c272SStefano Zampini   if (pcbddc->recompute_topography) {
4099566063dSJacob Faibussowitsch     PetscCall(PCBDDCScalingReset_Deluxe_Solvers(deluxe_ctx));
41072b8c272SStefano Zampini   }
41134a97f8cSStefano Zampini 
412b1b3d7a2SStefano Zampini   /* Compute data structures to solve sequential problems */
4139566063dSJacob Faibussowitsch   PetscCall(PCBDDCScalingSetUp_Deluxe_Private(pc));
4145db18549SStefano Zampini 
415b1b3d7a2SStefano Zampini   /* diagonal scaling on interface dofs not contained in cc */
416d62866d3SStefano Zampini   if (sub_schurs->is_vertices || sub_schurs->is_dir) {
4171b968477SStefano Zampini     PetscInt n_com,n_dir;
4181b968477SStefano Zampini     n_com = 0;
419d62866d3SStefano Zampini     if (sub_schurs->is_vertices) {
4209566063dSJacob Faibussowitsch       PetscCall(ISGetLocalSize(sub_schurs->is_vertices,&n_com));
4211b968477SStefano Zampini     }
4221b968477SStefano Zampini     n_dir = 0;
423d62866d3SStefano Zampini     if (sub_schurs->is_dir) {
4249566063dSJacob Faibussowitsch       PetscCall(ISGetLocalSize(sub_schurs->is_dir,&n_dir));
4251b968477SStefano Zampini     }
42672b8c272SStefano Zampini     if (!deluxe_ctx->n_simple) {
4271b968477SStefano Zampini       deluxe_ctx->n_simple = n_dir + n_com;
4289566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(deluxe_ctx->n_simple,&deluxe_ctx->idx_simple_B));
429d62866d3SStefano Zampini       if (sub_schurs->is_vertices) {
4309bb4a8caSStefano Zampini         PetscInt       nmap;
4319bb4a8caSStefano Zampini         const PetscInt *idxs;
4321b968477SStefano Zampini 
4339566063dSJacob Faibussowitsch         PetscCall(ISGetIndices(sub_schurs->is_vertices,&idxs));
4349566063dSJacob Faibussowitsch         PetscCall(ISGlobalToLocalMappingApply(pcis->BtoNmap,IS_GTOLM_DROP,n_com,idxs,&nmap,deluxe_ctx->idx_simple_B));
43563a3b9bcSJacob Faibussowitsch         PetscCheck(nmap == n_com,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Error when mapping simply scaled dofs (is_vertices)! %" PetscInt_FMT " != %" PetscInt_FMT,nmap,n_com);
4369566063dSJacob Faibussowitsch         PetscCall(ISRestoreIndices(sub_schurs->is_vertices,&idxs));
4371b968477SStefano Zampini       }
438d62866d3SStefano Zampini       if (sub_schurs->is_dir) {
4391b968477SStefano Zampini         PetscInt       nmap;
4401b968477SStefano Zampini         const PetscInt *idxs;
4411b968477SStefano Zampini 
4429566063dSJacob Faibussowitsch         PetscCall(ISGetIndices(sub_schurs->is_dir,&idxs));
4439566063dSJacob Faibussowitsch         PetscCall(ISGlobalToLocalMappingApply(pcis->BtoNmap,IS_GTOLM_DROP,n_dir,idxs,&nmap,deluxe_ctx->idx_simple_B+n_com));
44463a3b9bcSJacob Faibussowitsch         PetscCheck(nmap == n_dir,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Error when mapping simply scaled dofs (sub_schurs->is_dir)! %" PetscInt_FMT " != %" PetscInt_FMT,nmap,n_dir);
4459566063dSJacob Faibussowitsch         PetscCall(ISRestoreIndices(sub_schurs->is_dir,&idxs));
4461b968477SStefano Zampini       }
4479566063dSJacob Faibussowitsch       PetscCall(PetscSortInt(deluxe_ctx->n_simple,deluxe_ctx->idx_simple_B));
4489bb4a8caSStefano Zampini     } else {
4492472a847SBarry Smith       PetscCheck(deluxe_ctx->n_simple == n_dir + n_com,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Number of simply scaled dofs %" PetscInt_FMT " is different from the previous one computed %" PetscInt_FMT,n_dir + n_com,deluxe_ctx->n_simple);
45072b8c272SStefano Zampini     }
45172b8c272SStefano Zampini   } else {
452b1b3d7a2SStefano Zampini     deluxe_ctx->n_simple = 0;
4530a545947SLisandro Dalcin     deluxe_ctx->idx_simple_B = NULL;
454b1b3d7a2SStefano Zampini   }
45534a97f8cSStefano Zampini   PetscFunctionReturn(0);
45634a97f8cSStefano Zampini }
45734a97f8cSStefano Zampini 
4585a95e1ceSStefano Zampini static PetscErrorCode PCBDDCScalingSetUp_Deluxe_Private(PC pc)
4595db18549SStefano Zampini {
4605db18549SStefano Zampini   PC_BDDC                *pcbddc=(PC_BDDC*)pc->data;
4615db18549SStefano Zampini   PCBDDCDeluxeScaling    deluxe_ctx=pcbddc->deluxe_ctx;
462b96c3477SStefano Zampini   PCBDDCSubSchurs        sub_schurs = pcbddc->sub_schurs;
46372b8c272SStefano Zampini   PetscScalar            *matdata,*matdata2;
46472b8c272SStefano Zampini   PetscInt               i,max_subset_size,cum,cum2;
46572b8c272SStefano Zampini   const PetscInt         *idxs;
46672b8c272SStefano Zampini   PetscBool              newsetup = PETSC_FALSE;
4675db18549SStefano Zampini 
4685db18549SStefano Zampini   PetscFunctionBegin;
46928b400f6SJacob Faibussowitsch   PetscCheck(sub_schurs,PetscObjectComm((PetscObject)pc),PETSC_ERR_PLIB,"Missing PCBDDCSubSchurs");
47086bfa4cfSStefano Zampini   if (!sub_schurs->n_subs) PetscFunctionReturn(0);
4719221af80SStefano Zampini 
47272b8c272SStefano Zampini   /* Allocate arrays for subproblems */
47372b8c272SStefano Zampini   if (!deluxe_ctx->seq_n) {
47472b8c272SStefano Zampini     deluxe_ctx->seq_n = sub_schurs->n_subs;
4759566063dSJacob Faibussowitsch     PetscCall(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));
47672b8c272SStefano Zampini     newsetup = PETSC_TRUE;
47763a3b9bcSJacob Faibussowitsch   } else PetscCheck(deluxe_ctx->seq_n == sub_schurs->n_subs,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Number of deluxe subproblems %" PetscInt_FMT " is different from the sub_schurs %" PetscInt_FMT,deluxe_ctx->seq_n,sub_schurs->n_subs);
4786080607fSStefano Zampini 
4792f41f7d1SStefano Zampini   /* the change of basis is just a reference to sub_schurs->change (if any) */
4802f41f7d1SStefano Zampini   deluxe_ctx->change         = sub_schurs->change;
48172b8c272SStefano Zampini   deluxe_ctx->change_with_qr = sub_schurs->change_with_qr;
4825db18549SStefano Zampini 
48372b8c272SStefano Zampini   /* Create objects for deluxe */
48472b8c272SStefano Zampini   max_subset_size = 0;
48572b8c272SStefano Zampini   for (i=0;i<sub_schurs->n_subs;i++) {
48672b8c272SStefano Zampini     PetscInt subset_size;
4879566063dSJacob Faibussowitsch     PetscCall(ISGetLocalSize(sub_schurs->is_subs[i],&subset_size));
48872b8c272SStefano Zampini     max_subset_size = PetscMax(subset_size,max_subset_size);
48972b8c272SStefano Zampini   }
49072b8c272SStefano Zampini   if (newsetup) {
4919566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(2*max_subset_size,&deluxe_ctx->workspace));
49272b8c272SStefano Zampini   }
49372b8c272SStefano Zampini   cum = cum2 = 0;
4949566063dSJacob Faibussowitsch   PetscCall(ISGetIndices(sub_schurs->is_Ej_all,&idxs));
4959566063dSJacob Faibussowitsch   PetscCall(MatSeqAIJGetArray(sub_schurs->S_Ej_all,&matdata));
4969566063dSJacob Faibussowitsch   PetscCall(MatSeqAIJGetArray(sub_schurs->sum_S_Ej_all,&matdata2));
49772b8c272SStefano Zampini   for (i=0;i<deluxe_ctx->seq_n;i++) {
49872b8c272SStefano Zampini     PetscInt     subset_size;
499925dfd53SStefano Zampini 
5009566063dSJacob Faibussowitsch     PetscCall(ISGetLocalSize(sub_schurs->is_subs[i],&subset_size));
50172b8c272SStefano Zampini     if (newsetup) {
50272b8c272SStefano Zampini       IS  sub;
50372b8c272SStefano Zampini       /* work vectors */
5049566063dSJacob Faibussowitsch       PetscCall(VecCreateSeqWithArray(PETSC_COMM_SELF,1,subset_size,deluxe_ctx->workspace,&deluxe_ctx->seq_work1[i]));
5059566063dSJacob Faibussowitsch       PetscCall(VecCreateSeqWithArray(PETSC_COMM_SELF,1,subset_size,deluxe_ctx->workspace+subset_size,&deluxe_ctx->seq_work2[i]));
50672b8c272SStefano Zampini 
50772b8c272SStefano Zampini       /* scatters */
5089566063dSJacob Faibussowitsch       PetscCall(ISCreateGeneral(PETSC_COMM_SELF,subset_size,idxs+cum,PETSC_COPY_VALUES,&sub));
5099566063dSJacob Faibussowitsch       PetscCall(VecScatterCreate(pcbddc->work_scaling,sub,deluxe_ctx->seq_work1[i],NULL,&deluxe_ctx->seq_scctx[i]));
5109566063dSJacob Faibussowitsch       PetscCall(ISDestroy(&sub));
511839e9adbSstefano_zampini     }
51272b8c272SStefano Zampini 
51372b8c272SStefano Zampini     /* S_E_j */
5149566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&deluxe_ctx->seq_mat[i]));
5159566063dSJacob Faibussowitsch     PetscCall(MatCreateSeqDense(PETSC_COMM_SELF,subset_size,subset_size,matdata+cum2,&deluxe_ctx->seq_mat[i]));
516839e9adbSstefano_zampini 
51772b8c272SStefano Zampini     /* \sum_k S^k_E_j */
5189566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&deluxe_ctx->seq_mat_inv_sum[i]));
5199566063dSJacob Faibussowitsch     PetscCall(MatCreateSeqDense(PETSC_COMM_SELF,subset_size,subset_size,matdata2+cum2,&deluxe_ctx->seq_mat_inv_sum[i]));
5209566063dSJacob Faibussowitsch     PetscCall(MatSetOption(deluxe_ctx->seq_mat_inv_sum[i],MAT_SPD,sub_schurs->is_posdef));
5219566063dSJacob Faibussowitsch     PetscCall(MatSetOption(deluxe_ctx->seq_mat_inv_sum[i],MAT_HERMITIAN,sub_schurs->is_hermitian));
5223829e7e6SStefano Zampini     if (sub_schurs->is_hermitian) {
5239566063dSJacob Faibussowitsch       PetscCall(MatCholeskyFactor(deluxe_ctx->seq_mat_inv_sum[i],NULL,NULL));
524d62866d3SStefano Zampini     } else {
5259566063dSJacob Faibussowitsch       PetscCall(MatLUFactor(deluxe_ctx->seq_mat_inv_sum[i],NULL,NULL,NULL));
526d62866d3SStefano Zampini     }
527839e9adbSstefano_zampini     if (pcbddc->deluxe_singlemat) {
528839e9adbSstefano_zampini       Mat X,Y;
5293829e7e6SStefano Zampini       if (!sub_schurs->is_hermitian) {
5309566063dSJacob Faibussowitsch         PetscCall(MatTranspose(deluxe_ctx->seq_mat[i],MAT_INITIAL_MATRIX,&X));
531839e9adbSstefano_zampini       } else {
5329566063dSJacob Faibussowitsch         PetscCall(PetscObjectReference((PetscObject)deluxe_ctx->seq_mat[i]));
533839e9adbSstefano_zampini         X    = deluxe_ctx->seq_mat[i];
534839e9adbSstefano_zampini       }
5359566063dSJacob Faibussowitsch       PetscCall(MatDuplicate(X,MAT_DO_NOT_COPY_VALUES,&Y));
5363829e7e6SStefano Zampini       if (!sub_schurs->is_hermitian) {
5379566063dSJacob Faibussowitsch         PetscCall(PCBDDCMatTransposeMatSolve_SeqDense(deluxe_ctx->seq_mat_inv_sum[i],X,Y));
538839e9adbSstefano_zampini       } else {
5399566063dSJacob Faibussowitsch         PetscCall(MatMatSolve(deluxe_ctx->seq_mat_inv_sum[i],X,Y));
540839e9adbSstefano_zampini       }
541729c86d3Sstefano_zampini 
5429566063dSJacob Faibussowitsch       PetscCall(MatDestroy(&deluxe_ctx->seq_mat_inv_sum[i]));
5439566063dSJacob Faibussowitsch       PetscCall(MatDestroy(&deluxe_ctx->seq_mat[i]));
5449566063dSJacob Faibussowitsch       PetscCall(MatDestroy(&X));
5452f41f7d1SStefano Zampini       if (deluxe_ctx->change) {
5462f41f7d1SStefano Zampini         Mat C,CY;
54728b400f6SJacob Faibussowitsch         PetscCheck(deluxe_ctx->change_with_qr,PETSC_COMM_SELF,PETSC_ERR_SUP,"Only QR based change of basis");
5489566063dSJacob Faibussowitsch         PetscCall(KSPGetOperators(deluxe_ctx->change[i],&C,NULL));
5499566063dSJacob Faibussowitsch         PetscCall(MatMatMult(C,Y,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&CY));
5509566063dSJacob Faibussowitsch         PetscCall(MatMatTransposeMult(CY,C,MAT_REUSE_MATRIX,PETSC_DEFAULT,&Y));
5519566063dSJacob Faibussowitsch         PetscCall(MatDestroy(&CY));
5529566063dSJacob Faibussowitsch         PetscCall(MatProductClear(Y)); /* clear internal matproduct structure of Y since CY is destroyed */
5532f41f7d1SStefano Zampini       }
5549566063dSJacob Faibussowitsch       PetscCall(MatTranspose(Y,MAT_INPLACE_MATRIX,&Y));
555839e9adbSstefano_zampini       deluxe_ctx->seq_mat[i] = Y;
55672b8c272SStefano Zampini     }
55772b8c272SStefano Zampini     cum += subset_size;
55872b8c272SStefano Zampini     cum2 += subset_size*subset_size;
55972b8c272SStefano Zampini   }
5609566063dSJacob Faibussowitsch   PetscCall(ISRestoreIndices(sub_schurs->is_Ej_all,&idxs));
5619566063dSJacob Faibussowitsch   PetscCall(MatSeqAIJRestoreArray(sub_schurs->S_Ej_all,&matdata));
5629566063dSJacob Faibussowitsch   PetscCall(MatSeqAIJRestoreArray(sub_schurs->sum_S_Ej_all,&matdata2));
5632f41f7d1SStefano Zampini   if (pcbddc->deluxe_singlemat) {
5642f41f7d1SStefano Zampini     deluxe_ctx->change         = NULL;
5652f41f7d1SStefano Zampini     deluxe_ctx->change_with_qr = PETSC_FALSE;
5662f41f7d1SStefano Zampini   }
5675db18549SStefano Zampini 
56872b8c272SStefano Zampini   if (deluxe_ctx->change && !deluxe_ctx->change_with_qr) {
56972b8c272SStefano Zampini     for (i=0;i<deluxe_ctx->seq_n;i++) {
57072b8c272SStefano Zampini       if (newsetup) {
57172b8c272SStefano Zampini         PC pc;
57272b8c272SStefano Zampini 
5739566063dSJacob Faibussowitsch         PetscCall(KSPGetPC(deluxe_ctx->change[i],&pc));
5749566063dSJacob Faibussowitsch         PetscCall(PCSetType(pc,PCLU));
5759566063dSJacob Faibussowitsch         PetscCall(KSPSetFromOptions(deluxe_ctx->change[i]));
57672b8c272SStefano Zampini       }
5779566063dSJacob Faibussowitsch       PetscCall(KSPSetUp(deluxe_ctx->change[i]));
57872b8c272SStefano Zampini     }
57916e386b8SStefano Zampini   }
5805db18549SStefano Zampini   PetscFunctionReturn(0);
5815db18549SStefano Zampini }
582