xref: /petsc/src/ksp/pc/impls/bddc/bddcscalingbasic.c (revision 3829e7e6d9b23a878abd312e0c1d82bde59983dd)
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
43*3829e7e6SStefano Zampini   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Only LU factor supported");
44839e9adbSstefano_zampini 
45839e9adbSstefano_zampini   ierr = MatDenseRestoreArray(B,&b);CHKERRQ(ierr);
46839e9adbSstefano_zampini   ierr = MatDenseRestoreArray(X,&x);CHKERRQ(ierr);
47839e9adbSstefano_zampini   ierr = PetscLogFlops(nrhs*(2.0*m*m - m));CHKERRQ(ierr);
48839e9adbSstefano_zampini   PetscFunctionReturn(0);
49839e9adbSstefano_zampini }
50839e9adbSstefano_zampini 
51674ae819SStefano Zampini static PetscErrorCode PCBDDCScalingExtension_Basic(PC pc, Vec local_interface_vector, Vec global_vector)
52674ae819SStefano Zampini {
53674ae819SStefano Zampini   PC_IS*         pcis = (PC_IS*)pc->data;
54674ae819SStefano Zampini   PC_BDDC*       pcbddc = (PC_BDDC*)pc->data;
55674ae819SStefano Zampini   PetscErrorCode ierr;
56674ae819SStefano Zampini 
57674ae819SStefano Zampini   PetscFunctionBegin;
58674ae819SStefano Zampini   /* Apply partition of unity */
59674ae819SStefano Zampini   ierr = VecPointwiseMult(pcbddc->work_scaling,pcis->D,local_interface_vector);CHKERRQ(ierr);
60674ae819SStefano Zampini   ierr = VecSet(global_vector,0.0);CHKERRQ(ierr);
61674ae819SStefano Zampini   ierr = VecScatterBegin(pcis->global_to_B,pcbddc->work_scaling,global_vector,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
62674ae819SStefano Zampini   ierr = VecScatterEnd(pcis->global_to_B,pcbddc->work_scaling,global_vector,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
63674ae819SStefano Zampini   PetscFunctionReturn(0);
64674ae819SStefano Zampini }
65674ae819SStefano Zampini 
66674ae819SStefano Zampini static PetscErrorCode PCBDDCScalingExtension_Deluxe(PC pc, Vec x, Vec y)
67674ae819SStefano Zampini {
68674ae819SStefano Zampini   PC_IS*              pcis=(PC_IS*)pc->data;
69674ae819SStefano Zampini   PC_BDDC*            pcbddc=(PC_BDDC*)pc->data;
70674ae819SStefano Zampini   PCBDDCDeluxeScaling deluxe_ctx = pcbddc->deluxe_ctx;
71674ae819SStefano Zampini   PetscErrorCode      ierr;
72674ae819SStefano Zampini 
73674ae819SStefano Zampini   PetscFunctionBegin;
745a95e1ceSStefano Zampini   ierr = VecSet(pcbddc->work_scaling,0.0);CHKERRQ(ierr);
755a95e1ceSStefano Zampini   ierr = VecSet(y,0.0);CHKERRQ(ierr);
765a95e1ceSStefano Zampini   if (deluxe_ctx->n_simple) { /* scale deluxe vertices using diagonal scaling */
775a95e1ceSStefano Zampini     PetscInt          i;
782b095fd8SStefano Zampini     const PetscScalar *array_x,*array_D;
792b095fd8SStefano Zampini     PetscScalar       *array;
802b095fd8SStefano Zampini     ierr = VecGetArrayRead(x,&array_x);CHKERRQ(ierr);
812b095fd8SStefano Zampini     ierr = VecGetArrayRead(pcis->D,&array_D);CHKERRQ(ierr);
82674ae819SStefano Zampini     ierr = VecGetArray(pcbddc->work_scaling,&array);CHKERRQ(ierr);
83674ae819SStefano Zampini     for (i=0;i<deluxe_ctx->n_simple;i++) {
84674ae819SStefano Zampini       array[deluxe_ctx->idx_simple_B[i]] = array_x[deluxe_ctx->idx_simple_B[i]]*array_D[deluxe_ctx->idx_simple_B[i]];
85674ae819SStefano Zampini     }
86674ae819SStefano Zampini     ierr = VecRestoreArray(pcbddc->work_scaling,&array);CHKERRQ(ierr);
872b095fd8SStefano Zampini     ierr = VecRestoreArrayRead(pcis->D,&array_D);CHKERRQ(ierr);
882b095fd8SStefano Zampini     ierr = VecRestoreArrayRead(x,&array_x);CHKERRQ(ierr);
8934a97f8cSStefano Zampini   }
90ac632422SStefano Zampini   /* sequential part : all problems and Schur applications collapsed into a single matrix vector multiplication or a matvec and a solve */
915a95e1ceSStefano Zampini   if (deluxe_ctx->seq_mat) {
9272b8c272SStefano Zampini     PetscInt i;
9372b8c272SStefano Zampini     for (i=0;i<deluxe_ctx->seq_n;i++) {
9416e386b8SStefano Zampini       if (deluxe_ctx->change) {
9572b8c272SStefano Zampini         ierr = VecScatterBegin(deluxe_ctx->seq_scctx[i],x,deluxe_ctx->seq_work2[i],INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
9672b8c272SStefano Zampini         ierr = VecScatterEnd(deluxe_ctx->seq_scctx[i],x,deluxe_ctx->seq_work2[i],INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
9772b8c272SStefano Zampini         if (deluxe_ctx->change_with_qr) {
9872b8c272SStefano Zampini           Mat change;
9972b8c272SStefano Zampini 
10072b8c272SStefano Zampini           ierr = KSPGetOperators(deluxe_ctx->change[i],&change,NULL);CHKERRQ(ierr);
10172b8c272SStefano Zampini           ierr = MatMultTranspose(change,deluxe_ctx->seq_work2[i],deluxe_ctx->seq_work1[i]);CHKERRQ(ierr);
10272b8c272SStefano Zampini         } else {
10372b8c272SStefano Zampini           ierr = KSPSolve(deluxe_ctx->change[i],deluxe_ctx->seq_work2[i],deluxe_ctx->seq_work1[i]);CHKERRQ(ierr);
10416e386b8SStefano Zampini         }
10572b8c272SStefano Zampini       } else {
10672b8c272SStefano Zampini         ierr = VecScatterBegin(deluxe_ctx->seq_scctx[i],x,deluxe_ctx->seq_work1[i],INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
10772b8c272SStefano Zampini         ierr = VecScatterEnd(deluxe_ctx->seq_scctx[i],x,deluxe_ctx->seq_work1[i],INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
10872b8c272SStefano Zampini       }
10972b8c272SStefano Zampini       ierr = MatMultTranspose(deluxe_ctx->seq_mat[i],deluxe_ctx->seq_work1[i],deluxe_ctx->seq_work2[i]);CHKERRQ(ierr);
11072b8c272SStefano Zampini       if (deluxe_ctx->seq_mat_inv_sum[i]) {
11172b8c272SStefano Zampini         PetscScalar *x;
11272b8c272SStefano Zampini 
11372b8c272SStefano Zampini         ierr = VecGetArray(deluxe_ctx->seq_work2[i],&x);CHKERRQ(ierr);
11472b8c272SStefano Zampini         ierr = VecPlaceArray(deluxe_ctx->seq_work1[i],x);CHKERRQ(ierr);
11572b8c272SStefano Zampini         ierr = VecRestoreArray(deluxe_ctx->seq_work2[i],&x);CHKERRQ(ierr);
11672b8c272SStefano Zampini         ierr = MatSolveTranspose(deluxe_ctx->seq_mat_inv_sum[i],deluxe_ctx->seq_work1[i],deluxe_ctx->seq_work2[i]);CHKERRQ(ierr);
11772b8c272SStefano Zampini         ierr = VecResetArray(deluxe_ctx->seq_work1[i]);CHKERRQ(ierr);
118ac632422SStefano Zampini       }
11916e386b8SStefano Zampini       if (deluxe_ctx->change) {
12016e386b8SStefano Zampini         Mat change;
12116e386b8SStefano Zampini 
12272b8c272SStefano Zampini         ierr = KSPGetOperators(deluxe_ctx->change[i],&change,NULL);CHKERRQ(ierr);
12372b8c272SStefano Zampini         ierr = MatMult(change,deluxe_ctx->seq_work2[i],deluxe_ctx->seq_work1[i]);CHKERRQ(ierr);
12472b8c272SStefano Zampini         ierr = VecScatterBegin(deluxe_ctx->seq_scctx[i],deluxe_ctx->seq_work1[i],pcbddc->work_scaling,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
12572b8c272SStefano Zampini         ierr = VecScatterEnd(deluxe_ctx->seq_scctx[i],deluxe_ctx->seq_work1[i],pcbddc->work_scaling,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
12616e386b8SStefano Zampini       } else {
12772b8c272SStefano Zampini         ierr = VecScatterBegin(deluxe_ctx->seq_scctx[i],deluxe_ctx->seq_work2[i],pcbddc->work_scaling,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
12872b8c272SStefano Zampini         ierr = VecScatterEnd(deluxe_ctx->seq_scctx[i],deluxe_ctx->seq_work2[i],pcbddc->work_scaling,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
12972b8c272SStefano Zampini       }
130674ae819SStefano Zampini     }
13116e386b8SStefano Zampini   }
132674ae819SStefano Zampini   /* put local boundary part in global vector */
133674ae819SStefano Zampini   ierr = VecScatterBegin(pcis->global_to_B,pcbddc->work_scaling,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
134674ae819SStefano Zampini   ierr = VecScatterEnd(pcis->global_to_B,pcbddc->work_scaling,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
135674ae819SStefano Zampini   PetscFunctionReturn(0);
136674ae819SStefano Zampini }
137674ae819SStefano Zampini 
138674ae819SStefano Zampini PetscErrorCode PCBDDCScalingExtension(PC pc, Vec local_interface_vector, Vec global_vector)
139674ae819SStefano Zampini {
140674ae819SStefano Zampini   PC_BDDC        *pcbddc=(PC_BDDC*)pc->data;
141674ae819SStefano Zampini   PetscErrorCode ierr;
142674ae819SStefano Zampini 
143674ae819SStefano Zampini   PetscFunctionBegin;
144674ae819SStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
145674ae819SStefano Zampini   PetscValidHeaderSpecific(local_interface_vector,VEC_CLASSID,2);
146674ae819SStefano Zampini   PetscValidHeaderSpecific(global_vector,VEC_CLASSID,3);
1476c4ed002SBarry Smith   if (local_interface_vector == pcbddc->work_scaling) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Local vector cannot be pcbddc->work_scaling!\n");
148163d334eSBarry Smith   ierr = PetscUseMethod(pc,"PCBDDCScalingExtension_C",(PC,Vec,Vec),(pc,local_interface_vector,global_vector));CHKERRQ(ierr);
149674ae819SStefano Zampini   PetscFunctionReturn(0);
150674ae819SStefano Zampini }
151674ae819SStefano Zampini 
152674ae819SStefano Zampini static PetscErrorCode PCBDDCScalingRestriction_Basic(PC pc, Vec global_vector, Vec local_interface_vector)
153674ae819SStefano Zampini {
154674ae819SStefano Zampini   PetscErrorCode ierr;
155674ae819SStefano Zampini   PC_IS          *pcis = (PC_IS*)pc->data;
156674ae819SStefano Zampini 
157674ae819SStefano Zampini   PetscFunctionBegin;
158674ae819SStefano Zampini   ierr = VecScatterBegin(pcis->global_to_B,global_vector,local_interface_vector,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
159674ae819SStefano Zampini   ierr = VecScatterEnd(pcis->global_to_B,global_vector,local_interface_vector,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
160674ae819SStefano Zampini   /* Apply partition of unity */
161674ae819SStefano Zampini   ierr = VecPointwiseMult(local_interface_vector,pcis->D,local_interface_vector);CHKERRQ(ierr);
162674ae819SStefano Zampini   PetscFunctionReturn(0);
163674ae819SStefano Zampini }
164674ae819SStefano Zampini 
165674ae819SStefano Zampini static PetscErrorCode PCBDDCScalingRestriction_Deluxe(PC pc, Vec x, Vec y)
166674ae819SStefano Zampini {
167674ae819SStefano Zampini   PC_IS*              pcis=(PC_IS*)pc->data;
168674ae819SStefano Zampini   PC_BDDC*            pcbddc=(PC_BDDC*)pc->data;
169674ae819SStefano Zampini   PCBDDCDeluxeScaling deluxe_ctx = pcbddc->deluxe_ctx;
170674ae819SStefano Zampini   PetscErrorCode      ierr;
171674ae819SStefano Zampini 
172674ae819SStefano Zampini   PetscFunctionBegin;
173674ae819SStefano Zampini   /* get local boundary part of global vector */
174674ae819SStefano Zampini   ierr = VecScatterBegin(pcis->global_to_B,x,y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
175674ae819SStefano Zampini   ierr = VecScatterEnd(pcis->global_to_B,x,y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1765a95e1ceSStefano Zampini   if (deluxe_ctx->n_simple) { /* scale deluxe vertices using diagonal scaling */
1775a95e1ceSStefano Zampini     PetscInt          i;
1782b095fd8SStefano Zampini     PetscScalar       *array_y;
1792b095fd8SStefano Zampini     const PetscScalar *array_D;
180674ae819SStefano Zampini     ierr = VecGetArray(y,&array_y);CHKERRQ(ierr);
1812b095fd8SStefano Zampini     ierr = VecGetArrayRead(pcis->D,&array_D);CHKERRQ(ierr);
182674ae819SStefano Zampini     for (i=0;i<deluxe_ctx->n_simple;i++) {
183674ae819SStefano Zampini       array_y[deluxe_ctx->idx_simple_B[i]] *= array_D[deluxe_ctx->idx_simple_B[i]];
184674ae819SStefano Zampini     }
1852b095fd8SStefano Zampini     ierr = VecRestoreArrayRead(pcis->D,&array_D);CHKERRQ(ierr);
186674ae819SStefano Zampini     ierr = VecRestoreArray(y,&array_y);CHKERRQ(ierr);
18734a97f8cSStefano Zampini   }
188729c86d3Sstefano_zampini   /* sequential part : all problems and Schur applications collapsed into a single matrix vector multiplication or a matvec and a solve */
1895a95e1ceSStefano Zampini   if (deluxe_ctx->seq_mat) {
19072b8c272SStefano Zampini     PetscInt i;
19172b8c272SStefano Zampini     for (i=0;i<deluxe_ctx->seq_n;i++) {
19216e386b8SStefano Zampini       if (deluxe_ctx->change) {
19316e386b8SStefano Zampini         Mat change;
19416e386b8SStefano Zampini 
19572b8c272SStefano Zampini         ierr = VecScatterBegin(deluxe_ctx->seq_scctx[i],y,deluxe_ctx->seq_work2[i],INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
19672b8c272SStefano Zampini         ierr = VecScatterEnd(deluxe_ctx->seq_scctx[i],y,deluxe_ctx->seq_work2[i],INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
19772b8c272SStefano Zampini         ierr = KSPGetOperators(deluxe_ctx->change[i],&change,NULL);CHKERRQ(ierr);
19872b8c272SStefano Zampini         ierr = MatMultTranspose(change,deluxe_ctx->seq_work2[i],deluxe_ctx->seq_work1[i]);CHKERRQ(ierr);
19916e386b8SStefano Zampini       } else {
20072b8c272SStefano Zampini         ierr = VecScatterBegin(deluxe_ctx->seq_scctx[i],y,deluxe_ctx->seq_work1[i],INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
20172b8c272SStefano Zampini         ierr = VecScatterEnd(deluxe_ctx->seq_scctx[i],y,deluxe_ctx->seq_work1[i],INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
20216e386b8SStefano Zampini       }
20372b8c272SStefano Zampini       if (deluxe_ctx->seq_mat_inv_sum[i]) {
20472b8c272SStefano Zampini         PetscScalar *x;
20572b8c272SStefano Zampini 
20672b8c272SStefano Zampini         ierr = VecGetArray(deluxe_ctx->seq_work1[i],&x);CHKERRQ(ierr);
20772b8c272SStefano Zampini         ierr = VecPlaceArray(deluxe_ctx->seq_work2[i],x);CHKERRQ(ierr);
20872b8c272SStefano Zampini         ierr = VecRestoreArray(deluxe_ctx->seq_work1[i],&x);CHKERRQ(ierr);
20972b8c272SStefano Zampini         ierr = MatSolve(deluxe_ctx->seq_mat_inv_sum[i],deluxe_ctx->seq_work1[i],deluxe_ctx->seq_work2[i]);CHKERRQ(ierr);
21072b8c272SStefano Zampini         ierr = VecResetArray(deluxe_ctx->seq_work2[i]);CHKERRQ(ierr);
211ac632422SStefano Zampini       }
21272b8c272SStefano Zampini       ierr = MatMult(deluxe_ctx->seq_mat[i],deluxe_ctx->seq_work1[i],deluxe_ctx->seq_work2[i]);CHKERRQ(ierr);
21316e386b8SStefano Zampini       if (deluxe_ctx->change) {
21472b8c272SStefano Zampini         if (deluxe_ctx->change_with_qr) {
21572b8c272SStefano Zampini           Mat change;
21672b8c272SStefano Zampini 
21772b8c272SStefano Zampini           ierr = KSPGetOperators(deluxe_ctx->change[i],&change,NULL);CHKERRQ(ierr);
21872b8c272SStefano Zampini           ierr = MatMult(change,deluxe_ctx->seq_work2[i],deluxe_ctx->seq_work1[i]);CHKERRQ(ierr);
21972b8c272SStefano Zampini         } else {
22072b8c272SStefano Zampini           ierr = KSPSolveTranspose(deluxe_ctx->change[i],deluxe_ctx->seq_work2[i],deluxe_ctx->seq_work1[i]);CHKERRQ(ierr);
22116e386b8SStefano Zampini         }
22272b8c272SStefano Zampini         ierr = VecScatterBegin(deluxe_ctx->seq_scctx[i],deluxe_ctx->seq_work1[i],y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
22372b8c272SStefano Zampini         ierr = VecScatterEnd(deluxe_ctx->seq_scctx[i],deluxe_ctx->seq_work1[i],y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
22472b8c272SStefano Zampini       } else {
22572b8c272SStefano Zampini         ierr = VecScatterBegin(deluxe_ctx->seq_scctx[i],deluxe_ctx->seq_work2[i],y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
22672b8c272SStefano Zampini         ierr = VecScatterEnd(deluxe_ctx->seq_scctx[i],deluxe_ctx->seq_work2[i],y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
22772b8c272SStefano Zampini       }
22872b8c272SStefano Zampini     }
229674ae819SStefano Zampini   }
230674ae819SStefano Zampini   PetscFunctionReturn(0);
231674ae819SStefano Zampini }
232674ae819SStefano Zampini 
233674ae819SStefano Zampini PetscErrorCode PCBDDCScalingRestriction(PC pc, Vec global_vector, Vec local_interface_vector)
234674ae819SStefano Zampini {
235674ae819SStefano Zampini   PC_BDDC        *pcbddc=(PC_BDDC*)pc->data;
236674ae819SStefano Zampini   PetscErrorCode ierr;
237674ae819SStefano Zampini 
238674ae819SStefano Zampini   PetscFunctionBegin;
239674ae819SStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
240674ae819SStefano Zampini   PetscValidHeaderSpecific(global_vector,VEC_CLASSID,2);
241674ae819SStefano Zampini   PetscValidHeaderSpecific(local_interface_vector,VEC_CLASSID,3);
2426c4ed002SBarry Smith   if (local_interface_vector == pcbddc->work_scaling) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Local vector cannot be pcbddc->work_scaling!\n");
243163d334eSBarry Smith   ierr = PetscUseMethod(pc,"PCBDDCScalingRestriction_C",(PC,Vec,Vec),(pc,global_vector,local_interface_vector));CHKERRQ(ierr);
244674ae819SStefano Zampini   PetscFunctionReturn(0);
245674ae819SStefano Zampini }
246674ae819SStefano Zampini 
247674ae819SStefano Zampini PetscErrorCode PCBDDCScalingSetUp(PC pc)
248674ae819SStefano Zampini {
249674ae819SStefano Zampini   PC_IS*         pcis=(PC_IS*)pc->data;
250674ae819SStefano Zampini   PC_BDDC*       pcbddc=(PC_BDDC*)pc->data;
251674ae819SStefano Zampini   PetscErrorCode ierr;
252674ae819SStefano Zampini 
253674ae819SStefano Zampini   PetscFunctionBegin;
254674ae819SStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
255674ae819SStefano Zampini   /* create work vector for the operator */
25634a97f8cSStefano Zampini   ierr = VecDestroy(&pcbddc->work_scaling);CHKERRQ(ierr);
257674ae819SStefano Zampini   ierr = VecDuplicate(pcis->vec1_B,&pcbddc->work_scaling);CHKERRQ(ierr);
25834a97f8cSStefano Zampini   /* always rebuild pcis->D */
25928d874f6SStefano Zampini   if (pcis->use_stiffness_scaling) {
260674ae819SStefano Zampini     ierr = MatGetDiagonal(pcbddc->local_mat,pcis->vec1_N);CHKERRQ(ierr);
261674ae819SStefano Zampini     ierr = VecScatterBegin(pcis->N_to_B,pcis->vec1_N,pcis->D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
262674ae819SStefano Zampini     ierr = VecScatterEnd(pcis->N_to_B,pcis->vec1_N,pcis->D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
263674ae819SStefano Zampini   }
264674ae819SStefano Zampini   ierr = VecCopy(pcis->D,pcis->vec1_B);CHKERRQ(ierr);
265674ae819SStefano Zampini   ierr = VecSet(pcis->vec1_global,0.0);CHKERRQ(ierr);
266674ae819SStefano Zampini   ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_B,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
267674ae819SStefano Zampini   ierr = VecScatterEnd(pcis->global_to_B,pcis->vec1_B,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
268674ae819SStefano Zampini   ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_global,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
269674ae819SStefano Zampini   ierr = VecScatterEnd(pcis->global_to_B,pcis->vec1_global,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
270674ae819SStefano Zampini   ierr = VecPointwiseDivide(pcis->D,pcis->D,pcis->vec1_B);CHKERRQ(ierr);
271674ae819SStefano Zampini   /* now setup */
272681e7c04SStefano Zampini   if (pcbddc->use_deluxe_scaling) {
27334a97f8cSStefano Zampini     if (!pcbddc->deluxe_ctx) {
27434a97f8cSStefano Zampini       ierr = PCBDDCScalingCreate_Deluxe(pc);CHKERRQ(ierr);
27534a97f8cSStefano Zampini     }
27634a97f8cSStefano Zampini     ierr = PCBDDCScalingSetUp_Deluxe(pc);CHKERRQ(ierr);
277674ae819SStefano Zampini     ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCScalingRestriction_C",PCBDDCScalingRestriction_Deluxe);CHKERRQ(ierr);
278674ae819SStefano Zampini     ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCScalingExtension_C",PCBDDCScalingExtension_Deluxe);CHKERRQ(ierr);
279674ae819SStefano Zampini   } else {
280674ae819SStefano Zampini     ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCScalingRestriction_C",PCBDDCScalingRestriction_Basic);CHKERRQ(ierr);
281674ae819SStefano Zampini     ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCScalingExtension_C",PCBDDCScalingExtension_Basic);CHKERRQ(ierr);
282674ae819SStefano Zampini   }
28334a97f8cSStefano Zampini 
284674ae819SStefano Zampini   /* test */
285674ae819SStefano Zampini   if (pcbddc->dbg_flag) {
28672b8c272SStefano Zampini     Mat         B0_B = NULL;
28772b8c272SStefano Zampini     Vec         B0_Bv = NULL, B0_Bv2 = NULL;
28834a97f8cSStefano Zampini     Vec         vec2_global;
289674ae819SStefano Zampini     PetscViewer viewer = pcbddc->dbg_viewer;
29034a97f8cSStefano Zampini     PetscReal   error;
291674ae819SStefano Zampini 
292674ae819SStefano Zampini     /* extension -> from local to parallel */
29334a97f8cSStefano Zampini     ierr = VecSet(pcis->vec1_global,0.0);CHKERRQ(ierr);
29434a97f8cSStefano Zampini     ierr = VecSetRandom(pcis->vec1_B,NULL);CHKERRQ(ierr);
29534a97f8cSStefano Zampini     ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_B,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
29634a97f8cSStefano Zampini     ierr = VecScatterEnd(pcis->global_to_B,pcis->vec1_B,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
29734a97f8cSStefano Zampini     ierr = VecDuplicate(pcis->vec1_global,&vec2_global);CHKERRQ(ierr);
29834a97f8cSStefano Zampini     ierr = VecCopy(pcis->vec1_global,vec2_global);CHKERRQ(ierr);
299674ae819SStefano Zampini     ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_global,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
300674ae819SStefano Zampini     ierr = VecScatterEnd(pcis->global_to_B,pcis->vec1_global,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
30172b8c272SStefano Zampini     if (pcbddc->benign_n) {
30272b8c272SStefano Zampini       IS is_dummy;
30372b8c272SStefano Zampini 
30472b8c272SStefano Zampini       ierr = ISCreateStride(PETSC_COMM_SELF,pcbddc->benign_n,0,1,&is_dummy);CHKERRQ(ierr);
3057dae84e0SHong Zhang       ierr = MatCreateSubMatrix(pcbddc->benign_B0,is_dummy,pcis->is_B_local,MAT_INITIAL_MATRIX,&B0_B);CHKERRQ(ierr);
30672b8c272SStefano Zampini       ierr = ISDestroy(&is_dummy);CHKERRQ(ierr);
30772b8c272SStefano Zampini       ierr = MatCreateVecs(B0_B,NULL,&B0_Bv);CHKERRQ(ierr);
30872b8c272SStefano Zampini       ierr = VecDuplicate(B0_Bv,&B0_Bv2);CHKERRQ(ierr);
30972b8c272SStefano Zampini       ierr = MatMult(B0_B,pcis->vec1_B,B0_Bv);CHKERRQ(ierr);
31072b8c272SStefano Zampini     }
311674ae819SStefano Zampini     ierr = PCBDDCScalingExtension(pc,pcis->vec1_B,pcis->vec1_global);CHKERRQ(ierr);
31272b8c272SStefano Zampini     if (pcbddc->benign_saddle_point) {
31372b8c272SStefano Zampini       PetscReal errorl = 0.;
31472b8c272SStefano Zampini       ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_global,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
31572b8c272SStefano Zampini       ierr = VecScatterEnd(pcis->global_to_B,pcis->vec1_global,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
31672b8c272SStefano Zampini       if (pcbddc->benign_n) {
31772b8c272SStefano Zampini         ierr = MatMult(B0_B,pcis->vec1_B,B0_Bv2);CHKERRQ(ierr);
31872b8c272SStefano Zampini         ierr = VecAXPY(B0_Bv,-1.0,B0_Bv2);CHKERRQ(ierr);
31972b8c272SStefano Zampini         ierr = VecNorm(B0_Bv,NORM_INFINITY,&errorl);CHKERRQ(ierr);
32072b8c272SStefano Zampini       }
32172b8c272SStefano Zampini       ierr = MPI_Allreduce(&errorl,&error,1,MPIU_REAL,MPI_SUM,PetscObjectComm((PetscObject)pc));CHKERRQ(ierr);
32272b8c272SStefano Zampini       ierr = PetscViewerASCIIPrintf(viewer,"Error benign extension %1.14e\n",error);CHKERRQ(ierr);
32372b8c272SStefano Zampini     }
32434a97f8cSStefano Zampini     ierr = VecAXPY(pcis->vec1_global,-1.0,vec2_global);CHKERRQ(ierr);
32534a97f8cSStefano Zampini     ierr = VecNorm(pcis->vec1_global,NORM_INFINITY,&error);CHKERRQ(ierr);
326674ae819SStefano Zampini     ierr = PetscViewerASCIIPrintf(viewer,"Error scaling extension %1.14e\n",error);CHKERRQ(ierr);
32734a97f8cSStefano Zampini     ierr = VecDestroy(&vec2_global);CHKERRQ(ierr);
32834a97f8cSStefano Zampini 
329674ae819SStefano Zampini     /* restriction -> from parallel to local */
330674ae819SStefano Zampini     ierr = VecSet(pcis->vec1_global,0.0);CHKERRQ(ierr);
33134a97f8cSStefano Zampini     ierr = VecSetRandom(pcis->vec1_B,NULL);CHKERRQ(ierr);
332674ae819SStefano Zampini     ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_B,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
333674ae819SStefano Zampini     ierr = VecScatterEnd(pcis->global_to_B,pcis->vec1_B,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
33434a97f8cSStefano Zampini     ierr = PCBDDCScalingRestriction(pc,pcis->vec1_global,pcis->vec1_B);CHKERRQ(ierr);
33534a97f8cSStefano Zampini     ierr = VecScale(pcis->vec1_B,-1.0);CHKERRQ(ierr);
33634a97f8cSStefano Zampini     ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_B,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
33734a97f8cSStefano Zampini     ierr = VecScatterEnd(pcis->global_to_B,pcis->vec1_B,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
33834a97f8cSStefano Zampini     ierr = VecNorm(pcis->vec1_global,NORM_INFINITY,&error);CHKERRQ(ierr);
33934a97f8cSStefano Zampini     ierr = PetscViewerASCIIPrintf(viewer,"Error scaling restriction %1.14e\n",error);CHKERRQ(ierr);
34072b8c272SStefano Zampini     ierr = MatDestroy(&B0_B);CHKERRQ(ierr);
34172b8c272SStefano Zampini     ierr = VecDestroy(&B0_Bv);CHKERRQ(ierr);
34272b8c272SStefano Zampini     ierr = VecDestroy(&B0_Bv2);CHKERRQ(ierr);
343674ae819SStefano Zampini   }
344674ae819SStefano Zampini   PetscFunctionReturn(0);
345674ae819SStefano Zampini }
346674ae819SStefano Zampini 
347674ae819SStefano Zampini PetscErrorCode PCBDDCScalingDestroy(PC pc)
348674ae819SStefano Zampini {
349674ae819SStefano Zampini   PC_BDDC*       pcbddc=(PC_BDDC*)pc->data;
350674ae819SStefano Zampini   PetscErrorCode ierr;
351674ae819SStefano Zampini 
352674ae819SStefano Zampini   PetscFunctionBegin;
35334a97f8cSStefano Zampini   if (pcbddc->deluxe_ctx) {
35434a97f8cSStefano Zampini     ierr = PCBDDCScalingDestroy_Deluxe(pc);CHKERRQ(ierr);
355674ae819SStefano Zampini   }
356674ae819SStefano Zampini   ierr = VecDestroy(&pcbddc->work_scaling);CHKERRQ(ierr);
357674ae819SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCScalingRestriction_C",NULL);CHKERRQ(ierr);
358674ae819SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCScalingExtension_C",NULL);CHKERRQ(ierr);
359674ae819SStefano Zampini   PetscFunctionReturn(0);
360674ae819SStefano Zampini }
361674ae819SStefano Zampini 
36234a97f8cSStefano Zampini static PetscErrorCode PCBDDCScalingCreate_Deluxe(PC pc)
36334a97f8cSStefano Zampini {
36434a97f8cSStefano Zampini   PC_BDDC*            pcbddc=(PC_BDDC*)pc->data;
36534a97f8cSStefano Zampini   PCBDDCDeluxeScaling deluxe_ctx;
36634a97f8cSStefano Zampini   PetscErrorCode      ierr;
36734a97f8cSStefano Zampini 
36834a97f8cSStefano Zampini   PetscFunctionBegin;
36934a97f8cSStefano Zampini   ierr = PetscNew(&deluxe_ctx);CHKERRQ(ierr);
37034a97f8cSStefano Zampini   pcbddc->deluxe_ctx = deluxe_ctx;
37134a97f8cSStefano Zampini   PetscFunctionReturn(0);
37234a97f8cSStefano Zampini }
37334a97f8cSStefano Zampini 
37434a97f8cSStefano Zampini static PetscErrorCode PCBDDCScalingDestroy_Deluxe(PC pc)
37534a97f8cSStefano Zampini {
37634a97f8cSStefano Zampini   PC_BDDC*            pcbddc=(PC_BDDC*)pc->data;
37734a97f8cSStefano Zampini   PetscErrorCode      ierr;
37834a97f8cSStefano Zampini 
37934a97f8cSStefano Zampini   PetscFunctionBegin;
38034a97f8cSStefano Zampini   ierr = PCBDDCScalingReset_Deluxe_Solvers(pcbddc->deluxe_ctx);CHKERRQ(ierr);
38134a97f8cSStefano Zampini   ierr = PetscFree(pcbddc->deluxe_ctx);CHKERRQ(ierr);
38234a97f8cSStefano Zampini   PetscFunctionReturn(0);
38334a97f8cSStefano Zampini }
38434a97f8cSStefano Zampini 
38534a97f8cSStefano Zampini static PetscErrorCode PCBDDCScalingReset_Deluxe_Solvers(PCBDDCDeluxeScaling deluxe_ctx)
38634a97f8cSStefano Zampini {
38772b8c272SStefano Zampini   PetscInt       i;
38834a97f8cSStefano Zampini   PetscErrorCode ierr;
38934a97f8cSStefano Zampini 
39034a97f8cSStefano Zampini   PetscFunctionBegin;
39134a97f8cSStefano Zampini   ierr = PetscFree(deluxe_ctx->idx_simple_B);CHKERRQ(ierr);
39234a97f8cSStefano Zampini   deluxe_ctx->n_simple = 0;
39372b8c272SStefano Zampini   for (i=0;i<deluxe_ctx->seq_n;i++) {
39472b8c272SStefano Zampini     ierr = VecScatterDestroy(&deluxe_ctx->seq_scctx[i]);CHKERRQ(ierr);
39572b8c272SStefano Zampini     ierr = VecDestroy(&deluxe_ctx->seq_work1[i]);CHKERRQ(ierr);
39672b8c272SStefano Zampini     ierr = VecDestroy(&deluxe_ctx->seq_work2[i]);CHKERRQ(ierr);
39772b8c272SStefano Zampini     ierr = MatDestroy(&deluxe_ctx->seq_mat[i]);CHKERRQ(ierr);
39872b8c272SStefano Zampini     ierr = MatDestroy(&deluxe_ctx->seq_mat_inv_sum[i]);CHKERRQ(ierr);
39972b8c272SStefano Zampini   }
40072b8c272SStefano 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);
40172b8c272SStefano Zampini   ierr = PetscFree(deluxe_ctx->workspace);CHKERRQ(ierr);
40272b8c272SStefano Zampini   deluxe_ctx->seq_n = 0;
40334a97f8cSStefano Zampini   PetscFunctionReturn(0);
40434a97f8cSStefano Zampini }
40534a97f8cSStefano Zampini 
40634a97f8cSStefano Zampini static PetscErrorCode PCBDDCScalingSetUp_Deluxe(PC pc)
40734a97f8cSStefano Zampini {
40834a97f8cSStefano Zampini   PC_IS               *pcis=(PC_IS*)pc->data;
40934a97f8cSStefano Zampini   PC_BDDC             *pcbddc=(PC_BDDC*)pc->data;
41034a97f8cSStefano Zampini   PCBDDCDeluxeScaling deluxe_ctx=pcbddc->deluxe_ctx;
411b96c3477SStefano Zampini   PCBDDCSubSchurs     sub_schurs=pcbddc->sub_schurs;
41234a97f8cSStefano Zampini   PetscErrorCode      ierr;
41334a97f8cSStefano Zampini 
41434a97f8cSStefano Zampini   PetscFunctionBegin;
41572b8c272SStefano Zampini   /* reset data structures if the topology has changed */
41672b8c272SStefano Zampini   if (pcbddc->recompute_topography) {
41734a97f8cSStefano Zampini     ierr = PCBDDCScalingReset_Deluxe_Solvers(deluxe_ctx);CHKERRQ(ierr);
41872b8c272SStefano Zampini   }
41934a97f8cSStefano Zampini 
420b1b3d7a2SStefano Zampini   /* Compute data structures to solve sequential problems */
4215a95e1ceSStefano Zampini   ierr = PCBDDCScalingSetUp_Deluxe_Private(pc);CHKERRQ(ierr);
4225db18549SStefano Zampini 
423b1b3d7a2SStefano Zampini   /* diagonal scaling on interface dofs not contained in cc */
424d62866d3SStefano Zampini   if (sub_schurs->is_vertices || sub_schurs->is_dir) {
4251b968477SStefano Zampini     PetscInt n_com,n_dir;
4261b968477SStefano Zampini     n_com = 0;
427d62866d3SStefano Zampini     if (sub_schurs->is_vertices) {
428d62866d3SStefano Zampini       ierr = ISGetLocalSize(sub_schurs->is_vertices,&n_com);CHKERRQ(ierr);
4291b968477SStefano Zampini     }
4301b968477SStefano Zampini     n_dir = 0;
431d62866d3SStefano Zampini     if (sub_schurs->is_dir) {
432d62866d3SStefano Zampini       ierr = ISGetLocalSize(sub_schurs->is_dir,&n_dir);CHKERRQ(ierr);
4331b968477SStefano Zampini     }
43472b8c272SStefano Zampini     if (!deluxe_ctx->n_simple) {
4351b968477SStefano Zampini       deluxe_ctx->n_simple = n_dir + n_com;
4361b968477SStefano Zampini       ierr = PetscMalloc1(deluxe_ctx->n_simple,&deluxe_ctx->idx_simple_B);CHKERRQ(ierr);
437d62866d3SStefano Zampini       if (sub_schurs->is_vertices) {
4389bb4a8caSStefano Zampini         PetscInt       nmap;
4399bb4a8caSStefano Zampini         const PetscInt *idxs;
4401b968477SStefano Zampini 
441d62866d3SStefano Zampini         ierr = ISGetIndices(sub_schurs->is_vertices,&idxs);CHKERRQ(ierr);
4421b968477SStefano Zampini         ierr = ISGlobalToLocalMappingApply(pcis->BtoNmap,IS_GTOLM_DROP,n_com,idxs,&nmap,deluxe_ctx->idx_simple_B);CHKERRQ(ierr);
443af25d912SStefano 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);
444d62866d3SStefano Zampini         ierr = ISRestoreIndices(sub_schurs->is_vertices,&idxs);CHKERRQ(ierr);
4451b968477SStefano Zampini       }
446d62866d3SStefano Zampini       if (sub_schurs->is_dir) {
4471b968477SStefano Zampini         PetscInt       nmap;
4481b968477SStefano Zampini         const PetscInt *idxs;
4491b968477SStefano Zampini 
450d62866d3SStefano Zampini         ierr = ISGetIndices(sub_schurs->is_dir,&idxs);CHKERRQ(ierr);
4511b968477SStefano Zampini         ierr = ISGlobalToLocalMappingApply(pcis->BtoNmap,IS_GTOLM_DROP,n_dir,idxs,&nmap,deluxe_ctx->idx_simple_B+n_com);CHKERRQ(ierr);
452af25d912SStefano 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);
453d62866d3SStefano Zampini         ierr = ISRestoreIndices(sub_schurs->is_dir,&idxs);CHKERRQ(ierr);
4541b968477SStefano Zampini       }
4551b968477SStefano Zampini       ierr = PetscSortInt(deluxe_ctx->n_simple,deluxe_ctx->idx_simple_B);CHKERRQ(ierr);
4569bb4a8caSStefano Zampini     } else {
457af25d912SStefano 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);
45872b8c272SStefano Zampini     }
45972b8c272SStefano Zampini   } else {
460b1b3d7a2SStefano Zampini     deluxe_ctx->n_simple = 0;
4619bb4a8caSStefano Zampini     deluxe_ctx->idx_simple_B = 0;
462b1b3d7a2SStefano Zampini   }
46334a97f8cSStefano Zampini   PetscFunctionReturn(0);
46434a97f8cSStefano Zampini }
46534a97f8cSStefano Zampini 
4665a95e1ceSStefano Zampini static PetscErrorCode PCBDDCScalingSetUp_Deluxe_Private(PC pc)
4675db18549SStefano Zampini {
4685db18549SStefano Zampini   PC_BDDC                *pcbddc=(PC_BDDC*)pc->data;
4695db18549SStefano Zampini   PCBDDCDeluxeScaling    deluxe_ctx=pcbddc->deluxe_ctx;
470b96c3477SStefano Zampini   PCBDDCSubSchurs        sub_schurs = pcbddc->sub_schurs;
47172b8c272SStefano Zampini   PetscScalar            *matdata,*matdata2;
47272b8c272SStefano Zampini   PetscInt               i,max_subset_size,cum,cum2;
47372b8c272SStefano Zampini   const PetscInt         *idxs;
47472b8c272SStefano Zampini   PetscBool              newsetup = PETSC_FALSE;
4755db18549SStefano Zampini   PetscErrorCode         ierr;
4765db18549SStefano Zampini 
4775db18549SStefano Zampini   PetscFunctionBegin;
4785a95e1ceSStefano Zampini   if (!sub_schurs->n_subs) {
4799221af80SStefano Zampini     PetscFunctionReturn(0);
4809221af80SStefano Zampini   }
4819221af80SStefano Zampini 
48272b8c272SStefano Zampini   /* Allocate arrays for subproblems */
48372b8c272SStefano Zampini   if (!deluxe_ctx->seq_n) {
48472b8c272SStefano Zampini     deluxe_ctx->seq_n = sub_schurs->n_subs;
48572b8c272SStefano 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);
48672b8c272SStefano Zampini     newsetup = PETSC_TRUE;
48772b8c272SStefano Zampini   } else if (deluxe_ctx->seq_n != sub_schurs->n_subs) {
48872b8c272SStefano 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);
48972b8c272SStefano Zampini   }
4902f41f7d1SStefano Zampini   /* the change of basis is just a reference to sub_schurs->change (if any) */
4912f41f7d1SStefano Zampini   deluxe_ctx->change         = sub_schurs->change;
49272b8c272SStefano Zampini   deluxe_ctx->change_with_qr = sub_schurs->change_with_qr;
4935db18549SStefano Zampini 
49472b8c272SStefano Zampini   /* Create objects for deluxe */
49572b8c272SStefano Zampini   max_subset_size = 0;
49672b8c272SStefano Zampini   for (i=0;i<sub_schurs->n_subs;i++) {
49772b8c272SStefano Zampini     PetscInt subset_size;
49872b8c272SStefano Zampini     ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr);
49972b8c272SStefano Zampini     max_subset_size = PetscMax(subset_size,max_subset_size);
50072b8c272SStefano Zampini   }
50172b8c272SStefano Zampini   if (newsetup) {
50272b8c272SStefano Zampini     ierr = PetscMalloc1(2*max_subset_size,&deluxe_ctx->workspace);CHKERRQ(ierr);
50372b8c272SStefano Zampini   }
50472b8c272SStefano Zampini   cum = cum2 = 0;
50572b8c272SStefano Zampini   ierr = ISGetIndices(sub_schurs->is_Ej_all,&idxs);CHKERRQ(ierr);
50672b8c272SStefano Zampini   ierr = MatSeqAIJGetArray(sub_schurs->S_Ej_all,&matdata);CHKERRQ(ierr);
50772b8c272SStefano Zampini   ierr = MatSeqAIJGetArray(sub_schurs->sum_S_Ej_all,&matdata2);CHKERRQ(ierr);
50872b8c272SStefano Zampini   for (i=0;i<deluxe_ctx->seq_n;i++) {
50972b8c272SStefano Zampini     PetscInt     subset_size;
510925dfd53SStefano Zampini 
51172b8c272SStefano Zampini     ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr);
51272b8c272SStefano Zampini     if (newsetup) {
51372b8c272SStefano Zampini       IS  sub;
51472b8c272SStefano Zampini       /* work vectors */
51572b8c272SStefano Zampini       ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,1,subset_size,deluxe_ctx->workspace,&deluxe_ctx->seq_work1[i]);CHKERRQ(ierr);
51672b8c272SStefano Zampini       ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,1,subset_size,deluxe_ctx->workspace+subset_size,&deluxe_ctx->seq_work2[i]);CHKERRQ(ierr);
51772b8c272SStefano Zampini 
51872b8c272SStefano Zampini       /* scatters */
51972b8c272SStefano Zampini       ierr = ISCreateGeneral(PETSC_COMM_SELF,subset_size,idxs+cum,PETSC_COPY_VALUES,&sub);CHKERRQ(ierr);
52072b8c272SStefano Zampini       ierr = VecScatterCreate(pcbddc->work_scaling,sub,deluxe_ctx->seq_work1[i],NULL,&deluxe_ctx->seq_scctx[i]);CHKERRQ(ierr);
52172b8c272SStefano Zampini       ierr = ISDestroy(&sub);CHKERRQ(ierr);
522839e9adbSstefano_zampini     }
52372b8c272SStefano Zampini 
52472b8c272SStefano Zampini     /* S_E_j */
525839e9adbSstefano_zampini     ierr = MatDestroy(&deluxe_ctx->seq_mat[i]);CHKERRQ(ierr);
52672b8c272SStefano Zampini     ierr = MatCreateSeqDense(PETSC_COMM_SELF,subset_size,subset_size,matdata+cum2,&deluxe_ctx->seq_mat[i]);CHKERRQ(ierr);
527839e9adbSstefano_zampini 
52872b8c272SStefano Zampini     /* \sum_k S^k_E_j */
52972b8c272SStefano Zampini     ierr = MatDestroy(&deluxe_ctx->seq_mat_inv_sum[i]);CHKERRQ(ierr);
53072b8c272SStefano Zampini     ierr = MatCreateSeqDense(PETSC_COMM_SELF,subset_size,subset_size,matdata2+cum2,&deluxe_ctx->seq_mat_inv_sum[i]);CHKERRQ(ierr);
531*3829e7e6SStefano Zampini     ierr = MatSetOption(deluxe_ctx->seq_mat_inv_sum[i],MAT_SPD,sub_schurs->is_posdef);CHKERRQ(ierr);
532*3829e7e6SStefano Zampini     ierr = MatSetOption(deluxe_ctx->seq_mat_inv_sum[i],MAT_HERMITIAN,sub_schurs->is_hermitian);CHKERRQ(ierr);
533*3829e7e6SStefano Zampini     if (sub_schurs->is_hermitian) {
53472b8c272SStefano Zampini       ierr = MatCholeskyFactor(deluxe_ctx->seq_mat_inv_sum[i],NULL,NULL);CHKERRQ(ierr);
535d62866d3SStefano Zampini     } else {
53672b8c272SStefano Zampini       ierr = MatLUFactor(deluxe_ctx->seq_mat_inv_sum[i],NULL,NULL,NULL);CHKERRQ(ierr);
537d62866d3SStefano Zampini     }
538839e9adbSstefano_zampini     if (pcbddc->deluxe_singlemat) {
539839e9adbSstefano_zampini       Mat X,Y;
540*3829e7e6SStefano Zampini       if (!sub_schurs->is_hermitian) {
541839e9adbSstefano_zampini         ierr = MatTranspose(deluxe_ctx->seq_mat[i],MAT_INITIAL_MATRIX,&X);CHKERRQ(ierr);
542839e9adbSstefano_zampini       } else {
543839e9adbSstefano_zampini         ierr = PetscObjectReference((PetscObject)deluxe_ctx->seq_mat[i]);CHKERRQ(ierr);
544839e9adbSstefano_zampini         X    = deluxe_ctx->seq_mat[i];
545839e9adbSstefano_zampini       }
546839e9adbSstefano_zampini       ierr = MatDuplicate(X,MAT_DO_NOT_COPY_VALUES,&Y);CHKERRQ(ierr);
547*3829e7e6SStefano Zampini       if (!sub_schurs->is_hermitian) {
548839e9adbSstefano_zampini         ierr = PCBDDCMatTransposeMatSolve_SeqDense(deluxe_ctx->seq_mat_inv_sum[i],X,Y);CHKERRQ(ierr);
549839e9adbSstefano_zampini       } else {
550839e9adbSstefano_zampini         ierr = MatMatSolve(deluxe_ctx->seq_mat_inv_sum[i],X,Y);CHKERRQ(ierr);
551839e9adbSstefano_zampini       }
552729c86d3Sstefano_zampini 
553839e9adbSstefano_zampini       ierr = MatDestroy(&deluxe_ctx->seq_mat_inv_sum[i]);CHKERRQ(ierr);
554839e9adbSstefano_zampini       ierr = MatDestroy(&deluxe_ctx->seq_mat[i]);CHKERRQ(ierr);
555839e9adbSstefano_zampini       ierr = MatDestroy(&X);CHKERRQ(ierr);
5562f41f7d1SStefano Zampini       if (deluxe_ctx->change) {
5572f41f7d1SStefano Zampini         Mat C,CY;
5582f41f7d1SStefano Zampini 
5592f41f7d1SStefano Zampini         if (!deluxe_ctx->change_with_qr) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Only QR based change of basis");
5602f41f7d1SStefano Zampini         ierr = KSPGetOperators(deluxe_ctx->change[i],&C,NULL);CHKERRQ(ierr);
5612f41f7d1SStefano Zampini         ierr = MatMatMult(C,Y,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&CY);CHKERRQ(ierr);
5622f41f7d1SStefano Zampini         ierr = MatMatTransposeMult(CY,C,MAT_REUSE_MATRIX,PETSC_DEFAULT,&Y);CHKERRQ(ierr);
5632f41f7d1SStefano Zampini         ierr = MatDestroy(&CY);CHKERRQ(ierr);
5642f41f7d1SStefano Zampini       }
56503dfb2d7SStefano Zampini       ierr = MatTranspose(Y,MAT_INPLACE_MATRIX,&Y);CHKERRQ(ierr);
566839e9adbSstefano_zampini       deluxe_ctx->seq_mat[i] = Y;
56772b8c272SStefano Zampini     }
56872b8c272SStefano Zampini     cum += subset_size;
56972b8c272SStefano Zampini     cum2 += subset_size*subset_size;
57072b8c272SStefano Zampini   }
57172b8c272SStefano Zampini   ierr = ISRestoreIndices(sub_schurs->is_Ej_all,&idxs);CHKERRQ(ierr);
57272b8c272SStefano Zampini   ierr = MatSeqAIJRestoreArray(sub_schurs->S_Ej_all,&matdata);CHKERRQ(ierr);
57372b8c272SStefano Zampini   ierr = MatSeqAIJRestoreArray(sub_schurs->sum_S_Ej_all,&matdata2);CHKERRQ(ierr);
5742f41f7d1SStefano Zampini   if (pcbddc->deluxe_singlemat) {
5752f41f7d1SStefano Zampini     deluxe_ctx->change         = NULL;
5762f41f7d1SStefano Zampini     deluxe_ctx->change_with_qr = PETSC_FALSE;
5772f41f7d1SStefano Zampini   }
5785db18549SStefano Zampini 
57972b8c272SStefano Zampini   if (deluxe_ctx->change && !deluxe_ctx->change_with_qr) {
58072b8c272SStefano Zampini     for (i=0;i<deluxe_ctx->seq_n;i++) {
58172b8c272SStefano Zampini       if (newsetup) {
58272b8c272SStefano Zampini         PC pc;
58372b8c272SStefano Zampini 
58472b8c272SStefano Zampini         ierr = KSPGetPC(deluxe_ctx->change[i],&pc);CHKERRQ(ierr);
58572b8c272SStefano Zampini         ierr = PCSetType(pc,PCLU);CHKERRQ(ierr);
58672b8c272SStefano Zampini         ierr = KSPSetFromOptions(deluxe_ctx->change[i]);CHKERRQ(ierr);
58772b8c272SStefano Zampini       }
58872b8c272SStefano Zampini       ierr = KSPSetUp(deluxe_ctx->change[i]);CHKERRQ(ierr);
58972b8c272SStefano Zampini     }
59016e386b8SStefano Zampini   }
5915db18549SStefano Zampini   PetscFunctionReturn(0);
5925db18549SStefano Zampini }
593