xref: /petsc/src/ksp/pc/impls/bddc/bddcscalingbasic.c (revision 03dfb2d78eb3e713a8329a5c7dca01d2e8c2c7cf)
1ab5c6b0cSJed Brown #include <../src/ksp/pc/impls/bddc/bddc.h>
2ab5c6b0cSJed Brown #include <../src/ksp/pc/impls/bddc/bddcprivate.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;
16839e9adbSstefano_zampini   PetscErrorCode ierr;
17839e9adbSstefano_zampini   PetscScalar    *b,*x;
18839e9adbSstefano_zampini   PetscInt       n;
19839e9adbSstefano_zampini   PetscBLASInt   nrhs,info,m;
20839e9adbSstefano_zampini   PetscBool      flg;
21839e9adbSstefano_zampini 
22839e9adbSstefano_zampini   PetscFunctionBegin;
23839e9adbSstefano_zampini   ierr = PetscBLASIntCast(A->rmap->n,&m);CHKERRQ(ierr);
24839e9adbSstefano_zampini   ierr = PetscObjectTypeCompareAny((PetscObject)B,&flg,MATSEQDENSE,MATMPIDENSE,NULL);CHKERRQ(ierr);
25839e9adbSstefano_zampini   if (!flg) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONG,"Matrix B must be MATDENSE matrix");
26839e9adbSstefano_zampini   ierr = PetscObjectTypeCompareAny((PetscObject)X,&flg,MATSEQDENSE,MATMPIDENSE,NULL);CHKERRQ(ierr);
27839e9adbSstefano_zampini   if (!flg) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONG,"Matrix X must be MATDENSE matrix");
28839e9adbSstefano_zampini 
29839e9adbSstefano_zampini   ierr = MatGetSize(B,NULL,&n);CHKERRQ(ierr);
30839e9adbSstefano_zampini   ierr = PetscBLASIntCast(n,&nrhs);CHKERRQ(ierr);
31839e9adbSstefano_zampini   ierr = MatDenseGetArray(B,&b);CHKERRQ(ierr);
32839e9adbSstefano_zampini   ierr = MatDenseGetArray(X,&x);CHKERRQ(ierr);
33839e9adbSstefano_zampini 
34839e9adbSstefano_zampini   ierr = PetscMemcpy(x,b,m*nrhs*sizeof(PetscScalar));CHKERRQ(ierr);
35839e9adbSstefano_zampini 
36839e9adbSstefano_zampini   if (A->factortype == MAT_FACTOR_LU) {
37839e9adbSstefano_zampini #if defined(PETSC_MISSING_LAPACK_GETRS)
38839e9adbSstefano_zampini     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"GETRS - Lapack routine is unavailable.");
39839e9adbSstefano_zampini #else
40839e9adbSstefano_zampini     PetscStackCallBLAS("LAPACKgetrs",LAPACKgetrs_("T",&m,&nrhs,mat->v,&mat->lda,mat->pivots,x,&m,&info));
41839e9adbSstefano_zampini     if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"GETRS - Bad solve");
42839e9adbSstefano_zampini #endif
43839e9adbSstefano_zampini   } else if (A->factortype == MAT_FACTOR_CHOLESKY) {
44839e9adbSstefano_zampini #if defined(PETSC_MISSING_LAPACK_POTRS)
45839e9adbSstefano_zampini     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"POTRS - Lapack routine is unavailable.");
46839e9adbSstefano_zampini #else
47839e9adbSstefano_zampini     PetscStackCallBLAS("LAPACKpotrs",LAPACKpotrs_("L",&m,&nrhs,mat->v,&mat->lda,x,&m,&info));
48839e9adbSstefano_zampini     if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"POTRS Bad solve");
49839e9adbSstefano_zampini #endif
50839e9adbSstefano_zampini   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Matrix must be factored to solve");
51839e9adbSstefano_zampini 
52839e9adbSstefano_zampini   ierr = MatDenseRestoreArray(B,&b);CHKERRQ(ierr);
53839e9adbSstefano_zampini   ierr = MatDenseRestoreArray(X,&x);CHKERRQ(ierr);
54839e9adbSstefano_zampini   ierr = PetscLogFlops(nrhs*(2.0*m*m - m));CHKERRQ(ierr);
55839e9adbSstefano_zampini   PetscFunctionReturn(0);
56839e9adbSstefano_zampini }
57839e9adbSstefano_zampini 
58839e9adbSstefano_zampini 
59674ae819SStefano Zampini static PetscErrorCode PCBDDCScalingExtension_Basic(PC pc, Vec local_interface_vector, Vec global_vector)
60674ae819SStefano Zampini {
61674ae819SStefano Zampini   PC_IS*         pcis = (PC_IS*)pc->data;
62674ae819SStefano Zampini   PC_BDDC*       pcbddc = (PC_BDDC*)pc->data;
63674ae819SStefano Zampini   PetscErrorCode ierr;
64674ae819SStefano Zampini 
65674ae819SStefano Zampini   PetscFunctionBegin;
66674ae819SStefano Zampini   /* Apply partition of unity */
67674ae819SStefano Zampini   ierr = VecPointwiseMult(pcbddc->work_scaling,pcis->D,local_interface_vector);CHKERRQ(ierr);
68674ae819SStefano Zampini   ierr = VecSet(global_vector,0.0);CHKERRQ(ierr);
69674ae819SStefano Zampini   ierr = VecScatterBegin(pcis->global_to_B,pcbddc->work_scaling,global_vector,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
70674ae819SStefano Zampini   ierr = VecScatterEnd(pcis->global_to_B,pcbddc->work_scaling,global_vector,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
71674ae819SStefano Zampini   PetscFunctionReturn(0);
72674ae819SStefano Zampini }
73674ae819SStefano Zampini 
74674ae819SStefano Zampini static PetscErrorCode PCBDDCScalingExtension_Deluxe(PC pc, Vec x, Vec y)
75674ae819SStefano Zampini {
76674ae819SStefano Zampini   PC_IS*              pcis=(PC_IS*)pc->data;
77674ae819SStefano Zampini   PC_BDDC*            pcbddc=(PC_BDDC*)pc->data;
78674ae819SStefano Zampini   PCBDDCDeluxeScaling deluxe_ctx = pcbddc->deluxe_ctx;
79674ae819SStefano Zampini   PetscErrorCode      ierr;
80674ae819SStefano Zampini 
81674ae819SStefano Zampini   PetscFunctionBegin;
825a95e1ceSStefano Zampini   ierr = VecSet(pcbddc->work_scaling,0.0);CHKERRQ(ierr);
835a95e1ceSStefano Zampini   ierr = VecSet(y,0.0);CHKERRQ(ierr);
845a95e1ceSStefano Zampini   if (deluxe_ctx->n_simple) { /* scale deluxe vertices using diagonal scaling */
855a95e1ceSStefano Zampini     PetscInt          i;
862b095fd8SStefano Zampini     const PetscScalar *array_x,*array_D;
872b095fd8SStefano Zampini     PetscScalar       *array;
882b095fd8SStefano Zampini     ierr = VecGetArrayRead(x,&array_x);CHKERRQ(ierr);
892b095fd8SStefano Zampini     ierr = VecGetArrayRead(pcis->D,&array_D);CHKERRQ(ierr);
90674ae819SStefano Zampini     ierr = VecGetArray(pcbddc->work_scaling,&array);CHKERRQ(ierr);
91674ae819SStefano Zampini     for (i=0;i<deluxe_ctx->n_simple;i++) {
92674ae819SStefano Zampini       array[deluxe_ctx->idx_simple_B[i]] = array_x[deluxe_ctx->idx_simple_B[i]]*array_D[deluxe_ctx->idx_simple_B[i]];
93674ae819SStefano Zampini     }
94674ae819SStefano Zampini     ierr = VecRestoreArray(pcbddc->work_scaling,&array);CHKERRQ(ierr);
952b095fd8SStefano Zampini     ierr = VecRestoreArrayRead(pcis->D,&array_D);CHKERRQ(ierr);
962b095fd8SStefano Zampini     ierr = VecRestoreArrayRead(x,&array_x);CHKERRQ(ierr);
9734a97f8cSStefano Zampini   }
98ac632422SStefano Zampini   /* sequential part : all problems and Schur applications collapsed into a single matrix vector multiplication or a matvec and a solve */
995a95e1ceSStefano Zampini   if (deluxe_ctx->seq_mat) {
10072b8c272SStefano Zampini     PetscInt i;
10172b8c272SStefano Zampini     for (i=0;i<deluxe_ctx->seq_n;i++) {
10216e386b8SStefano Zampini       if (deluxe_ctx->change) {
10372b8c272SStefano Zampini         ierr = VecScatterBegin(deluxe_ctx->seq_scctx[i],x,deluxe_ctx->seq_work2[i],INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
10472b8c272SStefano Zampini         ierr = VecScatterEnd(deluxe_ctx->seq_scctx[i],x,deluxe_ctx->seq_work2[i],INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
10572b8c272SStefano Zampini         if (deluxe_ctx->change_with_qr) {
10672b8c272SStefano Zampini           Mat change;
10772b8c272SStefano Zampini 
10872b8c272SStefano Zampini           ierr = KSPGetOperators(deluxe_ctx->change[i],&change,NULL);CHKERRQ(ierr);
10972b8c272SStefano Zampini           ierr = MatMultTranspose(change,deluxe_ctx->seq_work2[i],deluxe_ctx->seq_work1[i]);CHKERRQ(ierr);
11072b8c272SStefano Zampini         } else {
11172b8c272SStefano Zampini           ierr = KSPSolve(deluxe_ctx->change[i],deluxe_ctx->seq_work2[i],deluxe_ctx->seq_work1[i]);CHKERRQ(ierr);
11216e386b8SStefano Zampini         }
11372b8c272SStefano Zampini       } else {
11472b8c272SStefano Zampini         ierr = VecScatterBegin(deluxe_ctx->seq_scctx[i],x,deluxe_ctx->seq_work1[i],INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
11572b8c272SStefano Zampini         ierr = VecScatterEnd(deluxe_ctx->seq_scctx[i],x,deluxe_ctx->seq_work1[i],INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
11672b8c272SStefano Zampini       }
11772b8c272SStefano Zampini       ierr = MatMultTranspose(deluxe_ctx->seq_mat[i],deluxe_ctx->seq_work1[i],deluxe_ctx->seq_work2[i]);CHKERRQ(ierr);
11872b8c272SStefano Zampini       if (deluxe_ctx->seq_mat_inv_sum[i]) {
11972b8c272SStefano Zampini         PetscScalar *x;
12072b8c272SStefano Zampini 
12172b8c272SStefano Zampini         ierr = VecGetArray(deluxe_ctx->seq_work2[i],&x);CHKERRQ(ierr);
12272b8c272SStefano Zampini         ierr = VecPlaceArray(deluxe_ctx->seq_work1[i],x);CHKERRQ(ierr);
12372b8c272SStefano Zampini         ierr = VecRestoreArray(deluxe_ctx->seq_work2[i],&x);CHKERRQ(ierr);
12472b8c272SStefano Zampini         ierr = MatSolveTranspose(deluxe_ctx->seq_mat_inv_sum[i],deluxe_ctx->seq_work1[i],deluxe_ctx->seq_work2[i]);CHKERRQ(ierr);
12572b8c272SStefano Zampini         ierr = VecResetArray(deluxe_ctx->seq_work1[i]);CHKERRQ(ierr);
126ac632422SStefano Zampini       }
12716e386b8SStefano Zampini       if (deluxe_ctx->change) {
12816e386b8SStefano Zampini         Mat change;
12916e386b8SStefano Zampini 
13072b8c272SStefano Zampini         ierr = KSPGetOperators(deluxe_ctx->change[i],&change,NULL);CHKERRQ(ierr);
13172b8c272SStefano Zampini         ierr = MatMult(change,deluxe_ctx->seq_work2[i],deluxe_ctx->seq_work1[i]);CHKERRQ(ierr);
13272b8c272SStefano Zampini         ierr = VecScatterBegin(deluxe_ctx->seq_scctx[i],deluxe_ctx->seq_work1[i],pcbddc->work_scaling,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
13372b8c272SStefano Zampini         ierr = VecScatterEnd(deluxe_ctx->seq_scctx[i],deluxe_ctx->seq_work1[i],pcbddc->work_scaling,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
13416e386b8SStefano Zampini       } else {
13572b8c272SStefano Zampini         ierr = VecScatterBegin(deluxe_ctx->seq_scctx[i],deluxe_ctx->seq_work2[i],pcbddc->work_scaling,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
13672b8c272SStefano Zampini         ierr = VecScatterEnd(deluxe_ctx->seq_scctx[i],deluxe_ctx->seq_work2[i],pcbddc->work_scaling,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
13772b8c272SStefano Zampini       }
138674ae819SStefano Zampini     }
13916e386b8SStefano Zampini   }
140674ae819SStefano Zampini   /* put local boundary part in global vector */
141674ae819SStefano Zampini   ierr = VecScatterBegin(pcis->global_to_B,pcbddc->work_scaling,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
142674ae819SStefano Zampini   ierr = VecScatterEnd(pcis->global_to_B,pcbddc->work_scaling,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
143674ae819SStefano Zampini   PetscFunctionReturn(0);
144674ae819SStefano Zampini }
145674ae819SStefano Zampini 
146674ae819SStefano Zampini PetscErrorCode PCBDDCScalingExtension(PC pc, Vec local_interface_vector, Vec global_vector)
147674ae819SStefano Zampini {
148674ae819SStefano Zampini   PC_BDDC        *pcbddc=(PC_BDDC*)pc->data;
149674ae819SStefano Zampini   PetscErrorCode ierr;
150674ae819SStefano Zampini 
151674ae819SStefano Zampini   PetscFunctionBegin;
152674ae819SStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
153674ae819SStefano Zampini   PetscValidHeaderSpecific(local_interface_vector,VEC_CLASSID,2);
154674ae819SStefano Zampini   PetscValidHeaderSpecific(global_vector,VEC_CLASSID,3);
1556c4ed002SBarry Smith   if (local_interface_vector == pcbddc->work_scaling) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Local vector cannot be pcbddc->work_scaling!\n");
156163d334eSBarry Smith   ierr = PetscUseMethod(pc,"PCBDDCScalingExtension_C",(PC,Vec,Vec),(pc,local_interface_vector,global_vector));CHKERRQ(ierr);
157674ae819SStefano Zampini   PetscFunctionReturn(0);
158674ae819SStefano Zampini }
159674ae819SStefano Zampini 
160674ae819SStefano Zampini static PetscErrorCode PCBDDCScalingRestriction_Basic(PC pc, Vec global_vector, Vec local_interface_vector)
161674ae819SStefano Zampini {
162674ae819SStefano Zampini   PetscErrorCode ierr;
163674ae819SStefano Zampini   PC_IS          *pcis = (PC_IS*)pc->data;
164674ae819SStefano Zampini 
165674ae819SStefano Zampini   PetscFunctionBegin;
166674ae819SStefano Zampini   ierr = VecScatterBegin(pcis->global_to_B,global_vector,local_interface_vector,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
167674ae819SStefano Zampini   ierr = VecScatterEnd(pcis->global_to_B,global_vector,local_interface_vector,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
168674ae819SStefano Zampini   /* Apply partition of unity */
169674ae819SStefano Zampini   ierr = VecPointwiseMult(local_interface_vector,pcis->D,local_interface_vector);CHKERRQ(ierr);
170674ae819SStefano Zampini   PetscFunctionReturn(0);
171674ae819SStefano Zampini }
172674ae819SStefano Zampini 
173674ae819SStefano Zampini static PetscErrorCode PCBDDCScalingRestriction_Deluxe(PC pc, Vec x, Vec y)
174674ae819SStefano Zampini {
175674ae819SStefano Zampini   PC_IS*              pcis=(PC_IS*)pc->data;
176674ae819SStefano Zampini   PC_BDDC*            pcbddc=(PC_BDDC*)pc->data;
177674ae819SStefano Zampini   PCBDDCDeluxeScaling deluxe_ctx = pcbddc->deluxe_ctx;
178674ae819SStefano Zampini   PetscErrorCode      ierr;
179674ae819SStefano Zampini 
180674ae819SStefano Zampini   PetscFunctionBegin;
181674ae819SStefano Zampini   /* get local boundary part of global vector */
182674ae819SStefano Zampini   ierr = VecScatterBegin(pcis->global_to_B,x,y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
183674ae819SStefano Zampini   ierr = VecScatterEnd(pcis->global_to_B,x,y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1845a95e1ceSStefano Zampini   if (deluxe_ctx->n_simple) { /* scale deluxe vertices using diagonal scaling */
1855a95e1ceSStefano Zampini     PetscInt          i;
1862b095fd8SStefano Zampini     PetscScalar       *array_y;
1872b095fd8SStefano Zampini     const PetscScalar *array_D;
188674ae819SStefano Zampini     ierr = VecGetArray(y,&array_y);CHKERRQ(ierr);
1892b095fd8SStefano Zampini     ierr = VecGetArrayRead(pcis->D,&array_D);CHKERRQ(ierr);
190674ae819SStefano Zampini     for (i=0;i<deluxe_ctx->n_simple;i++) {
191674ae819SStefano Zampini       array_y[deluxe_ctx->idx_simple_B[i]] *= array_D[deluxe_ctx->idx_simple_B[i]];
192674ae819SStefano Zampini     }
1932b095fd8SStefano Zampini     ierr = VecRestoreArrayRead(pcis->D,&array_D);CHKERRQ(ierr);
194674ae819SStefano Zampini     ierr = VecRestoreArray(y,&array_y);CHKERRQ(ierr);
19534a97f8cSStefano Zampini   }
196729c86d3Sstefano_zampini   /* sequential part : all problems and Schur applications collapsed into a single matrix vector multiplication or a matvec and a solve */
1975a95e1ceSStefano Zampini   if (deluxe_ctx->seq_mat) {
19872b8c272SStefano Zampini     PetscInt i;
19972b8c272SStefano Zampini     for (i=0;i<deluxe_ctx->seq_n;i++) {
20016e386b8SStefano Zampini       if (deluxe_ctx->change) {
20116e386b8SStefano Zampini         Mat change;
20216e386b8SStefano Zampini 
20372b8c272SStefano Zampini         ierr = VecScatterBegin(deluxe_ctx->seq_scctx[i],y,deluxe_ctx->seq_work2[i],INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
20472b8c272SStefano Zampini         ierr = VecScatterEnd(deluxe_ctx->seq_scctx[i],y,deluxe_ctx->seq_work2[i],INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
20572b8c272SStefano Zampini         ierr = KSPGetOperators(deluxe_ctx->change[i],&change,NULL);CHKERRQ(ierr);
20672b8c272SStefano Zampini         ierr = MatMultTranspose(change,deluxe_ctx->seq_work2[i],deluxe_ctx->seq_work1[i]);CHKERRQ(ierr);
20716e386b8SStefano Zampini       } else {
20872b8c272SStefano Zampini         ierr = VecScatterBegin(deluxe_ctx->seq_scctx[i],y,deluxe_ctx->seq_work1[i],INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
20972b8c272SStefano Zampini         ierr = VecScatterEnd(deluxe_ctx->seq_scctx[i],y,deluxe_ctx->seq_work1[i],INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
21016e386b8SStefano Zampini       }
21172b8c272SStefano Zampini       if (deluxe_ctx->seq_mat_inv_sum[i]) {
21272b8c272SStefano Zampini         PetscScalar *x;
21372b8c272SStefano Zampini 
21472b8c272SStefano Zampini         ierr = VecGetArray(deluxe_ctx->seq_work1[i],&x);CHKERRQ(ierr);
21572b8c272SStefano Zampini         ierr = VecPlaceArray(deluxe_ctx->seq_work2[i],x);CHKERRQ(ierr);
21672b8c272SStefano Zampini         ierr = VecRestoreArray(deluxe_ctx->seq_work1[i],&x);CHKERRQ(ierr);
21772b8c272SStefano Zampini         ierr = MatSolve(deluxe_ctx->seq_mat_inv_sum[i],deluxe_ctx->seq_work1[i],deluxe_ctx->seq_work2[i]);CHKERRQ(ierr);
21872b8c272SStefano Zampini         ierr = VecResetArray(deluxe_ctx->seq_work2[i]);CHKERRQ(ierr);
219ac632422SStefano Zampini       }
22072b8c272SStefano Zampini       ierr = MatMult(deluxe_ctx->seq_mat[i],deluxe_ctx->seq_work1[i],deluxe_ctx->seq_work2[i]);CHKERRQ(ierr);
22116e386b8SStefano Zampini       if (deluxe_ctx->change) {
22272b8c272SStefano Zampini         if (deluxe_ctx->change_with_qr) {
22372b8c272SStefano Zampini           Mat change;
22472b8c272SStefano Zampini 
22572b8c272SStefano Zampini           ierr = KSPGetOperators(deluxe_ctx->change[i],&change,NULL);CHKERRQ(ierr);
22672b8c272SStefano Zampini           ierr = MatMult(change,deluxe_ctx->seq_work2[i],deluxe_ctx->seq_work1[i]);CHKERRQ(ierr);
22772b8c272SStefano Zampini         } else {
22872b8c272SStefano Zampini           ierr = KSPSolveTranspose(deluxe_ctx->change[i],deluxe_ctx->seq_work2[i],deluxe_ctx->seq_work1[i]);CHKERRQ(ierr);
22916e386b8SStefano Zampini         }
23072b8c272SStefano Zampini         ierr = VecScatterBegin(deluxe_ctx->seq_scctx[i],deluxe_ctx->seq_work1[i],y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
23172b8c272SStefano Zampini         ierr = VecScatterEnd(deluxe_ctx->seq_scctx[i],deluxe_ctx->seq_work1[i],y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
23272b8c272SStefano Zampini       } else {
23372b8c272SStefano Zampini         ierr = VecScatterBegin(deluxe_ctx->seq_scctx[i],deluxe_ctx->seq_work2[i],y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
23472b8c272SStefano Zampini         ierr = VecScatterEnd(deluxe_ctx->seq_scctx[i],deluxe_ctx->seq_work2[i],y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
23572b8c272SStefano Zampini       }
23672b8c272SStefano Zampini     }
237674ae819SStefano Zampini   }
238674ae819SStefano Zampini   PetscFunctionReturn(0);
239674ae819SStefano Zampini }
240674ae819SStefano Zampini 
241674ae819SStefano Zampini PetscErrorCode PCBDDCScalingRestriction(PC pc, Vec global_vector, Vec local_interface_vector)
242674ae819SStefano Zampini {
243674ae819SStefano Zampini   PC_BDDC        *pcbddc=(PC_BDDC*)pc->data;
244674ae819SStefano Zampini   PetscErrorCode ierr;
245674ae819SStefano Zampini 
246674ae819SStefano Zampini   PetscFunctionBegin;
247674ae819SStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
248674ae819SStefano Zampini   PetscValidHeaderSpecific(global_vector,VEC_CLASSID,2);
249674ae819SStefano Zampini   PetscValidHeaderSpecific(local_interface_vector,VEC_CLASSID,3);
2506c4ed002SBarry Smith   if (local_interface_vector == pcbddc->work_scaling) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Local vector cannot be pcbddc->work_scaling!\n");
251163d334eSBarry Smith   ierr = PetscUseMethod(pc,"PCBDDCScalingRestriction_C",(PC,Vec,Vec),(pc,global_vector,local_interface_vector));CHKERRQ(ierr);
252674ae819SStefano Zampini   PetscFunctionReturn(0);
253674ae819SStefano Zampini }
254674ae819SStefano Zampini 
255674ae819SStefano Zampini PetscErrorCode PCBDDCScalingSetUp(PC pc)
256674ae819SStefano Zampini {
257674ae819SStefano Zampini   PC_IS*         pcis=(PC_IS*)pc->data;
258674ae819SStefano Zampini   PC_BDDC*       pcbddc=(PC_BDDC*)pc->data;
259674ae819SStefano Zampini   PetscErrorCode ierr;
260674ae819SStefano Zampini 
261674ae819SStefano Zampini   PetscFunctionBegin;
262674ae819SStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
263674ae819SStefano Zampini   /* create work vector for the operator */
26434a97f8cSStefano Zampini   ierr = VecDestroy(&pcbddc->work_scaling);CHKERRQ(ierr);
265674ae819SStefano Zampini   ierr = VecDuplicate(pcis->vec1_B,&pcbddc->work_scaling);CHKERRQ(ierr);
26634a97f8cSStefano Zampini   /* always rebuild pcis->D */
26728d874f6SStefano Zampini   if (pcis->use_stiffness_scaling) {
268674ae819SStefano Zampini     ierr = MatGetDiagonal(pcbddc->local_mat,pcis->vec1_N);CHKERRQ(ierr);
269674ae819SStefano Zampini     ierr = VecScatterBegin(pcis->N_to_B,pcis->vec1_N,pcis->D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
270674ae819SStefano Zampini     ierr = VecScatterEnd(pcis->N_to_B,pcis->vec1_N,pcis->D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
271674ae819SStefano Zampini   }
272674ae819SStefano Zampini   ierr = VecCopy(pcis->D,pcis->vec1_B);CHKERRQ(ierr);
273674ae819SStefano Zampini   ierr = VecSet(pcis->vec1_global,0.0);CHKERRQ(ierr);
274674ae819SStefano Zampini   ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_B,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
275674ae819SStefano Zampini   ierr = VecScatterEnd(pcis->global_to_B,pcis->vec1_B,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
276674ae819SStefano Zampini   ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_global,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
277674ae819SStefano Zampini   ierr = VecScatterEnd(pcis->global_to_B,pcis->vec1_global,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
278674ae819SStefano Zampini   ierr = VecPointwiseDivide(pcis->D,pcis->D,pcis->vec1_B);CHKERRQ(ierr);
279674ae819SStefano Zampini   /* now setup */
280681e7c04SStefano Zampini   if (pcbddc->use_deluxe_scaling) {
28134a97f8cSStefano Zampini     if (!pcbddc->deluxe_ctx) {
28234a97f8cSStefano Zampini       ierr = PCBDDCScalingCreate_Deluxe(pc);CHKERRQ(ierr);
28334a97f8cSStefano Zampini     }
28434a97f8cSStefano Zampini     ierr = PCBDDCScalingSetUp_Deluxe(pc);CHKERRQ(ierr);
285674ae819SStefano Zampini     ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCScalingRestriction_C",PCBDDCScalingRestriction_Deluxe);CHKERRQ(ierr);
286674ae819SStefano Zampini     ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCScalingExtension_C",PCBDDCScalingExtension_Deluxe);CHKERRQ(ierr);
287674ae819SStefano Zampini   } else {
288674ae819SStefano Zampini     ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCScalingRestriction_C",PCBDDCScalingRestriction_Basic);CHKERRQ(ierr);
289674ae819SStefano Zampini     ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCScalingExtension_C",PCBDDCScalingExtension_Basic);CHKERRQ(ierr);
290674ae819SStefano Zampini   }
29134a97f8cSStefano Zampini 
292674ae819SStefano Zampini   /* test */
293674ae819SStefano Zampini   if (pcbddc->dbg_flag) {
29472b8c272SStefano Zampini     Mat         B0_B = NULL;
29572b8c272SStefano Zampini     Vec         B0_Bv = NULL, B0_Bv2 = NULL;
29634a97f8cSStefano Zampini     Vec         vec2_global;
297674ae819SStefano Zampini     PetscViewer viewer = pcbddc->dbg_viewer;
29834a97f8cSStefano Zampini     PetscReal   error;
299674ae819SStefano Zampini 
300674ae819SStefano Zampini     /* extension -> from local to parallel */
30134a97f8cSStefano Zampini     ierr = VecSet(pcis->vec1_global,0.0);CHKERRQ(ierr);
30234a97f8cSStefano Zampini     ierr = VecSetRandom(pcis->vec1_B,NULL);CHKERRQ(ierr);
30334a97f8cSStefano Zampini     ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_B,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
30434a97f8cSStefano Zampini     ierr = VecScatterEnd(pcis->global_to_B,pcis->vec1_B,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
30534a97f8cSStefano Zampini     ierr = VecDuplicate(pcis->vec1_global,&vec2_global);CHKERRQ(ierr);
30634a97f8cSStefano Zampini     ierr = VecCopy(pcis->vec1_global,vec2_global);CHKERRQ(ierr);
307674ae819SStefano Zampini     ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_global,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
308674ae819SStefano Zampini     ierr = VecScatterEnd(pcis->global_to_B,pcis->vec1_global,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
30972b8c272SStefano Zampini     if (pcbddc->benign_n) {
31072b8c272SStefano Zampini       IS is_dummy;
31172b8c272SStefano Zampini 
31272b8c272SStefano Zampini       ierr = ISCreateStride(PETSC_COMM_SELF,pcbddc->benign_n,0,1,&is_dummy);CHKERRQ(ierr);
3137dae84e0SHong Zhang       ierr = MatCreateSubMatrix(pcbddc->benign_B0,is_dummy,pcis->is_B_local,MAT_INITIAL_MATRIX,&B0_B);CHKERRQ(ierr);
31472b8c272SStefano Zampini       ierr = ISDestroy(&is_dummy);CHKERRQ(ierr);
31572b8c272SStefano Zampini       ierr = MatCreateVecs(B0_B,NULL,&B0_Bv);CHKERRQ(ierr);
31672b8c272SStefano Zampini       ierr = VecDuplicate(B0_Bv,&B0_Bv2);CHKERRQ(ierr);
31772b8c272SStefano Zampini       ierr = MatMult(B0_B,pcis->vec1_B,B0_Bv);CHKERRQ(ierr);
31872b8c272SStefano Zampini     }
319674ae819SStefano Zampini     ierr = PCBDDCScalingExtension(pc,pcis->vec1_B,pcis->vec1_global);CHKERRQ(ierr);
32072b8c272SStefano Zampini     if (pcbddc->benign_saddle_point) {
32172b8c272SStefano Zampini       PetscReal errorl = 0.;
32272b8c272SStefano Zampini       ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_global,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
32372b8c272SStefano Zampini       ierr = VecScatterEnd(pcis->global_to_B,pcis->vec1_global,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
32472b8c272SStefano Zampini       if (pcbddc->benign_n) {
32572b8c272SStefano Zampini         ierr = MatMult(B0_B,pcis->vec1_B,B0_Bv2);CHKERRQ(ierr);
32672b8c272SStefano Zampini         ierr = VecAXPY(B0_Bv,-1.0,B0_Bv2);CHKERRQ(ierr);
32772b8c272SStefano Zampini         ierr = VecNorm(B0_Bv,NORM_INFINITY,&errorl);CHKERRQ(ierr);
32872b8c272SStefano Zampini       }
32972b8c272SStefano Zampini       ierr = MPI_Allreduce(&errorl,&error,1,MPIU_REAL,MPI_SUM,PetscObjectComm((PetscObject)pc));CHKERRQ(ierr);
33072b8c272SStefano Zampini       ierr = PetscViewerASCIIPrintf(viewer,"Error benign extension %1.14e\n",error);CHKERRQ(ierr);
33172b8c272SStefano Zampini     }
33234a97f8cSStefano Zampini     ierr = VecAXPY(pcis->vec1_global,-1.0,vec2_global);CHKERRQ(ierr);
33334a97f8cSStefano Zampini     ierr = VecNorm(pcis->vec1_global,NORM_INFINITY,&error);CHKERRQ(ierr);
334674ae819SStefano Zampini     ierr = PetscViewerASCIIPrintf(viewer,"Error scaling extension %1.14e\n",error);CHKERRQ(ierr);
33534a97f8cSStefano Zampini     ierr = VecDestroy(&vec2_global);CHKERRQ(ierr);
33634a97f8cSStefano Zampini 
337674ae819SStefano Zampini     /* restriction -> from parallel to local */
338674ae819SStefano Zampini     ierr = VecSet(pcis->vec1_global,0.0);CHKERRQ(ierr);
33934a97f8cSStefano Zampini     ierr = VecSetRandom(pcis->vec1_B,NULL);CHKERRQ(ierr);
340674ae819SStefano Zampini     ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_B,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
341674ae819SStefano Zampini     ierr = VecScatterEnd(pcis->global_to_B,pcis->vec1_B,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
34234a97f8cSStefano Zampini     ierr = PCBDDCScalingRestriction(pc,pcis->vec1_global,pcis->vec1_B);CHKERRQ(ierr);
34334a97f8cSStefano Zampini     ierr = VecScale(pcis->vec1_B,-1.0);CHKERRQ(ierr);
34434a97f8cSStefano Zampini     ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_B,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
34534a97f8cSStefano Zampini     ierr = VecScatterEnd(pcis->global_to_B,pcis->vec1_B,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
34634a97f8cSStefano Zampini     ierr = VecNorm(pcis->vec1_global,NORM_INFINITY,&error);CHKERRQ(ierr);
34734a97f8cSStefano Zampini     ierr = PetscViewerASCIIPrintf(viewer,"Error scaling restriction %1.14e\n",error);CHKERRQ(ierr);
34872b8c272SStefano Zampini     ierr = MatDestroy(&B0_B);CHKERRQ(ierr);
34972b8c272SStefano Zampini     ierr = VecDestroy(&B0_Bv);CHKERRQ(ierr);
35072b8c272SStefano Zampini     ierr = VecDestroy(&B0_Bv2);CHKERRQ(ierr);
351674ae819SStefano Zampini   }
352674ae819SStefano Zampini   PetscFunctionReturn(0);
353674ae819SStefano Zampini }
354674ae819SStefano Zampini 
355674ae819SStefano Zampini PetscErrorCode PCBDDCScalingDestroy(PC pc)
356674ae819SStefano Zampini {
357674ae819SStefano Zampini   PC_BDDC*       pcbddc=(PC_BDDC*)pc->data;
358674ae819SStefano Zampini   PetscErrorCode ierr;
359674ae819SStefano Zampini 
360674ae819SStefano Zampini   PetscFunctionBegin;
36134a97f8cSStefano Zampini   if (pcbddc->deluxe_ctx) {
36234a97f8cSStefano Zampini     ierr = PCBDDCScalingDestroy_Deluxe(pc);CHKERRQ(ierr);
363674ae819SStefano Zampini   }
364674ae819SStefano Zampini   ierr = VecDestroy(&pcbddc->work_scaling);CHKERRQ(ierr);
365674ae819SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCScalingRestriction_C",NULL);CHKERRQ(ierr);
366674ae819SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCScalingExtension_C",NULL);CHKERRQ(ierr);
367674ae819SStefano Zampini   PetscFunctionReturn(0);
368674ae819SStefano Zampini }
369674ae819SStefano Zampini 
37034a97f8cSStefano Zampini static PetscErrorCode PCBDDCScalingCreate_Deluxe(PC pc)
37134a97f8cSStefano Zampini {
37234a97f8cSStefano Zampini   PC_BDDC*            pcbddc=(PC_BDDC*)pc->data;
37334a97f8cSStefano Zampini   PCBDDCDeluxeScaling deluxe_ctx;
37434a97f8cSStefano Zampini   PetscErrorCode      ierr;
37534a97f8cSStefano Zampini 
37634a97f8cSStefano Zampini   PetscFunctionBegin;
37734a97f8cSStefano Zampini   ierr = PetscNew(&deluxe_ctx);CHKERRQ(ierr);
37834a97f8cSStefano Zampini   pcbddc->deluxe_ctx = deluxe_ctx;
37934a97f8cSStefano Zampini   PetscFunctionReturn(0);
38034a97f8cSStefano Zampini }
38134a97f8cSStefano Zampini 
38234a97f8cSStefano Zampini static PetscErrorCode PCBDDCScalingDestroy_Deluxe(PC pc)
38334a97f8cSStefano Zampini {
38434a97f8cSStefano Zampini   PC_BDDC*            pcbddc=(PC_BDDC*)pc->data;
38534a97f8cSStefano Zampini   PetscErrorCode      ierr;
38634a97f8cSStefano Zampini 
38734a97f8cSStefano Zampini   PetscFunctionBegin;
38834a97f8cSStefano Zampini   ierr = PCBDDCScalingReset_Deluxe_Solvers(pcbddc->deluxe_ctx);CHKERRQ(ierr);
38934a97f8cSStefano Zampini   ierr = PetscFree(pcbddc->deluxe_ctx);CHKERRQ(ierr);
39034a97f8cSStefano Zampini   PetscFunctionReturn(0);
39134a97f8cSStefano Zampini }
39234a97f8cSStefano Zampini 
39334a97f8cSStefano Zampini static PetscErrorCode PCBDDCScalingReset_Deluxe_Solvers(PCBDDCDeluxeScaling deluxe_ctx)
39434a97f8cSStefano Zampini {
39572b8c272SStefano Zampini   PetscInt       i;
39634a97f8cSStefano Zampini   PetscErrorCode ierr;
39734a97f8cSStefano Zampini 
39834a97f8cSStefano Zampini   PetscFunctionBegin;
39934a97f8cSStefano Zampini   ierr = PetscFree(deluxe_ctx->idx_simple_B);CHKERRQ(ierr);
40034a97f8cSStefano Zampini   deluxe_ctx->n_simple = 0;
40172b8c272SStefano Zampini   for (i=0;i<deluxe_ctx->seq_n;i++) {
40272b8c272SStefano Zampini     ierr = VecScatterDestroy(&deluxe_ctx->seq_scctx[i]);CHKERRQ(ierr);
40372b8c272SStefano Zampini     ierr = VecDestroy(&deluxe_ctx->seq_work1[i]);CHKERRQ(ierr);
40472b8c272SStefano Zampini     ierr = VecDestroy(&deluxe_ctx->seq_work2[i]);CHKERRQ(ierr);
40572b8c272SStefano Zampini     ierr = MatDestroy(&deluxe_ctx->seq_mat[i]);CHKERRQ(ierr);
40672b8c272SStefano Zampini     ierr = MatDestroy(&deluxe_ctx->seq_mat_inv_sum[i]);CHKERRQ(ierr);
40772b8c272SStefano Zampini   }
40872b8c272SStefano 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);
40972b8c272SStefano Zampini   ierr = PetscFree(deluxe_ctx->workspace);CHKERRQ(ierr);
41072b8c272SStefano Zampini   deluxe_ctx->seq_n = 0;
41134a97f8cSStefano Zampini   PetscFunctionReturn(0);
41234a97f8cSStefano Zampini }
41334a97f8cSStefano Zampini 
41434a97f8cSStefano Zampini static PetscErrorCode PCBDDCScalingSetUp_Deluxe(PC pc)
41534a97f8cSStefano Zampini {
41634a97f8cSStefano Zampini   PC_IS               *pcis=(PC_IS*)pc->data;
41734a97f8cSStefano Zampini   PC_BDDC             *pcbddc=(PC_BDDC*)pc->data;
41834a97f8cSStefano Zampini   PCBDDCDeluxeScaling deluxe_ctx=pcbddc->deluxe_ctx;
419b96c3477SStefano Zampini   PCBDDCSubSchurs     sub_schurs=pcbddc->sub_schurs;
42034a97f8cSStefano Zampini   PetscErrorCode      ierr;
42134a97f8cSStefano Zampini 
42234a97f8cSStefano Zampini   PetscFunctionBegin;
42372b8c272SStefano Zampini   /* reset data structures if the topology has changed */
42472b8c272SStefano Zampini   if (pcbddc->recompute_topography) {
42534a97f8cSStefano Zampini     ierr = PCBDDCScalingReset_Deluxe_Solvers(deluxe_ctx);CHKERRQ(ierr);
42672b8c272SStefano Zampini   }
42734a97f8cSStefano Zampini 
428b1b3d7a2SStefano Zampini   /* Compute data structures to solve sequential problems */
4295a95e1ceSStefano Zampini   ierr = PCBDDCScalingSetUp_Deluxe_Private(pc);CHKERRQ(ierr);
4305db18549SStefano Zampini 
431b1b3d7a2SStefano Zampini   /* diagonal scaling on interface dofs not contained in cc */
432d62866d3SStefano Zampini   if (sub_schurs->is_vertices || sub_schurs->is_dir) {
4331b968477SStefano Zampini     PetscInt n_com,n_dir;
4341b968477SStefano Zampini     n_com = 0;
435d62866d3SStefano Zampini     if (sub_schurs->is_vertices) {
436d62866d3SStefano Zampini       ierr = ISGetLocalSize(sub_schurs->is_vertices,&n_com);CHKERRQ(ierr);
4371b968477SStefano Zampini     }
4381b968477SStefano Zampini     n_dir = 0;
439d62866d3SStefano Zampini     if (sub_schurs->is_dir) {
440d62866d3SStefano Zampini       ierr = ISGetLocalSize(sub_schurs->is_dir,&n_dir);CHKERRQ(ierr);
4411b968477SStefano Zampini     }
44272b8c272SStefano Zampini     if (!deluxe_ctx->n_simple) {
4431b968477SStefano Zampini       deluxe_ctx->n_simple = n_dir + n_com;
4441b968477SStefano Zampini       ierr = PetscMalloc1(deluxe_ctx->n_simple,&deluxe_ctx->idx_simple_B);CHKERRQ(ierr);
445d62866d3SStefano Zampini       if (sub_schurs->is_vertices) {
4469bb4a8caSStefano Zampini         PetscInt       nmap;
4479bb4a8caSStefano Zampini         const PetscInt *idxs;
4481b968477SStefano Zampini 
449d62866d3SStefano Zampini         ierr = ISGetIndices(sub_schurs->is_vertices,&idxs);CHKERRQ(ierr);
4501b968477SStefano Zampini         ierr = ISGlobalToLocalMappingApply(pcis->BtoNmap,IS_GTOLM_DROP,n_com,idxs,&nmap,deluxe_ctx->idx_simple_B);CHKERRQ(ierr);
451af25d912SStefano Zampini         if (nmap != n_com) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Error when mapping simply scaled dofs (is_vertices)! %d != %d",nmap,n_com);
452d62866d3SStefano Zampini         ierr = ISRestoreIndices(sub_schurs->is_vertices,&idxs);CHKERRQ(ierr);
4531b968477SStefano Zampini       }
454d62866d3SStefano Zampini       if (sub_schurs->is_dir) {
4551b968477SStefano Zampini         PetscInt       nmap;
4561b968477SStefano Zampini         const PetscInt *idxs;
4571b968477SStefano Zampini 
458d62866d3SStefano Zampini         ierr = ISGetIndices(sub_schurs->is_dir,&idxs);CHKERRQ(ierr);
4591b968477SStefano Zampini         ierr = ISGlobalToLocalMappingApply(pcis->BtoNmap,IS_GTOLM_DROP,n_dir,idxs,&nmap,deluxe_ctx->idx_simple_B+n_com);CHKERRQ(ierr);
460af25d912SStefano Zampini         if (nmap != n_dir) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Error when mapping simply scaled dofs (sub_schurs->is_dir)! %d != %d",nmap,n_dir);
461d62866d3SStefano Zampini         ierr = ISRestoreIndices(sub_schurs->is_dir,&idxs);CHKERRQ(ierr);
4621b968477SStefano Zampini       }
4631b968477SStefano Zampini       ierr = PetscSortInt(deluxe_ctx->n_simple,deluxe_ctx->idx_simple_B);CHKERRQ(ierr);
4649bb4a8caSStefano Zampini     } else {
465af25d912SStefano Zampini       if (deluxe_ctx->n_simple != n_dir + n_com) 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);
46672b8c272SStefano Zampini     }
46772b8c272SStefano Zampini   } else {
468b1b3d7a2SStefano Zampini     deluxe_ctx->n_simple = 0;
4699bb4a8caSStefano Zampini     deluxe_ctx->idx_simple_B = 0;
470b1b3d7a2SStefano Zampini   }
47134a97f8cSStefano Zampini   PetscFunctionReturn(0);
47234a97f8cSStefano Zampini }
47334a97f8cSStefano Zampini 
4745a95e1ceSStefano Zampini static PetscErrorCode PCBDDCScalingSetUp_Deluxe_Private(PC pc)
4755db18549SStefano Zampini {
4765db18549SStefano Zampini   PC_BDDC                *pcbddc=(PC_BDDC*)pc->data;
4775db18549SStefano Zampini   PCBDDCDeluxeScaling    deluxe_ctx=pcbddc->deluxe_ctx;
478b96c3477SStefano Zampini   PCBDDCSubSchurs        sub_schurs = pcbddc->sub_schurs;
47972b8c272SStefano Zampini   PetscScalar            *matdata,*matdata2;
48072b8c272SStefano Zampini   PetscInt               i,max_subset_size,cum,cum2;
48172b8c272SStefano Zampini   const PetscInt         *idxs;
48272b8c272SStefano Zampini   PetscBool              newsetup = PETSC_FALSE;
4835db18549SStefano Zampini   PetscErrorCode         ierr;
4845db18549SStefano Zampini 
4855db18549SStefano Zampini   PetscFunctionBegin;
4865a95e1ceSStefano Zampini   if (!sub_schurs->n_subs) {
4879221af80SStefano Zampini     PetscFunctionReturn(0);
4889221af80SStefano Zampini   }
4899221af80SStefano Zampini 
49072b8c272SStefano Zampini   /* Allocate arrays for subproblems */
49172b8c272SStefano Zampini   if (!deluxe_ctx->seq_n) {
49272b8c272SStefano Zampini     deluxe_ctx->seq_n = sub_schurs->n_subs;
49372b8c272SStefano 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);
49472b8c272SStefano Zampini     newsetup = PETSC_TRUE;
49572b8c272SStefano Zampini   } else if (deluxe_ctx->seq_n != sub_schurs->n_subs) {
49672b8c272SStefano 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);
49772b8c272SStefano Zampini   }
4982f41f7d1SStefano Zampini   /* the change of basis is just a reference to sub_schurs->change (if any) */
4992f41f7d1SStefano Zampini   deluxe_ctx->change         = sub_schurs->change;
50072b8c272SStefano Zampini   deluxe_ctx->change_with_qr = sub_schurs->change_with_qr;
5015db18549SStefano Zampini 
50272b8c272SStefano Zampini   /* Create objects for deluxe */
50372b8c272SStefano Zampini   max_subset_size = 0;
50472b8c272SStefano Zampini   for (i=0;i<sub_schurs->n_subs;i++) {
50572b8c272SStefano Zampini     PetscInt subset_size;
50672b8c272SStefano Zampini     ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr);
50772b8c272SStefano Zampini     max_subset_size = PetscMax(subset_size,max_subset_size);
50872b8c272SStefano Zampini   }
50972b8c272SStefano Zampini   if (newsetup) {
51072b8c272SStefano Zampini     ierr = PetscMalloc1(2*max_subset_size,&deluxe_ctx->workspace);CHKERRQ(ierr);
51172b8c272SStefano Zampini   }
51272b8c272SStefano Zampini   cum = cum2 = 0;
51372b8c272SStefano Zampini   ierr = ISGetIndices(sub_schurs->is_Ej_all,&idxs);CHKERRQ(ierr);
51472b8c272SStefano Zampini   ierr = MatSeqAIJGetArray(sub_schurs->S_Ej_all,&matdata);CHKERRQ(ierr);
51572b8c272SStefano Zampini   ierr = MatSeqAIJGetArray(sub_schurs->sum_S_Ej_all,&matdata2);CHKERRQ(ierr);
51672b8c272SStefano Zampini   for (i=0;i<deluxe_ctx->seq_n;i++) {
51772b8c272SStefano Zampini     PetscInt     subset_size;
518925dfd53SStefano Zampini 
51972b8c272SStefano Zampini     ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr);
52072b8c272SStefano Zampini     if (newsetup) {
52172b8c272SStefano Zampini       IS  sub;
52272b8c272SStefano Zampini       /* work vectors */
52372b8c272SStefano Zampini       ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,1,subset_size,deluxe_ctx->workspace,&deluxe_ctx->seq_work1[i]);CHKERRQ(ierr);
52472b8c272SStefano Zampini       ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,1,subset_size,deluxe_ctx->workspace+subset_size,&deluxe_ctx->seq_work2[i]);CHKERRQ(ierr);
52572b8c272SStefano Zampini 
52672b8c272SStefano Zampini       /* scatters */
52772b8c272SStefano Zampini       ierr = ISCreateGeneral(PETSC_COMM_SELF,subset_size,idxs+cum,PETSC_COPY_VALUES,&sub);CHKERRQ(ierr);
52872b8c272SStefano Zampini       ierr = VecScatterCreate(pcbddc->work_scaling,sub,deluxe_ctx->seq_work1[i],NULL,&deluxe_ctx->seq_scctx[i]);CHKERRQ(ierr);
52972b8c272SStefano Zampini       ierr = ISDestroy(&sub);CHKERRQ(ierr);
530839e9adbSstefano_zampini     }
53172b8c272SStefano Zampini 
53272b8c272SStefano Zampini     /* S_E_j */
533839e9adbSstefano_zampini     ierr = MatDestroy(&deluxe_ctx->seq_mat[i]);CHKERRQ(ierr);
53472b8c272SStefano Zampini     ierr = MatCreateSeqDense(PETSC_COMM_SELF,subset_size,subset_size,matdata+cum2,&deluxe_ctx->seq_mat[i]);CHKERRQ(ierr);
535839e9adbSstefano_zampini 
53672b8c272SStefano Zampini     /* \sum_k S^k_E_j */
53772b8c272SStefano Zampini     ierr = MatDestroy(&deluxe_ctx->seq_mat_inv_sum[i]);CHKERRQ(ierr);
53872b8c272SStefano Zampini     ierr = MatCreateSeqDense(PETSC_COMM_SELF,subset_size,subset_size,matdata2+cum2,&deluxe_ctx->seq_mat_inv_sum[i]);CHKERRQ(ierr);
539729c86d3Sstefano_zampini 
54072b8c272SStefano Zampini     if (sub_schurs->is_hermitian && sub_schurs->is_posdef) {
54172b8c272SStefano Zampini       ierr = MatCholeskyFactor(deluxe_ctx->seq_mat_inv_sum[i],NULL,NULL);CHKERRQ(ierr);
542d62866d3SStefano Zampini     } else {
54372b8c272SStefano Zampini       ierr = MatLUFactor(deluxe_ctx->seq_mat_inv_sum[i],NULL,NULL,NULL);CHKERRQ(ierr);
544d62866d3SStefano Zampini     }
545839e9adbSstefano_zampini     if (pcbddc->deluxe_singlemat) {
546839e9adbSstefano_zampini       Mat X,Y;
547839e9adbSstefano_zampini       if (!sub_schurs->is_hermitian || !sub_schurs->is_posdef) {
548839e9adbSstefano_zampini         ierr = MatTranspose(deluxe_ctx->seq_mat[i],MAT_INITIAL_MATRIX,&X);CHKERRQ(ierr);
549839e9adbSstefano_zampini       } else {
550839e9adbSstefano_zampini         ierr = PetscObjectReference((PetscObject)deluxe_ctx->seq_mat[i]);CHKERRQ(ierr);
551839e9adbSstefano_zampini         X    = deluxe_ctx->seq_mat[i];
552839e9adbSstefano_zampini       }
553839e9adbSstefano_zampini       ierr = MatDuplicate(X,MAT_DO_NOT_COPY_VALUES,&Y);CHKERRQ(ierr);
554839e9adbSstefano_zampini       if (!sub_schurs->is_hermitian || !sub_schurs->is_posdef) {
555839e9adbSstefano_zampini         ierr = PCBDDCMatTransposeMatSolve_SeqDense(deluxe_ctx->seq_mat_inv_sum[i],X,Y);CHKERRQ(ierr);
556839e9adbSstefano_zampini       } else {
557839e9adbSstefano_zampini         ierr = MatMatSolve(deluxe_ctx->seq_mat_inv_sum[i],X,Y);CHKERRQ(ierr);
558839e9adbSstefano_zampini       }
559729c86d3Sstefano_zampini 
560839e9adbSstefano_zampini       ierr = MatDestroy(&deluxe_ctx->seq_mat_inv_sum[i]);CHKERRQ(ierr);
561839e9adbSstefano_zampini       ierr = MatDestroy(&deluxe_ctx->seq_mat[i]);CHKERRQ(ierr);
562839e9adbSstefano_zampini       ierr = MatDestroy(&X);CHKERRQ(ierr);
5632f41f7d1SStefano Zampini       if (deluxe_ctx->change) {
5642f41f7d1SStefano Zampini         Mat C,CY;
5652f41f7d1SStefano Zampini 
5662f41f7d1SStefano Zampini         if (!deluxe_ctx->change_with_qr) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Only QR based change of basis");
5672f41f7d1SStefano Zampini         ierr = KSPGetOperators(deluxe_ctx->change[i],&C,NULL);CHKERRQ(ierr);
5682f41f7d1SStefano Zampini         ierr = MatMatMult(C,Y,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&CY);CHKERRQ(ierr);
5692f41f7d1SStefano Zampini         ierr = MatMatTransposeMult(CY,C,MAT_REUSE_MATRIX,PETSC_DEFAULT,&Y);CHKERRQ(ierr);
5702f41f7d1SStefano Zampini         ierr = MatDestroy(&CY);CHKERRQ(ierr);
5712f41f7d1SStefano Zampini       }
572*03dfb2d7SStefano Zampini       ierr = MatTranspose(Y,MAT_INPLACE_MATRIX,&Y);CHKERRQ(ierr);
573839e9adbSstefano_zampini       deluxe_ctx->seq_mat[i] = Y;
57472b8c272SStefano Zampini     }
57572b8c272SStefano Zampini     cum += subset_size;
57672b8c272SStefano Zampini     cum2 += subset_size*subset_size;
57772b8c272SStefano Zampini   }
57872b8c272SStefano Zampini   ierr = ISRestoreIndices(sub_schurs->is_Ej_all,&idxs);CHKERRQ(ierr);
57972b8c272SStefano Zampini   ierr = MatSeqAIJRestoreArray(sub_schurs->S_Ej_all,&matdata);CHKERRQ(ierr);
58072b8c272SStefano Zampini   ierr = MatSeqAIJRestoreArray(sub_schurs->sum_S_Ej_all,&matdata2);CHKERRQ(ierr);
5812f41f7d1SStefano Zampini   if (pcbddc->deluxe_singlemat) {
5822f41f7d1SStefano Zampini     deluxe_ctx->change         = NULL;
5832f41f7d1SStefano Zampini     deluxe_ctx->change_with_qr = PETSC_FALSE;
5842f41f7d1SStefano Zampini   }
5855db18549SStefano Zampini 
58672b8c272SStefano Zampini   if (deluxe_ctx->change && !deluxe_ctx->change_with_qr) {
58772b8c272SStefano Zampini     for (i=0;i<deluxe_ctx->seq_n;i++) {
58872b8c272SStefano Zampini       if (newsetup) {
58972b8c272SStefano Zampini         PC pc;
59072b8c272SStefano Zampini 
59172b8c272SStefano Zampini         ierr = KSPGetPC(deluxe_ctx->change[i],&pc);CHKERRQ(ierr);
59272b8c272SStefano Zampini         ierr = PCSetType(pc,PCLU);CHKERRQ(ierr);
59372b8c272SStefano Zampini         ierr = KSPSetFromOptions(deluxe_ctx->change[i]);CHKERRQ(ierr);
59472b8c272SStefano Zampini       }
59572b8c272SStefano Zampini       ierr = KSPSetUp(deluxe_ctx->change[i]);CHKERRQ(ierr);
59672b8c272SStefano Zampini     }
59716e386b8SStefano Zampini   }
5985db18549SStefano Zampini   PetscFunctionReturn(0);
5995db18549SStefano Zampini }
600