xref: /petsc/src/ksp/pc/impls/bddc/bddcscalingbasic.c (revision 729c86d3a64c57e0da2e90e7cd2fc81b1248e478)
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 
13674ae819SStefano Zampini #undef __FUNCT__
14839e9adbSstefano_zampini #define __FUNCT__ "PCBDDCMatTransposeMatSolve_SeqDense"
15839e9adbSstefano_zampini static PetscErrorCode PCBDDCMatTransposeMatSolve_SeqDense(Mat A,Mat B,Mat X)
16839e9adbSstefano_zampini {
17839e9adbSstefano_zampini   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
18839e9adbSstefano_zampini   PetscErrorCode ierr;
19839e9adbSstefano_zampini   PetscScalar    *b,*x;
20839e9adbSstefano_zampini   PetscInt       n;
21839e9adbSstefano_zampini   PetscBLASInt   nrhs,info,m;
22839e9adbSstefano_zampini   PetscBool      flg;
23839e9adbSstefano_zampini 
24839e9adbSstefano_zampini   PetscFunctionBegin;
25839e9adbSstefano_zampini   ierr = PetscBLASIntCast(A->rmap->n,&m);CHKERRQ(ierr);
26839e9adbSstefano_zampini   ierr = PetscObjectTypeCompareAny((PetscObject)B,&flg,MATSEQDENSE,MATMPIDENSE,NULL);CHKERRQ(ierr);
27839e9adbSstefano_zampini   if (!flg) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONG,"Matrix B must be MATDENSE matrix");
28839e9adbSstefano_zampini   ierr = PetscObjectTypeCompareAny((PetscObject)X,&flg,MATSEQDENSE,MATMPIDENSE,NULL);CHKERRQ(ierr);
29839e9adbSstefano_zampini   if (!flg) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONG,"Matrix X must be MATDENSE matrix");
30839e9adbSstefano_zampini 
31839e9adbSstefano_zampini   ierr = MatGetSize(B,NULL,&n);CHKERRQ(ierr);
32839e9adbSstefano_zampini   ierr = PetscBLASIntCast(n,&nrhs);CHKERRQ(ierr);
33839e9adbSstefano_zampini   ierr = MatDenseGetArray(B,&b);CHKERRQ(ierr);
34839e9adbSstefano_zampini   ierr = MatDenseGetArray(X,&x);CHKERRQ(ierr);
35839e9adbSstefano_zampini 
36839e9adbSstefano_zampini   ierr = PetscMemcpy(x,b,m*nrhs*sizeof(PetscScalar));CHKERRQ(ierr);
37839e9adbSstefano_zampini 
38839e9adbSstefano_zampini   if (A->factortype == MAT_FACTOR_LU) {
39839e9adbSstefano_zampini #if defined(PETSC_MISSING_LAPACK_GETRS)
40839e9adbSstefano_zampini     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"GETRS - Lapack routine is unavailable.");
41839e9adbSstefano_zampini #else
42839e9adbSstefano_zampini     PetscStackCallBLAS("LAPACKgetrs",LAPACKgetrs_("T",&m,&nrhs,mat->v,&mat->lda,mat->pivots,x,&m,&info));
43839e9adbSstefano_zampini     if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"GETRS - Bad solve");
44839e9adbSstefano_zampini #endif
45839e9adbSstefano_zampini   } else if (A->factortype == MAT_FACTOR_CHOLESKY) {
46839e9adbSstefano_zampini #if defined(PETSC_MISSING_LAPACK_POTRS)
47839e9adbSstefano_zampini     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"POTRS - Lapack routine is unavailable.");
48839e9adbSstefano_zampini #else
49839e9adbSstefano_zampini     PetscStackCallBLAS("LAPACKpotrs",LAPACKpotrs_("L",&m,&nrhs,mat->v,&mat->lda,x,&m,&info));
50839e9adbSstefano_zampini     if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"POTRS Bad solve");
51839e9adbSstefano_zampini #endif
52839e9adbSstefano_zampini   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Matrix must be factored to solve");
53839e9adbSstefano_zampini 
54839e9adbSstefano_zampini   ierr = MatDenseRestoreArray(B,&b);CHKERRQ(ierr);
55839e9adbSstefano_zampini   ierr = MatDenseRestoreArray(X,&x);CHKERRQ(ierr);
56839e9adbSstefano_zampini   ierr = PetscLogFlops(nrhs*(2.0*m*m - m));CHKERRQ(ierr);
57839e9adbSstefano_zampini   PetscFunctionReturn(0);
58839e9adbSstefano_zampini }
59839e9adbSstefano_zampini 
60839e9adbSstefano_zampini 
61839e9adbSstefano_zampini #undef __FUNCT__
62674ae819SStefano Zampini #define __FUNCT__ "PCBDDCScalingExtension_Basic"
63674ae819SStefano Zampini static PetscErrorCode PCBDDCScalingExtension_Basic(PC pc, Vec local_interface_vector, Vec global_vector)
64674ae819SStefano Zampini {
65674ae819SStefano Zampini   PC_IS*         pcis = (PC_IS*)pc->data;
66674ae819SStefano Zampini   PC_BDDC*       pcbddc = (PC_BDDC*)pc->data;
67674ae819SStefano Zampini   PetscErrorCode ierr;
68674ae819SStefano Zampini 
69674ae819SStefano Zampini   PetscFunctionBegin;
70674ae819SStefano Zampini   /* Apply partition of unity */
71674ae819SStefano Zampini   ierr = VecPointwiseMult(pcbddc->work_scaling,pcis->D,local_interface_vector);CHKERRQ(ierr);
72674ae819SStefano Zampini   ierr = VecSet(global_vector,0.0);CHKERRQ(ierr);
73674ae819SStefano Zampini   ierr = VecScatterBegin(pcis->global_to_B,pcbddc->work_scaling,global_vector,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
74674ae819SStefano Zampini   ierr = VecScatterEnd(pcis->global_to_B,pcbddc->work_scaling,global_vector,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
75674ae819SStefano Zampini   PetscFunctionReturn(0);
76674ae819SStefano Zampini }
77674ae819SStefano Zampini 
78674ae819SStefano Zampini #undef __FUNCT__
79674ae819SStefano Zampini #define __FUNCT__ "PCBDDCScalingExtension_Deluxe"
80674ae819SStefano Zampini static PetscErrorCode PCBDDCScalingExtension_Deluxe(PC pc, Vec x, Vec y)
81674ae819SStefano Zampini {
82674ae819SStefano Zampini   PC_IS*              pcis=(PC_IS*)pc->data;
83674ae819SStefano Zampini   PC_BDDC*            pcbddc=(PC_BDDC*)pc->data;
84674ae819SStefano Zampini   PCBDDCDeluxeScaling deluxe_ctx = pcbddc->deluxe_ctx;
85674ae819SStefano Zampini   PetscErrorCode      ierr;
86674ae819SStefano Zampini 
87674ae819SStefano Zampini   PetscFunctionBegin;
885a95e1ceSStefano Zampini   ierr = VecSet(pcbddc->work_scaling,0.0);CHKERRQ(ierr);
895a95e1ceSStefano Zampini   ierr = VecSet(y,0.0);CHKERRQ(ierr);
905a95e1ceSStefano Zampini   if (deluxe_ctx->n_simple) { /* scale deluxe vertices using diagonal scaling */
915a95e1ceSStefano Zampini     PetscInt          i;
922b095fd8SStefano Zampini     const PetscScalar *array_x,*array_D;
932b095fd8SStefano Zampini     PetscScalar       *array;
942b095fd8SStefano Zampini     ierr = VecGetArrayRead(x,&array_x);CHKERRQ(ierr);
952b095fd8SStefano Zampini     ierr = VecGetArrayRead(pcis->D,&array_D);CHKERRQ(ierr);
96674ae819SStefano Zampini     ierr = VecGetArray(pcbddc->work_scaling,&array);CHKERRQ(ierr);
97674ae819SStefano Zampini     for (i=0;i<deluxe_ctx->n_simple;i++) {
98674ae819SStefano Zampini       array[deluxe_ctx->idx_simple_B[i]] = array_x[deluxe_ctx->idx_simple_B[i]]*array_D[deluxe_ctx->idx_simple_B[i]];
99674ae819SStefano Zampini     }
100674ae819SStefano Zampini     ierr = VecRestoreArray(pcbddc->work_scaling,&array);CHKERRQ(ierr);
1012b095fd8SStefano Zampini     ierr = VecRestoreArrayRead(pcis->D,&array_D);CHKERRQ(ierr);
1022b095fd8SStefano Zampini     ierr = VecRestoreArrayRead(x,&array_x);CHKERRQ(ierr);
10334a97f8cSStefano Zampini   }
104ac632422SStefano Zampini   /* sequential part : all problems and Schur applications collapsed into a single matrix vector multiplication or a matvec and a solve */
1055a95e1ceSStefano Zampini   if (deluxe_ctx->seq_mat) {
10672b8c272SStefano Zampini     PetscInt i;
10772b8c272SStefano Zampini     for (i=0;i<deluxe_ctx->seq_n;i++) {
10816e386b8SStefano Zampini       if (deluxe_ctx->change) {
10972b8c272SStefano Zampini         ierr = VecScatterBegin(deluxe_ctx->seq_scctx[i],x,deluxe_ctx->seq_work2[i],INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
11072b8c272SStefano Zampini         ierr = VecScatterEnd(deluxe_ctx->seq_scctx[i],x,deluxe_ctx->seq_work2[i],INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
11172b8c272SStefano Zampini         if (deluxe_ctx->change_with_qr) {
11272b8c272SStefano Zampini           Mat change;
11372b8c272SStefano Zampini 
11472b8c272SStefano Zampini           ierr = KSPGetOperators(deluxe_ctx->change[i],&change,NULL);CHKERRQ(ierr);
11572b8c272SStefano Zampini           ierr = MatMultTranspose(change,deluxe_ctx->seq_work2[i],deluxe_ctx->seq_work1[i]);CHKERRQ(ierr);
11672b8c272SStefano Zampini         } else {
11772b8c272SStefano Zampini           ierr = KSPSolve(deluxe_ctx->change[i],deluxe_ctx->seq_work2[i],deluxe_ctx->seq_work1[i]);CHKERRQ(ierr);
11816e386b8SStefano Zampini         }
11972b8c272SStefano Zampini       } else {
12072b8c272SStefano Zampini         ierr = VecScatterBegin(deluxe_ctx->seq_scctx[i],x,deluxe_ctx->seq_work1[i],INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
12172b8c272SStefano Zampini         ierr = VecScatterEnd(deluxe_ctx->seq_scctx[i],x,deluxe_ctx->seq_work1[i],INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
12272b8c272SStefano Zampini       }
12372b8c272SStefano Zampini       ierr = MatMultTranspose(deluxe_ctx->seq_mat[i],deluxe_ctx->seq_work1[i],deluxe_ctx->seq_work2[i]);CHKERRQ(ierr);
12472b8c272SStefano Zampini       if (deluxe_ctx->seq_mat_inv_sum[i]) {
12572b8c272SStefano Zampini         PetscScalar *x;
12672b8c272SStefano Zampini 
12772b8c272SStefano Zampini         ierr = VecGetArray(deluxe_ctx->seq_work2[i],&x);CHKERRQ(ierr);
12872b8c272SStefano Zampini         ierr = VecPlaceArray(deluxe_ctx->seq_work1[i],x);CHKERRQ(ierr);
12972b8c272SStefano Zampini         ierr = VecRestoreArray(deluxe_ctx->seq_work2[i],&x);CHKERRQ(ierr);
13072b8c272SStefano Zampini         ierr = MatSolveTranspose(deluxe_ctx->seq_mat_inv_sum[i],deluxe_ctx->seq_work1[i],deluxe_ctx->seq_work2[i]);CHKERRQ(ierr);
13172b8c272SStefano Zampini         ierr = VecResetArray(deluxe_ctx->seq_work1[i]);CHKERRQ(ierr);
132ac632422SStefano Zampini       }
13316e386b8SStefano Zampini       if (deluxe_ctx->change) {
13416e386b8SStefano Zampini         Mat change;
13516e386b8SStefano Zampini 
13672b8c272SStefano Zampini         ierr = KSPGetOperators(deluxe_ctx->change[i],&change,NULL);CHKERRQ(ierr);
13772b8c272SStefano Zampini         ierr = MatMult(change,deluxe_ctx->seq_work2[i],deluxe_ctx->seq_work1[i]);CHKERRQ(ierr);
13872b8c272SStefano Zampini         ierr = VecScatterBegin(deluxe_ctx->seq_scctx[i],deluxe_ctx->seq_work1[i],pcbddc->work_scaling,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
13972b8c272SStefano Zampini         ierr = VecScatterEnd(deluxe_ctx->seq_scctx[i],deluxe_ctx->seq_work1[i],pcbddc->work_scaling,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
14016e386b8SStefano Zampini       } else {
14172b8c272SStefano Zampini         ierr = VecScatterBegin(deluxe_ctx->seq_scctx[i],deluxe_ctx->seq_work2[i],pcbddc->work_scaling,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
14272b8c272SStefano Zampini         ierr = VecScatterEnd(deluxe_ctx->seq_scctx[i],deluxe_ctx->seq_work2[i],pcbddc->work_scaling,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
14372b8c272SStefano Zampini       }
144674ae819SStefano Zampini     }
14516e386b8SStefano Zampini   }
146674ae819SStefano Zampini   /* put local boundary part in global vector */
147674ae819SStefano Zampini   ierr = VecScatterBegin(pcis->global_to_B,pcbddc->work_scaling,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
148674ae819SStefano Zampini   ierr = VecScatterEnd(pcis->global_to_B,pcbddc->work_scaling,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
149674ae819SStefano Zampini   PetscFunctionReturn(0);
150674ae819SStefano Zampini }
151674ae819SStefano Zampini 
152674ae819SStefano Zampini #undef __FUNCT__
153674ae819SStefano Zampini #define __FUNCT__ "PCBDDCScalingExtension"
154674ae819SStefano Zampini PetscErrorCode PCBDDCScalingExtension(PC pc, Vec local_interface_vector, Vec global_vector)
155674ae819SStefano Zampini {
156674ae819SStefano Zampini   PC_BDDC        *pcbddc=(PC_BDDC*)pc->data;
157674ae819SStefano Zampini   PetscErrorCode ierr;
158674ae819SStefano Zampini 
159674ae819SStefano Zampini   PetscFunctionBegin;
160674ae819SStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
161674ae819SStefano Zampini   PetscValidHeaderSpecific(local_interface_vector,VEC_CLASSID,2);
162674ae819SStefano Zampini   PetscValidHeaderSpecific(global_vector,VEC_CLASSID,3);
1636c4ed002SBarry Smith   if (local_interface_vector == pcbddc->work_scaling) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Local vector cannot be pcbddc->work_scaling!\n");
164163d334eSBarry Smith   ierr = PetscUseMethod(pc,"PCBDDCScalingExtension_C",(PC,Vec,Vec),(pc,local_interface_vector,global_vector));CHKERRQ(ierr);
165674ae819SStefano Zampini   PetscFunctionReturn(0);
166674ae819SStefano Zampini }
167674ae819SStefano Zampini 
168674ae819SStefano Zampini #undef __FUNCT__
169674ae819SStefano Zampini #define __FUNCT__ "PCBDDCScalingRestriction_Basic"
170674ae819SStefano Zampini static PetscErrorCode PCBDDCScalingRestriction_Basic(PC pc, Vec global_vector, Vec local_interface_vector)
171674ae819SStefano Zampini {
172674ae819SStefano Zampini   PetscErrorCode ierr;
173674ae819SStefano Zampini   PC_IS          *pcis = (PC_IS*)pc->data;
174674ae819SStefano Zampini 
175674ae819SStefano Zampini   PetscFunctionBegin;
176674ae819SStefano Zampini   ierr = VecScatterBegin(pcis->global_to_B,global_vector,local_interface_vector,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
177674ae819SStefano Zampini   ierr = VecScatterEnd(pcis->global_to_B,global_vector,local_interface_vector,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
178674ae819SStefano Zampini   /* Apply partition of unity */
179674ae819SStefano Zampini   ierr = VecPointwiseMult(local_interface_vector,pcis->D,local_interface_vector);CHKERRQ(ierr);
180674ae819SStefano Zampini   PetscFunctionReturn(0);
181674ae819SStefano Zampini }
182674ae819SStefano Zampini 
183674ae819SStefano Zampini #undef __FUNCT__
184674ae819SStefano Zampini #define __FUNCT__ "PCBDDCScalingRestriction_Deluxe"
185674ae819SStefano Zampini static PetscErrorCode PCBDDCScalingRestriction_Deluxe(PC pc, Vec x, Vec y)
186674ae819SStefano Zampini {
187674ae819SStefano Zampini   PC_IS*              pcis=(PC_IS*)pc->data;
188674ae819SStefano Zampini   PC_BDDC*            pcbddc=(PC_BDDC*)pc->data;
189674ae819SStefano Zampini   PCBDDCDeluxeScaling deluxe_ctx = pcbddc->deluxe_ctx;
190674ae819SStefano Zampini   PetscErrorCode      ierr;
191674ae819SStefano Zampini 
192674ae819SStefano Zampini   PetscFunctionBegin;
193674ae819SStefano Zampini   /* get local boundary part of global vector */
194674ae819SStefano Zampini   ierr = VecScatterBegin(pcis->global_to_B,x,y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
195674ae819SStefano Zampini   ierr = VecScatterEnd(pcis->global_to_B,x,y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1965a95e1ceSStefano Zampini   if (deluxe_ctx->n_simple) { /* scale deluxe vertices using diagonal scaling */
1975a95e1ceSStefano Zampini     PetscInt          i;
1982b095fd8SStefano Zampini     PetscScalar       *array_y;
1992b095fd8SStefano Zampini     const PetscScalar *array_D;
200674ae819SStefano Zampini     ierr = VecGetArray(y,&array_y);CHKERRQ(ierr);
2012b095fd8SStefano Zampini     ierr = VecGetArrayRead(pcis->D,&array_D);CHKERRQ(ierr);
202674ae819SStefano Zampini     for (i=0;i<deluxe_ctx->n_simple;i++) {
203674ae819SStefano Zampini       array_y[deluxe_ctx->idx_simple_B[i]] *= array_D[deluxe_ctx->idx_simple_B[i]];
204674ae819SStefano Zampini     }
2052b095fd8SStefano Zampini     ierr = VecRestoreArrayRead(pcis->D,&array_D);CHKERRQ(ierr);
206674ae819SStefano Zampini     ierr = VecRestoreArray(y,&array_y);CHKERRQ(ierr);
20734a97f8cSStefano Zampini   }
208*729c86d3Sstefano_zampini   /* sequential part : all problems and Schur applications collapsed into a single matrix vector multiplication or a matvec and a solve */
2095a95e1ceSStefano Zampini   if (deluxe_ctx->seq_mat) {
21072b8c272SStefano Zampini     PetscInt i;
21172b8c272SStefano Zampini     for (i=0;i<deluxe_ctx->seq_n;i++) {
21216e386b8SStefano Zampini       if (deluxe_ctx->change) {
21316e386b8SStefano Zampini         Mat change;
21416e386b8SStefano Zampini 
21572b8c272SStefano Zampini         ierr = VecScatterBegin(deluxe_ctx->seq_scctx[i],y,deluxe_ctx->seq_work2[i],INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
21672b8c272SStefano Zampini         ierr = VecScatterEnd(deluxe_ctx->seq_scctx[i],y,deluxe_ctx->seq_work2[i],INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
21772b8c272SStefano Zampini         ierr = KSPGetOperators(deluxe_ctx->change[i],&change,NULL);CHKERRQ(ierr);
21872b8c272SStefano Zampini         ierr = MatMultTranspose(change,deluxe_ctx->seq_work2[i],deluxe_ctx->seq_work1[i]);CHKERRQ(ierr);
21916e386b8SStefano Zampini       } else {
22072b8c272SStefano Zampini         ierr = VecScatterBegin(deluxe_ctx->seq_scctx[i],y,deluxe_ctx->seq_work1[i],INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
22172b8c272SStefano Zampini         ierr = VecScatterEnd(deluxe_ctx->seq_scctx[i],y,deluxe_ctx->seq_work1[i],INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
22216e386b8SStefano Zampini       }
22372b8c272SStefano Zampini       if (deluxe_ctx->seq_mat_inv_sum[i]) {
22472b8c272SStefano Zampini         PetscScalar *x;
22572b8c272SStefano Zampini 
22672b8c272SStefano Zampini         ierr = VecGetArray(deluxe_ctx->seq_work1[i],&x);CHKERRQ(ierr);
22772b8c272SStefano Zampini         ierr = VecPlaceArray(deluxe_ctx->seq_work2[i],x);CHKERRQ(ierr);
22872b8c272SStefano Zampini         ierr = VecRestoreArray(deluxe_ctx->seq_work1[i],&x);CHKERRQ(ierr);
22972b8c272SStefano Zampini         ierr = MatSolve(deluxe_ctx->seq_mat_inv_sum[i],deluxe_ctx->seq_work1[i],deluxe_ctx->seq_work2[i]);CHKERRQ(ierr);
23072b8c272SStefano Zampini         ierr = VecResetArray(deluxe_ctx->seq_work2[i]);CHKERRQ(ierr);
231ac632422SStefano Zampini       }
23272b8c272SStefano Zampini       ierr = MatMult(deluxe_ctx->seq_mat[i],deluxe_ctx->seq_work1[i],deluxe_ctx->seq_work2[i]);CHKERRQ(ierr);
23316e386b8SStefano Zampini       if (deluxe_ctx->change) {
23472b8c272SStefano Zampini         if (deluxe_ctx->change_with_qr) {
23572b8c272SStefano Zampini           Mat change;
23672b8c272SStefano Zampini 
23772b8c272SStefano Zampini           ierr = KSPGetOperators(deluxe_ctx->change[i],&change,NULL);CHKERRQ(ierr);
23872b8c272SStefano Zampini           ierr = MatMult(change,deluxe_ctx->seq_work2[i],deluxe_ctx->seq_work1[i]);CHKERRQ(ierr);
23972b8c272SStefano Zampini         } else {
24072b8c272SStefano Zampini           ierr = KSPSolveTranspose(deluxe_ctx->change[i],deluxe_ctx->seq_work2[i],deluxe_ctx->seq_work1[i]);CHKERRQ(ierr);
24116e386b8SStefano Zampini         }
24272b8c272SStefano Zampini         ierr = VecScatterBegin(deluxe_ctx->seq_scctx[i],deluxe_ctx->seq_work1[i],y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
24372b8c272SStefano Zampini         ierr = VecScatterEnd(deluxe_ctx->seq_scctx[i],deluxe_ctx->seq_work1[i],y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
24472b8c272SStefano Zampini       } else {
24572b8c272SStefano Zampini         ierr = VecScatterBegin(deluxe_ctx->seq_scctx[i],deluxe_ctx->seq_work2[i],y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
24672b8c272SStefano Zampini         ierr = VecScatterEnd(deluxe_ctx->seq_scctx[i],deluxe_ctx->seq_work2[i],y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
24772b8c272SStefano Zampini       }
24872b8c272SStefano Zampini     }
249674ae819SStefano Zampini   }
250674ae819SStefano Zampini   PetscFunctionReturn(0);
251674ae819SStefano Zampini }
252674ae819SStefano Zampini 
253674ae819SStefano Zampini #undef __FUNCT__
254674ae819SStefano Zampini #define __FUNCT__ "PCBDDCScalingRestriction"
255674ae819SStefano Zampini PetscErrorCode PCBDDCScalingRestriction(PC pc, Vec global_vector, Vec local_interface_vector)
256674ae819SStefano Zampini {
257674ae819SStefano Zampini   PC_BDDC        *pcbddc=(PC_BDDC*)pc->data;
258674ae819SStefano Zampini   PetscErrorCode ierr;
259674ae819SStefano Zampini 
260674ae819SStefano Zampini   PetscFunctionBegin;
261674ae819SStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
262674ae819SStefano Zampini   PetscValidHeaderSpecific(global_vector,VEC_CLASSID,2);
263674ae819SStefano Zampini   PetscValidHeaderSpecific(local_interface_vector,VEC_CLASSID,3);
2646c4ed002SBarry Smith   if (local_interface_vector == pcbddc->work_scaling) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Local vector cannot be pcbddc->work_scaling!\n");
265163d334eSBarry Smith   ierr = PetscUseMethod(pc,"PCBDDCScalingRestriction_C",(PC,Vec,Vec),(pc,global_vector,local_interface_vector));CHKERRQ(ierr);
266674ae819SStefano Zampini   PetscFunctionReturn(0);
267674ae819SStefano Zampini }
268674ae819SStefano Zampini 
269674ae819SStefano Zampini #undef __FUNCT__
270674ae819SStefano Zampini #define __FUNCT__ "PCBDDCScalingSetUp"
271674ae819SStefano Zampini PetscErrorCode PCBDDCScalingSetUp(PC pc)
272674ae819SStefano Zampini {
273674ae819SStefano Zampini   PC_IS*         pcis=(PC_IS*)pc->data;
274674ae819SStefano Zampini   PC_BDDC*       pcbddc=(PC_BDDC*)pc->data;
275674ae819SStefano Zampini   PetscErrorCode ierr;
276674ae819SStefano Zampini 
277674ae819SStefano Zampini   PetscFunctionBegin;
278674ae819SStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
279674ae819SStefano Zampini   /* create work vector for the operator */
28034a97f8cSStefano Zampini   ierr = VecDestroy(&pcbddc->work_scaling);CHKERRQ(ierr);
281674ae819SStefano Zampini   ierr = VecDuplicate(pcis->vec1_B,&pcbddc->work_scaling);CHKERRQ(ierr);
28234a97f8cSStefano Zampini   /* always rebuild pcis->D */
28328d874f6SStefano Zampini   if (pcis->use_stiffness_scaling) {
284674ae819SStefano Zampini     ierr = MatGetDiagonal(pcbddc->local_mat,pcis->vec1_N);CHKERRQ(ierr);
285674ae819SStefano Zampini     ierr = VecScatterBegin(pcis->N_to_B,pcis->vec1_N,pcis->D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
286674ae819SStefano Zampini     ierr = VecScatterEnd(pcis->N_to_B,pcis->vec1_N,pcis->D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
287674ae819SStefano Zampini   }
288674ae819SStefano Zampini   ierr = VecCopy(pcis->D,pcis->vec1_B);CHKERRQ(ierr);
289674ae819SStefano Zampini   ierr = VecSet(pcis->vec1_global,0.0);CHKERRQ(ierr);
290674ae819SStefano Zampini   ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_B,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
291674ae819SStefano Zampini   ierr = VecScatterEnd(pcis->global_to_B,pcis->vec1_B,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
292674ae819SStefano Zampini   ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_global,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
293674ae819SStefano Zampini   ierr = VecScatterEnd(pcis->global_to_B,pcis->vec1_global,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
294674ae819SStefano Zampini   ierr = VecPointwiseDivide(pcis->D,pcis->D,pcis->vec1_B);CHKERRQ(ierr);
295674ae819SStefano Zampini   /* now setup */
296681e7c04SStefano Zampini   if (pcbddc->use_deluxe_scaling) {
29734a97f8cSStefano Zampini     if (!pcbddc->deluxe_ctx) {
29834a97f8cSStefano Zampini       ierr = PCBDDCScalingCreate_Deluxe(pc);CHKERRQ(ierr);
29934a97f8cSStefano Zampini     }
30034a97f8cSStefano Zampini     ierr = PCBDDCScalingSetUp_Deluxe(pc);CHKERRQ(ierr);
301674ae819SStefano Zampini     ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCScalingRestriction_C",PCBDDCScalingRestriction_Deluxe);CHKERRQ(ierr);
302674ae819SStefano Zampini     ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCScalingExtension_C",PCBDDCScalingExtension_Deluxe);CHKERRQ(ierr);
303674ae819SStefano Zampini   } else {
304674ae819SStefano Zampini     ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCScalingRestriction_C",PCBDDCScalingRestriction_Basic);CHKERRQ(ierr);
305674ae819SStefano Zampini     ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCScalingExtension_C",PCBDDCScalingExtension_Basic);CHKERRQ(ierr);
306674ae819SStefano Zampini   }
30734a97f8cSStefano Zampini 
308674ae819SStefano Zampini   /* test */
309674ae819SStefano Zampini   if (pcbddc->dbg_flag) {
31072b8c272SStefano Zampini     Mat         B0_B = NULL;
31172b8c272SStefano Zampini     Vec         B0_Bv = NULL, B0_Bv2 = NULL;
31234a97f8cSStefano Zampini     Vec         vec2_global;
313674ae819SStefano Zampini     PetscViewer viewer = pcbddc->dbg_viewer;
31434a97f8cSStefano Zampini     PetscReal   error;
315674ae819SStefano Zampini 
316674ae819SStefano Zampini     /* extension -> from local to parallel */
31734a97f8cSStefano Zampini     ierr = VecSet(pcis->vec1_global,0.0);CHKERRQ(ierr);
31834a97f8cSStefano Zampini     ierr = VecSetRandom(pcis->vec1_B,NULL);CHKERRQ(ierr);
31934a97f8cSStefano Zampini     ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_B,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
32034a97f8cSStefano Zampini     ierr = VecScatterEnd(pcis->global_to_B,pcis->vec1_B,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
32134a97f8cSStefano Zampini     ierr = VecDuplicate(pcis->vec1_global,&vec2_global);CHKERRQ(ierr);
32234a97f8cSStefano Zampini     ierr = VecCopy(pcis->vec1_global,vec2_global);CHKERRQ(ierr);
323674ae819SStefano Zampini     ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_global,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
324674ae819SStefano Zampini     ierr = VecScatterEnd(pcis->global_to_B,pcis->vec1_global,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
32572b8c272SStefano Zampini     if (pcbddc->benign_n) {
32672b8c272SStefano Zampini       IS is_dummy;
32772b8c272SStefano Zampini 
32872b8c272SStefano Zampini       ierr = ISCreateStride(PETSC_COMM_SELF,pcbddc->benign_n,0,1,&is_dummy);CHKERRQ(ierr);
32972b8c272SStefano Zampini       ierr = MatGetSubMatrix(pcbddc->benign_B0,is_dummy,pcis->is_B_local,MAT_INITIAL_MATRIX,&B0_B);CHKERRQ(ierr);
33072b8c272SStefano Zampini       ierr = ISDestroy(&is_dummy);CHKERRQ(ierr);
33172b8c272SStefano Zampini       ierr = MatCreateVecs(B0_B,NULL,&B0_Bv);CHKERRQ(ierr);
33272b8c272SStefano Zampini       ierr = VecDuplicate(B0_Bv,&B0_Bv2);CHKERRQ(ierr);
33372b8c272SStefano Zampini       ierr = MatMult(B0_B,pcis->vec1_B,B0_Bv);CHKERRQ(ierr);
33472b8c272SStefano Zampini     }
335674ae819SStefano Zampini     ierr = PCBDDCScalingExtension(pc,pcis->vec1_B,pcis->vec1_global);CHKERRQ(ierr);
33672b8c272SStefano Zampini     if (pcbddc->benign_saddle_point) {
33772b8c272SStefano Zampini       PetscReal errorl = 0.;
33872b8c272SStefano Zampini       ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_global,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
33972b8c272SStefano Zampini       ierr = VecScatterEnd(pcis->global_to_B,pcis->vec1_global,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
34072b8c272SStefano Zampini       if (pcbddc->benign_n) {
34172b8c272SStefano Zampini         ierr = MatMult(B0_B,pcis->vec1_B,B0_Bv2);CHKERRQ(ierr);
34272b8c272SStefano Zampini         ierr = VecAXPY(B0_Bv,-1.0,B0_Bv2);CHKERRQ(ierr);
34372b8c272SStefano Zampini         ierr = VecNorm(B0_Bv,NORM_INFINITY,&errorl);CHKERRQ(ierr);
34472b8c272SStefano Zampini       }
34572b8c272SStefano Zampini       ierr = MPI_Allreduce(&errorl,&error,1,MPIU_REAL,MPI_SUM,PetscObjectComm((PetscObject)pc));CHKERRQ(ierr);
34672b8c272SStefano Zampini       ierr = PetscViewerASCIIPrintf(viewer,"Error benign extension %1.14e\n",error);CHKERRQ(ierr);
34772b8c272SStefano Zampini     }
34834a97f8cSStefano Zampini     ierr = VecAXPY(pcis->vec1_global,-1.0,vec2_global);CHKERRQ(ierr);
34934a97f8cSStefano Zampini     ierr = VecNorm(pcis->vec1_global,NORM_INFINITY,&error);CHKERRQ(ierr);
350674ae819SStefano Zampini     ierr = PetscViewerASCIIPrintf(viewer,"Error scaling extension %1.14e\n",error);CHKERRQ(ierr);
35134a97f8cSStefano Zampini     ierr = VecDestroy(&vec2_global);CHKERRQ(ierr);
35234a97f8cSStefano Zampini 
353674ae819SStefano Zampini     /* restriction -> from parallel to local */
354674ae819SStefano Zampini     ierr = VecSet(pcis->vec1_global,0.0);CHKERRQ(ierr);
35534a97f8cSStefano Zampini     ierr = VecSetRandom(pcis->vec1_B,NULL);CHKERRQ(ierr);
356674ae819SStefano Zampini     ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_B,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
357674ae819SStefano Zampini     ierr = VecScatterEnd(pcis->global_to_B,pcis->vec1_B,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
35834a97f8cSStefano Zampini     ierr = PCBDDCScalingRestriction(pc,pcis->vec1_global,pcis->vec1_B);CHKERRQ(ierr);
35934a97f8cSStefano Zampini     ierr = VecScale(pcis->vec1_B,-1.0);CHKERRQ(ierr);
36034a97f8cSStefano Zampini     ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_B,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
36134a97f8cSStefano Zampini     ierr = VecScatterEnd(pcis->global_to_B,pcis->vec1_B,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
36234a97f8cSStefano Zampini     ierr = VecNorm(pcis->vec1_global,NORM_INFINITY,&error);CHKERRQ(ierr);
36334a97f8cSStefano Zampini     ierr = PetscViewerASCIIPrintf(viewer,"Error scaling restriction %1.14e\n",error);CHKERRQ(ierr);
36472b8c272SStefano Zampini     ierr = MatDestroy(&B0_B);CHKERRQ(ierr);
36572b8c272SStefano Zampini     ierr = VecDestroy(&B0_Bv);CHKERRQ(ierr);
36672b8c272SStefano Zampini     ierr = VecDestroy(&B0_Bv2);CHKERRQ(ierr);
367674ae819SStefano Zampini   }
368674ae819SStefano Zampini   PetscFunctionReturn(0);
369674ae819SStefano Zampini }
370674ae819SStefano Zampini 
371674ae819SStefano Zampini #undef __FUNCT__
372674ae819SStefano Zampini #define __FUNCT__ "PCBDDCScalingDestroy"
373674ae819SStefano Zampini PetscErrorCode PCBDDCScalingDestroy(PC pc)
374674ae819SStefano Zampini {
375674ae819SStefano Zampini   PC_BDDC*       pcbddc=(PC_BDDC*)pc->data;
376674ae819SStefano Zampini   PetscErrorCode ierr;
377674ae819SStefano Zampini 
378674ae819SStefano Zampini   PetscFunctionBegin;
37934a97f8cSStefano Zampini   if (pcbddc->deluxe_ctx) {
38034a97f8cSStefano Zampini     ierr = PCBDDCScalingDestroy_Deluxe(pc);CHKERRQ(ierr);
381674ae819SStefano Zampini   }
382674ae819SStefano Zampini   ierr = VecDestroy(&pcbddc->work_scaling);CHKERRQ(ierr);
383674ae819SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCScalingRestriction_C",NULL);CHKERRQ(ierr);
384674ae819SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCScalingExtension_C",NULL);CHKERRQ(ierr);
385674ae819SStefano Zampini   PetscFunctionReturn(0);
386674ae819SStefano Zampini }
387674ae819SStefano Zampini 
38834a97f8cSStefano Zampini #undef __FUNCT__
38934a97f8cSStefano Zampini #define __FUNCT__ "PCBDDCScalingCreate_Deluxe"
39034a97f8cSStefano Zampini static PetscErrorCode PCBDDCScalingCreate_Deluxe(PC pc)
39134a97f8cSStefano Zampini {
39234a97f8cSStefano Zampini   PC_BDDC*            pcbddc=(PC_BDDC*)pc->data;
39334a97f8cSStefano Zampini   PCBDDCDeluxeScaling deluxe_ctx;
39434a97f8cSStefano Zampini   PetscErrorCode      ierr;
39534a97f8cSStefano Zampini 
39634a97f8cSStefano Zampini   PetscFunctionBegin;
39734a97f8cSStefano Zampini   ierr = PetscNew(&deluxe_ctx);CHKERRQ(ierr);
39834a97f8cSStefano Zampini   pcbddc->deluxe_ctx = deluxe_ctx;
39934a97f8cSStefano Zampini   PetscFunctionReturn(0);
40034a97f8cSStefano Zampini }
40134a97f8cSStefano Zampini 
40234a97f8cSStefano Zampini #undef __FUNCT__
40334a97f8cSStefano Zampini #define __FUNCT__ "PCBDDCScalingDestroy_Deluxe"
40434a97f8cSStefano Zampini static PetscErrorCode PCBDDCScalingDestroy_Deluxe(PC pc)
40534a97f8cSStefano Zampini {
40634a97f8cSStefano Zampini   PC_BDDC*            pcbddc=(PC_BDDC*)pc->data;
40734a97f8cSStefano Zampini   PetscErrorCode      ierr;
40834a97f8cSStefano Zampini 
40934a97f8cSStefano Zampini   PetscFunctionBegin;
41034a97f8cSStefano Zampini   ierr = PCBDDCScalingReset_Deluxe_Solvers(pcbddc->deluxe_ctx);CHKERRQ(ierr);
41134a97f8cSStefano Zampini   ierr = PetscFree(pcbddc->deluxe_ctx);CHKERRQ(ierr);
41234a97f8cSStefano Zampini   PetscFunctionReturn(0);
41334a97f8cSStefano Zampini }
41434a97f8cSStefano Zampini 
41534a97f8cSStefano Zampini #undef __FUNCT__
41634a97f8cSStefano Zampini #define __FUNCT__ "PCBDDCScalingReset_Deluxe_Solvers"
41734a97f8cSStefano Zampini static PetscErrorCode PCBDDCScalingReset_Deluxe_Solvers(PCBDDCDeluxeScaling deluxe_ctx)
41834a97f8cSStefano Zampini {
41972b8c272SStefano Zampini   PetscInt       i;
42034a97f8cSStefano Zampini   PetscErrorCode ierr;
42134a97f8cSStefano Zampini 
42234a97f8cSStefano Zampini   PetscFunctionBegin;
42334a97f8cSStefano Zampini   ierr = PetscFree(deluxe_ctx->idx_simple_B);CHKERRQ(ierr);
42434a97f8cSStefano Zampini   deluxe_ctx->n_simple = 0;
42572b8c272SStefano Zampini   for (i=0;i<deluxe_ctx->seq_n;i++) {
42672b8c272SStefano Zampini     ierr = VecScatterDestroy(&deluxe_ctx->seq_scctx[i]);CHKERRQ(ierr);
42772b8c272SStefano Zampini     ierr = VecDestroy(&deluxe_ctx->seq_work1[i]);CHKERRQ(ierr);
42872b8c272SStefano Zampini     ierr = VecDestroy(&deluxe_ctx->seq_work2[i]);CHKERRQ(ierr);
42972b8c272SStefano Zampini     ierr = MatDestroy(&deluxe_ctx->seq_mat[i]);CHKERRQ(ierr);
43072b8c272SStefano Zampini     ierr = MatDestroy(&deluxe_ctx->seq_mat_inv_sum[i]);CHKERRQ(ierr);
43172b8c272SStefano Zampini   }
43272b8c272SStefano 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);
43372b8c272SStefano Zampini   ierr = PetscFree(deluxe_ctx->workspace);CHKERRQ(ierr);
43472b8c272SStefano Zampini   deluxe_ctx->seq_n = 0;
43534a97f8cSStefano Zampini   PetscFunctionReturn(0);
43634a97f8cSStefano Zampini }
43734a97f8cSStefano Zampini 
43834a97f8cSStefano Zampini #undef __FUNCT__
43934a97f8cSStefano Zampini #define __FUNCT__ "PCBDDCScalingSetUp_Deluxe"
44034a97f8cSStefano Zampini static PetscErrorCode PCBDDCScalingSetUp_Deluxe(PC pc)
44134a97f8cSStefano Zampini {
44234a97f8cSStefano Zampini   PC_IS               *pcis=(PC_IS*)pc->data;
44334a97f8cSStefano Zampini   PC_BDDC             *pcbddc=(PC_BDDC*)pc->data;
44434a97f8cSStefano Zampini   PCBDDCDeluxeScaling deluxe_ctx=pcbddc->deluxe_ctx;
445b96c3477SStefano Zampini   PCBDDCSubSchurs     sub_schurs=pcbddc->sub_schurs;
44634a97f8cSStefano Zampini   PetscErrorCode      ierr;
44734a97f8cSStefano Zampini 
44834a97f8cSStefano Zampini   PetscFunctionBegin;
44972b8c272SStefano Zampini   /* reset data structures if the topology has changed */
45072b8c272SStefano Zampini   if (pcbddc->recompute_topography) {
45134a97f8cSStefano Zampini     ierr = PCBDDCScalingReset_Deluxe_Solvers(deluxe_ctx);CHKERRQ(ierr);
45272b8c272SStefano Zampini   }
45334a97f8cSStefano Zampini 
454b1b3d7a2SStefano Zampini   /* Compute data structures to solve sequential problems */
4555a95e1ceSStefano Zampini   ierr = PCBDDCScalingSetUp_Deluxe_Private(pc);CHKERRQ(ierr);
4565db18549SStefano Zampini 
457b1b3d7a2SStefano Zampini   /* diagonal scaling on interface dofs not contained in cc */
458d62866d3SStefano Zampini   if (sub_schurs->is_vertices || sub_schurs->is_dir) {
4591b968477SStefano Zampini     PetscInt n_com,n_dir;
4601b968477SStefano Zampini     n_com = 0;
461d62866d3SStefano Zampini     if (sub_schurs->is_vertices) {
462d62866d3SStefano Zampini       ierr = ISGetLocalSize(sub_schurs->is_vertices,&n_com);CHKERRQ(ierr);
4631b968477SStefano Zampini     }
4641b968477SStefano Zampini     n_dir = 0;
465d62866d3SStefano Zampini     if (sub_schurs->is_dir) {
466d62866d3SStefano Zampini       ierr = ISGetLocalSize(sub_schurs->is_dir,&n_dir);CHKERRQ(ierr);
4671b968477SStefano Zampini     }
46872b8c272SStefano Zampini     if (!deluxe_ctx->n_simple) {
4691b968477SStefano Zampini       deluxe_ctx->n_simple = n_dir + n_com;
4701b968477SStefano Zampini       ierr = PetscMalloc1(deluxe_ctx->n_simple,&deluxe_ctx->idx_simple_B);CHKERRQ(ierr);
471d62866d3SStefano Zampini       if (sub_schurs->is_vertices) {
4729bb4a8caSStefano Zampini         PetscInt       nmap;
4739bb4a8caSStefano Zampini         const PetscInt *idxs;
4741b968477SStefano Zampini 
475d62866d3SStefano Zampini         ierr = ISGetIndices(sub_schurs->is_vertices,&idxs);CHKERRQ(ierr);
4761b968477SStefano Zampini         ierr = ISGlobalToLocalMappingApply(pcis->BtoNmap,IS_GTOLM_DROP,n_com,idxs,&nmap,deluxe_ctx->idx_simple_B);CHKERRQ(ierr);
477af25d912SStefano 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);
478d62866d3SStefano Zampini         ierr = ISRestoreIndices(sub_schurs->is_vertices,&idxs);CHKERRQ(ierr);
4791b968477SStefano Zampini       }
480d62866d3SStefano Zampini       if (sub_schurs->is_dir) {
4811b968477SStefano Zampini         PetscInt       nmap;
4821b968477SStefano Zampini         const PetscInt *idxs;
4831b968477SStefano Zampini 
484d62866d3SStefano Zampini         ierr = ISGetIndices(sub_schurs->is_dir,&idxs);CHKERRQ(ierr);
4851b968477SStefano Zampini         ierr = ISGlobalToLocalMappingApply(pcis->BtoNmap,IS_GTOLM_DROP,n_dir,idxs,&nmap,deluxe_ctx->idx_simple_B+n_com);CHKERRQ(ierr);
486af25d912SStefano 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);
487d62866d3SStefano Zampini         ierr = ISRestoreIndices(sub_schurs->is_dir,&idxs);CHKERRQ(ierr);
4881b968477SStefano Zampini       }
4891b968477SStefano Zampini       ierr = PetscSortInt(deluxe_ctx->n_simple,deluxe_ctx->idx_simple_B);CHKERRQ(ierr);
4909bb4a8caSStefano Zampini     } else {
491af25d912SStefano 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);
49272b8c272SStefano Zampini     }
49372b8c272SStefano Zampini   } else {
494b1b3d7a2SStefano Zampini     deluxe_ctx->n_simple = 0;
4959bb4a8caSStefano Zampini     deluxe_ctx->idx_simple_B = 0;
496b1b3d7a2SStefano Zampini   }
49734a97f8cSStefano Zampini   PetscFunctionReturn(0);
49834a97f8cSStefano Zampini }
49934a97f8cSStefano Zampini 
50034a97f8cSStefano Zampini #undef __FUNCT__
5015a95e1ceSStefano Zampini #define __FUNCT__ "PCBDDCScalingSetUp_Deluxe_Private"
5025a95e1ceSStefano Zampini static PetscErrorCode PCBDDCScalingSetUp_Deluxe_Private(PC pc)
5035db18549SStefano Zampini {
5045db18549SStefano Zampini   PC_BDDC                *pcbddc=(PC_BDDC*)pc->data;
5055db18549SStefano Zampini   PCBDDCDeluxeScaling    deluxe_ctx=pcbddc->deluxe_ctx;
506b96c3477SStefano Zampini   PCBDDCSubSchurs        sub_schurs = pcbddc->sub_schurs;
50772b8c272SStefano Zampini   PetscScalar            *matdata,*matdata2;
50872b8c272SStefano Zampini   PetscInt               i,max_subset_size,cum,cum2;
50972b8c272SStefano Zampini   const PetscInt         *idxs;
51072b8c272SStefano Zampini   PetscBool              newsetup = PETSC_FALSE;
5115db18549SStefano Zampini   PetscErrorCode         ierr;
5125db18549SStefano Zampini 
5135db18549SStefano Zampini   PetscFunctionBegin;
5145a95e1ceSStefano Zampini   if (!sub_schurs->n_subs) {
5159221af80SStefano Zampini     PetscFunctionReturn(0);
5169221af80SStefano Zampini   }
5179221af80SStefano Zampini 
51872b8c272SStefano Zampini   /* Allocate arrays for subproblems */
51972b8c272SStefano Zampini   if (!deluxe_ctx->seq_n) {
52072b8c272SStefano Zampini     deluxe_ctx->seq_n = sub_schurs->n_subs;
52172b8c272SStefano 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);
52272b8c272SStefano Zampini     newsetup = PETSC_TRUE;
52372b8c272SStefano Zampini   } else if (deluxe_ctx->seq_n != sub_schurs->n_subs) {
52472b8c272SStefano 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);
52572b8c272SStefano Zampini   }
52672b8c272SStefano Zampini   deluxe_ctx->change_with_qr = sub_schurs->change_with_qr;
5275db18549SStefano Zampini 
52872b8c272SStefano Zampini   /* Create objects for deluxe */
52972b8c272SStefano Zampini   max_subset_size = 0;
53072b8c272SStefano Zampini   for (i=0;i<sub_schurs->n_subs;i++) {
53172b8c272SStefano Zampini     PetscInt subset_size;
53272b8c272SStefano Zampini     ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr);
53372b8c272SStefano Zampini     max_subset_size = PetscMax(subset_size,max_subset_size);
53472b8c272SStefano Zampini   }
53572b8c272SStefano Zampini   if (newsetup) {
53672b8c272SStefano Zampini     ierr = PetscMalloc1(2*max_subset_size,&deluxe_ctx->workspace);CHKERRQ(ierr);
53772b8c272SStefano Zampini   }
53872b8c272SStefano Zampini   cum = cum2 = 0;
53972b8c272SStefano Zampini   ierr = ISGetIndices(sub_schurs->is_Ej_all,&idxs);CHKERRQ(ierr);
54072b8c272SStefano Zampini   ierr = MatSeqAIJGetArray(sub_schurs->S_Ej_all,&matdata);CHKERRQ(ierr);
54172b8c272SStefano Zampini   ierr = MatSeqAIJGetArray(sub_schurs->sum_S_Ej_all,&matdata2);CHKERRQ(ierr);
54272b8c272SStefano Zampini   for (i=0;i<deluxe_ctx->seq_n;i++) {
54372b8c272SStefano Zampini     PetscInt     subset_size;
544925dfd53SStefano Zampini 
54572b8c272SStefano Zampini     ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr);
54672b8c272SStefano Zampini     if (newsetup) {
54772b8c272SStefano Zampini       IS  sub;
54872b8c272SStefano Zampini       /* work vectors */
54972b8c272SStefano Zampini       ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,1,subset_size,deluxe_ctx->workspace,&deluxe_ctx->seq_work1[i]);CHKERRQ(ierr);
55072b8c272SStefano Zampini       ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,1,subset_size,deluxe_ctx->workspace+subset_size,&deluxe_ctx->seq_work2[i]);CHKERRQ(ierr);
55172b8c272SStefano Zampini 
55272b8c272SStefano Zampini       /* scatters */
55372b8c272SStefano Zampini       ierr = ISCreateGeneral(PETSC_COMM_SELF,subset_size,idxs+cum,PETSC_COPY_VALUES,&sub);CHKERRQ(ierr);
55472b8c272SStefano Zampini       ierr = VecScatterCreate(pcbddc->work_scaling,sub,deluxe_ctx->seq_work1[i],NULL,&deluxe_ctx->seq_scctx[i]);CHKERRQ(ierr);
55572b8c272SStefano Zampini       ierr = ISDestroy(&sub);CHKERRQ(ierr);
556839e9adbSstefano_zampini     }
55772b8c272SStefano Zampini 
55872b8c272SStefano Zampini     /* S_E_j */
559839e9adbSstefano_zampini     ierr = MatDestroy(&deluxe_ctx->seq_mat[i]);CHKERRQ(ierr);
56072b8c272SStefano Zampini     ierr = MatCreateSeqDense(PETSC_COMM_SELF,subset_size,subset_size,matdata+cum2,&deluxe_ctx->seq_mat[i]);CHKERRQ(ierr);
561839e9adbSstefano_zampini 
56272b8c272SStefano Zampini     /* \sum_k S^k_E_j */
56372b8c272SStefano Zampini     ierr = MatDestroy(&deluxe_ctx->seq_mat_inv_sum[i]);CHKERRQ(ierr);
56472b8c272SStefano Zampini     ierr = MatCreateSeqDense(PETSC_COMM_SELF,subset_size,subset_size,matdata2+cum2,&deluxe_ctx->seq_mat_inv_sum[i]);CHKERRQ(ierr);
565*729c86d3Sstefano_zampini 
56672b8c272SStefano Zampini     if (sub_schurs->is_hermitian && sub_schurs->is_posdef) {
56772b8c272SStefano Zampini       ierr = MatCholeskyFactor(deluxe_ctx->seq_mat_inv_sum[i],NULL,NULL);CHKERRQ(ierr);
568d62866d3SStefano Zampini     } else {
56972b8c272SStefano Zampini       ierr = MatLUFactor(deluxe_ctx->seq_mat_inv_sum[i],NULL,NULL,NULL);CHKERRQ(ierr);
570d62866d3SStefano Zampini     }
571839e9adbSstefano_zampini     if (pcbddc->deluxe_singlemat) {
572839e9adbSstefano_zampini       Mat X,Y;
573839e9adbSstefano_zampini 
574839e9adbSstefano_zampini       if (!sub_schurs->is_hermitian || !sub_schurs->is_posdef) {
575839e9adbSstefano_zampini         ierr = MatTranspose(deluxe_ctx->seq_mat[i],MAT_INITIAL_MATRIX,&X);CHKERRQ(ierr);
576839e9adbSstefano_zampini       } else {
577839e9adbSstefano_zampini         ierr = PetscObjectReference((PetscObject)deluxe_ctx->seq_mat[i]);CHKERRQ(ierr);
578839e9adbSstefano_zampini         X    = deluxe_ctx->seq_mat[i];
579839e9adbSstefano_zampini       }
580839e9adbSstefano_zampini       ierr = MatDuplicate(X,MAT_DO_NOT_COPY_VALUES,&Y);CHKERRQ(ierr);
581839e9adbSstefano_zampini       if (!sub_schurs->is_hermitian || !sub_schurs->is_posdef) {
582839e9adbSstefano_zampini         ierr = PCBDDCMatTransposeMatSolve_SeqDense(deluxe_ctx->seq_mat_inv_sum[i],X,Y);CHKERRQ(ierr);
583839e9adbSstefano_zampini       } else {
584839e9adbSstefano_zampini         ierr = MatMatSolve(deluxe_ctx->seq_mat_inv_sum[i],X,Y);CHKERRQ(ierr);
585839e9adbSstefano_zampini       }
586*729c86d3Sstefano_zampini 
587839e9adbSstefano_zampini       ierr = MatDestroy(&deluxe_ctx->seq_mat_inv_sum[i]);CHKERRQ(ierr);
588839e9adbSstefano_zampini       ierr = MatDestroy(&deluxe_ctx->seq_mat[i]);CHKERRQ(ierr);
589839e9adbSstefano_zampini       ierr = MatDestroy(&X);CHKERRQ(ierr);
590*729c86d3Sstefano_zampini       ierr = MatTranspose(Y,MAT_REUSE_MATRIX,&Y);CHKERRQ(ierr);
591839e9adbSstefano_zampini       deluxe_ctx->seq_mat[i] = Y;
59272b8c272SStefano Zampini     }
59372b8c272SStefano Zampini     cum += subset_size;
59472b8c272SStefano Zampini     cum2 += subset_size*subset_size;
59572b8c272SStefano Zampini   }
59672b8c272SStefano Zampini   ierr = ISRestoreIndices(sub_schurs->is_Ej_all,&idxs);CHKERRQ(ierr);
59772b8c272SStefano Zampini   ierr = MatSeqAIJRestoreArray(sub_schurs->S_Ej_all,&matdata);CHKERRQ(ierr);
59872b8c272SStefano Zampini   ierr = MatSeqAIJRestoreArray(sub_schurs->sum_S_Ej_all,&matdata2);CHKERRQ(ierr);
5995db18549SStefano Zampini 
60072b8c272SStefano Zampini   /* the change of basis is just a reference to sub_schurs->change (if any) */
60116e386b8SStefano Zampini   deluxe_ctx->change = sub_schurs->change;
60272b8c272SStefano Zampini   if (deluxe_ctx->change && !deluxe_ctx->change_with_qr) {
60372b8c272SStefano Zampini     for (i=0;i<deluxe_ctx->seq_n;i++) {
60472b8c272SStefano Zampini       if (newsetup) {
60572b8c272SStefano Zampini         PC pc;
60672b8c272SStefano Zampini 
60772b8c272SStefano Zampini         ierr = KSPGetPC(deluxe_ctx->change[i],&pc);CHKERRQ(ierr);
60872b8c272SStefano Zampini         ierr = PCSetType(pc,PCLU);CHKERRQ(ierr);
60972b8c272SStefano Zampini         ierr = KSPSetFromOptions(deluxe_ctx->change[i]);CHKERRQ(ierr);
61072b8c272SStefano Zampini       }
61172b8c272SStefano Zampini       ierr = KSPSetUp(deluxe_ctx->change[i]);CHKERRQ(ierr);
61272b8c272SStefano Zampini     }
61316e386b8SStefano Zampini   }
6145db18549SStefano Zampini   PetscFunctionReturn(0);
6155db18549SStefano Zampini }
616