xref: /petsc/src/ksp/pc/impls/bddc/bddcscalingbasic.c (revision 2c71b3e237ead271e4f3aa1505f92bf476e3413d)
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;
171683a169SBarry Smith   const PetscScalar *b;
181683a169SBarry Smith   PetscScalar       *x;
19839e9adbSstefano_zampini   PetscInt          n;
20839e9adbSstefano_zampini   PetscBLASInt      nrhs,info,m;
21839e9adbSstefano_zampini   PetscBool         flg;
22839e9adbSstefano_zampini 
23839e9adbSstefano_zampini   PetscFunctionBegin;
24839e9adbSstefano_zampini   ierr = PetscBLASIntCast(A->rmap->n,&m);CHKERRQ(ierr);
25839e9adbSstefano_zampini   ierr = PetscObjectTypeCompareAny((PetscObject)B,&flg,MATSEQDENSE,MATMPIDENSE,NULL);CHKERRQ(ierr);
26*2c71b3e2SJacob Faibussowitsch   PetscCheckFalse(!flg,PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONG,"Matrix B must be MATDENSE matrix");
27839e9adbSstefano_zampini   ierr = PetscObjectTypeCompareAny((PetscObject)X,&flg,MATSEQDENSE,MATMPIDENSE,NULL);CHKERRQ(ierr);
28*2c71b3e2SJacob Faibussowitsch   PetscCheckFalse(!flg,PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONG,"Matrix X must be MATDENSE matrix");
29839e9adbSstefano_zampini 
30839e9adbSstefano_zampini   ierr = MatGetSize(B,NULL,&n);CHKERRQ(ierr);
31839e9adbSstefano_zampini   ierr = PetscBLASIntCast(n,&nrhs);CHKERRQ(ierr);
321683a169SBarry Smith   ierr = MatDenseGetArrayRead(B,&b);CHKERRQ(ierr);
33839e9adbSstefano_zampini   ierr = MatDenseGetArray(X,&x);CHKERRQ(ierr);
34580bdb30SBarry Smith   ierr = PetscArraycpy(x,b,m*nrhs);CHKERRQ(ierr);
351683a169SBarry Smith   ierr = MatDenseRestoreArrayRead(B,&b);CHKERRQ(ierr);
36839e9adbSstefano_zampini 
37839e9adbSstefano_zampini   if (A->factortype == MAT_FACTOR_LU) {
38839e9adbSstefano_zampini     PetscStackCallBLAS("LAPACKgetrs",LAPACKgetrs_("T",&m,&nrhs,mat->v,&mat->lda,mat->pivots,x,&m,&info));
39*2c71b3e2SJacob Faibussowitsch     PetscCheckFalse(info,PETSC_COMM_SELF,PETSC_ERR_LIB,"GETRS - Bad solve");
403829e7e6SStefano Zampini   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Only LU factor supported");
41839e9adbSstefano_zampini 
42839e9adbSstefano_zampini   ierr = MatDenseRestoreArray(X,&x);CHKERRQ(ierr);
43839e9adbSstefano_zampini   ierr = PetscLogFlops(nrhs*(2.0*m*m - m));CHKERRQ(ierr);
44839e9adbSstefano_zampini   PetscFunctionReturn(0);
45839e9adbSstefano_zampini }
46839e9adbSstefano_zampini 
47674ae819SStefano Zampini static PetscErrorCode PCBDDCScalingExtension_Basic(PC pc, Vec local_interface_vector, Vec global_vector)
48674ae819SStefano Zampini {
49674ae819SStefano Zampini   PC_IS*         pcis = (PC_IS*)pc->data;
50674ae819SStefano Zampini   PC_BDDC*       pcbddc = (PC_BDDC*)pc->data;
51674ae819SStefano Zampini   PetscErrorCode ierr;
52674ae819SStefano Zampini 
53674ae819SStefano Zampini   PetscFunctionBegin;
54674ae819SStefano Zampini   /* Apply partition of unity */
55674ae819SStefano Zampini   ierr = VecPointwiseMult(pcbddc->work_scaling,pcis->D,local_interface_vector);CHKERRQ(ierr);
56674ae819SStefano Zampini   ierr = VecSet(global_vector,0.0);CHKERRQ(ierr);
57674ae819SStefano Zampini   ierr = VecScatterBegin(pcis->global_to_B,pcbddc->work_scaling,global_vector,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
58674ae819SStefano Zampini   ierr = VecScatterEnd(pcis->global_to_B,pcbddc->work_scaling,global_vector,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
59674ae819SStefano Zampini   PetscFunctionReturn(0);
60674ae819SStefano Zampini }
61674ae819SStefano Zampini 
62674ae819SStefano Zampini static PetscErrorCode PCBDDCScalingExtension_Deluxe(PC pc, Vec x, Vec y)
63674ae819SStefano Zampini {
64674ae819SStefano Zampini   PC_IS*              pcis=(PC_IS*)pc->data;
65674ae819SStefano Zampini   PC_BDDC*            pcbddc=(PC_BDDC*)pc->data;
66674ae819SStefano Zampini   PCBDDCDeluxeScaling deluxe_ctx = pcbddc->deluxe_ctx;
67674ae819SStefano Zampini   PetscErrorCode      ierr;
68674ae819SStefano Zampini 
69674ae819SStefano Zampini   PetscFunctionBegin;
705a95e1ceSStefano Zampini   ierr = VecSet(pcbddc->work_scaling,0.0);CHKERRQ(ierr);
715a95e1ceSStefano Zampini   ierr = VecSet(y,0.0);CHKERRQ(ierr);
725a95e1ceSStefano Zampini   if (deluxe_ctx->n_simple) { /* scale deluxe vertices using diagonal scaling */
735a95e1ceSStefano Zampini     PetscInt          i;
742b095fd8SStefano Zampini     const PetscScalar *array_x,*array_D;
752b095fd8SStefano Zampini     PetscScalar       *array;
762b095fd8SStefano Zampini     ierr = VecGetArrayRead(x,&array_x);CHKERRQ(ierr);
772b095fd8SStefano Zampini     ierr = VecGetArrayRead(pcis->D,&array_D);CHKERRQ(ierr);
78674ae819SStefano Zampini     ierr = VecGetArray(pcbddc->work_scaling,&array);CHKERRQ(ierr);
79674ae819SStefano Zampini     for (i=0;i<deluxe_ctx->n_simple;i++) {
80674ae819SStefano Zampini       array[deluxe_ctx->idx_simple_B[i]] = array_x[deluxe_ctx->idx_simple_B[i]]*array_D[deluxe_ctx->idx_simple_B[i]];
81674ae819SStefano Zampini     }
82674ae819SStefano Zampini     ierr = VecRestoreArray(pcbddc->work_scaling,&array);CHKERRQ(ierr);
832b095fd8SStefano Zampini     ierr = VecRestoreArrayRead(pcis->D,&array_D);CHKERRQ(ierr);
842b095fd8SStefano Zampini     ierr = VecRestoreArrayRead(x,&array_x);CHKERRQ(ierr);
8534a97f8cSStefano Zampini   }
86ac632422SStefano Zampini   /* sequential part : all problems and Schur applications collapsed into a single matrix vector multiplication or a matvec and a solve */
875a95e1ceSStefano Zampini   if (deluxe_ctx->seq_mat) {
8872b8c272SStefano Zampini     PetscInt i;
8972b8c272SStefano Zampini     for (i=0;i<deluxe_ctx->seq_n;i++) {
9016e386b8SStefano Zampini       if (deluxe_ctx->change) {
9172b8c272SStefano Zampini         ierr = VecScatterBegin(deluxe_ctx->seq_scctx[i],x,deluxe_ctx->seq_work2[i],INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
9272b8c272SStefano Zampini         ierr = VecScatterEnd(deluxe_ctx->seq_scctx[i],x,deluxe_ctx->seq_work2[i],INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
9372b8c272SStefano Zampini         if (deluxe_ctx->change_with_qr) {
9472b8c272SStefano Zampini           Mat change;
9572b8c272SStefano Zampini 
9672b8c272SStefano Zampini           ierr = KSPGetOperators(deluxe_ctx->change[i],&change,NULL);CHKERRQ(ierr);
9772b8c272SStefano Zampini           ierr = MatMultTranspose(change,deluxe_ctx->seq_work2[i],deluxe_ctx->seq_work1[i]);CHKERRQ(ierr);
9872b8c272SStefano Zampini         } else {
9972b8c272SStefano Zampini           ierr = KSPSolve(deluxe_ctx->change[i],deluxe_ctx->seq_work2[i],deluxe_ctx->seq_work1[i]);CHKERRQ(ierr);
10016e386b8SStefano Zampini         }
10172b8c272SStefano Zampini       } else {
10272b8c272SStefano Zampini         ierr = VecScatterBegin(deluxe_ctx->seq_scctx[i],x,deluxe_ctx->seq_work1[i],INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
10372b8c272SStefano Zampini         ierr = VecScatterEnd(deluxe_ctx->seq_scctx[i],x,deluxe_ctx->seq_work1[i],INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
10472b8c272SStefano Zampini       }
10572b8c272SStefano Zampini       ierr = MatMultTranspose(deluxe_ctx->seq_mat[i],deluxe_ctx->seq_work1[i],deluxe_ctx->seq_work2[i]);CHKERRQ(ierr);
10672b8c272SStefano Zampini       if (deluxe_ctx->seq_mat_inv_sum[i]) {
10772b8c272SStefano Zampini         PetscScalar *x;
10872b8c272SStefano Zampini 
10972b8c272SStefano Zampini         ierr = VecGetArray(deluxe_ctx->seq_work2[i],&x);CHKERRQ(ierr);
11072b8c272SStefano Zampini         ierr = VecPlaceArray(deluxe_ctx->seq_work1[i],x);CHKERRQ(ierr);
11172b8c272SStefano Zampini         ierr = VecRestoreArray(deluxe_ctx->seq_work2[i],&x);CHKERRQ(ierr);
11272b8c272SStefano Zampini         ierr = MatSolveTranspose(deluxe_ctx->seq_mat_inv_sum[i],deluxe_ctx->seq_work1[i],deluxe_ctx->seq_work2[i]);CHKERRQ(ierr);
11372b8c272SStefano Zampini         ierr = VecResetArray(deluxe_ctx->seq_work1[i]);CHKERRQ(ierr);
114ac632422SStefano Zampini       }
11516e386b8SStefano Zampini       if (deluxe_ctx->change) {
11616e386b8SStefano Zampini         Mat change;
11716e386b8SStefano Zampini 
11872b8c272SStefano Zampini         ierr = KSPGetOperators(deluxe_ctx->change[i],&change,NULL);CHKERRQ(ierr);
11972b8c272SStefano Zampini         ierr = MatMult(change,deluxe_ctx->seq_work2[i],deluxe_ctx->seq_work1[i]);CHKERRQ(ierr);
12072b8c272SStefano Zampini         ierr = VecScatterBegin(deluxe_ctx->seq_scctx[i],deluxe_ctx->seq_work1[i],pcbddc->work_scaling,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
12172b8c272SStefano Zampini         ierr = VecScatterEnd(deluxe_ctx->seq_scctx[i],deluxe_ctx->seq_work1[i],pcbddc->work_scaling,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
12216e386b8SStefano Zampini       } else {
12372b8c272SStefano Zampini         ierr = VecScatterBegin(deluxe_ctx->seq_scctx[i],deluxe_ctx->seq_work2[i],pcbddc->work_scaling,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
12472b8c272SStefano Zampini         ierr = VecScatterEnd(deluxe_ctx->seq_scctx[i],deluxe_ctx->seq_work2[i],pcbddc->work_scaling,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
12572b8c272SStefano Zampini       }
126674ae819SStefano Zampini     }
12716e386b8SStefano Zampini   }
128674ae819SStefano Zampini   /* put local boundary part in global vector */
129674ae819SStefano Zampini   ierr = VecScatterBegin(pcis->global_to_B,pcbddc->work_scaling,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
130674ae819SStefano Zampini   ierr = VecScatterEnd(pcis->global_to_B,pcbddc->work_scaling,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
131674ae819SStefano Zampini   PetscFunctionReturn(0);
132674ae819SStefano Zampini }
133674ae819SStefano Zampini 
134674ae819SStefano Zampini PetscErrorCode PCBDDCScalingExtension(PC pc, Vec local_interface_vector, Vec global_vector)
135674ae819SStefano Zampini {
136674ae819SStefano Zampini   PC_BDDC        *pcbddc=(PC_BDDC*)pc->data;
137674ae819SStefano Zampini   PetscErrorCode ierr;
138674ae819SStefano Zampini 
139674ae819SStefano Zampini   PetscFunctionBegin;
140674ae819SStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
141674ae819SStefano Zampini   PetscValidHeaderSpecific(local_interface_vector,VEC_CLASSID,2);
142674ae819SStefano Zampini   PetscValidHeaderSpecific(global_vector,VEC_CLASSID,3);
143*2c71b3e2SJacob Faibussowitsch   PetscCheckFalse(local_interface_vector == pcbddc->work_scaling,PETSC_COMM_SELF,PETSC_ERR_SUP,"Local vector cannot be pcbddc->work_scaling!");
144163d334eSBarry Smith   ierr = PetscUseMethod(pc,"PCBDDCScalingExtension_C",(PC,Vec,Vec),(pc,local_interface_vector,global_vector));CHKERRQ(ierr);
145674ae819SStefano Zampini   PetscFunctionReturn(0);
146674ae819SStefano Zampini }
147674ae819SStefano Zampini 
148674ae819SStefano Zampini static PetscErrorCode PCBDDCScalingRestriction_Basic(PC pc, Vec global_vector, Vec local_interface_vector)
149674ae819SStefano Zampini {
150674ae819SStefano Zampini   PetscErrorCode ierr;
151674ae819SStefano Zampini   PC_IS          *pcis = (PC_IS*)pc->data;
152674ae819SStefano Zampini 
153674ae819SStefano Zampini   PetscFunctionBegin;
154674ae819SStefano Zampini   ierr = VecScatterBegin(pcis->global_to_B,global_vector,local_interface_vector,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
155674ae819SStefano Zampini   ierr = VecScatterEnd(pcis->global_to_B,global_vector,local_interface_vector,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
156674ae819SStefano Zampini   /* Apply partition of unity */
157674ae819SStefano Zampini   ierr = VecPointwiseMult(local_interface_vector,pcis->D,local_interface_vector);CHKERRQ(ierr);
158674ae819SStefano Zampini   PetscFunctionReturn(0);
159674ae819SStefano Zampini }
160674ae819SStefano Zampini 
161674ae819SStefano Zampini static PetscErrorCode PCBDDCScalingRestriction_Deluxe(PC pc, Vec x, Vec y)
162674ae819SStefano Zampini {
163674ae819SStefano Zampini   PC_IS*              pcis=(PC_IS*)pc->data;
164674ae819SStefano Zampini   PC_BDDC*            pcbddc=(PC_BDDC*)pc->data;
165674ae819SStefano Zampini   PCBDDCDeluxeScaling deluxe_ctx = pcbddc->deluxe_ctx;
166674ae819SStefano Zampini   PetscErrorCode      ierr;
167674ae819SStefano Zampini 
168674ae819SStefano Zampini   PetscFunctionBegin;
169674ae819SStefano Zampini   /* get local boundary part of global vector */
170674ae819SStefano Zampini   ierr = VecScatterBegin(pcis->global_to_B,x,y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
171674ae819SStefano Zampini   ierr = VecScatterEnd(pcis->global_to_B,x,y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1725a95e1ceSStefano Zampini   if (deluxe_ctx->n_simple) { /* scale deluxe vertices using diagonal scaling */
1735a95e1ceSStefano Zampini     PetscInt          i;
1742b095fd8SStefano Zampini     PetscScalar       *array_y;
1752b095fd8SStefano Zampini     const PetscScalar *array_D;
176674ae819SStefano Zampini     ierr = VecGetArray(y,&array_y);CHKERRQ(ierr);
1772b095fd8SStefano Zampini     ierr = VecGetArrayRead(pcis->D,&array_D);CHKERRQ(ierr);
178674ae819SStefano Zampini     for (i=0;i<deluxe_ctx->n_simple;i++) {
179674ae819SStefano Zampini       array_y[deluxe_ctx->idx_simple_B[i]] *= array_D[deluxe_ctx->idx_simple_B[i]];
180674ae819SStefano Zampini     }
1812b095fd8SStefano Zampini     ierr = VecRestoreArrayRead(pcis->D,&array_D);CHKERRQ(ierr);
182674ae819SStefano Zampini     ierr = VecRestoreArray(y,&array_y);CHKERRQ(ierr);
18334a97f8cSStefano Zampini   }
184729c86d3Sstefano_zampini   /* sequential part : all problems and Schur applications collapsed into a single matrix vector multiplication or a matvec and a solve */
1855a95e1ceSStefano Zampini   if (deluxe_ctx->seq_mat) {
18672b8c272SStefano Zampini     PetscInt i;
18772b8c272SStefano Zampini     for (i=0;i<deluxe_ctx->seq_n;i++) {
18816e386b8SStefano Zampini       if (deluxe_ctx->change) {
18916e386b8SStefano Zampini         Mat change;
19016e386b8SStefano Zampini 
19172b8c272SStefano Zampini         ierr = VecScatterBegin(deluxe_ctx->seq_scctx[i],y,deluxe_ctx->seq_work2[i],INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
19272b8c272SStefano Zampini         ierr = VecScatterEnd(deluxe_ctx->seq_scctx[i],y,deluxe_ctx->seq_work2[i],INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
19372b8c272SStefano Zampini         ierr = KSPGetOperators(deluxe_ctx->change[i],&change,NULL);CHKERRQ(ierr);
19472b8c272SStefano Zampini         ierr = MatMultTranspose(change,deluxe_ctx->seq_work2[i],deluxe_ctx->seq_work1[i]);CHKERRQ(ierr);
19516e386b8SStefano Zampini       } else {
19672b8c272SStefano Zampini         ierr = VecScatterBegin(deluxe_ctx->seq_scctx[i],y,deluxe_ctx->seq_work1[i],INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
19772b8c272SStefano Zampini         ierr = VecScatterEnd(deluxe_ctx->seq_scctx[i],y,deluxe_ctx->seq_work1[i],INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
19816e386b8SStefano Zampini       }
19972b8c272SStefano Zampini       if (deluxe_ctx->seq_mat_inv_sum[i]) {
20072b8c272SStefano Zampini         PetscScalar *x;
20172b8c272SStefano Zampini 
20272b8c272SStefano Zampini         ierr = VecGetArray(deluxe_ctx->seq_work1[i],&x);CHKERRQ(ierr);
20372b8c272SStefano Zampini         ierr = VecPlaceArray(deluxe_ctx->seq_work2[i],x);CHKERRQ(ierr);
20472b8c272SStefano Zampini         ierr = VecRestoreArray(deluxe_ctx->seq_work1[i],&x);CHKERRQ(ierr);
20572b8c272SStefano Zampini         ierr = MatSolve(deluxe_ctx->seq_mat_inv_sum[i],deluxe_ctx->seq_work1[i],deluxe_ctx->seq_work2[i]);CHKERRQ(ierr);
20672b8c272SStefano Zampini         ierr = VecResetArray(deluxe_ctx->seq_work2[i]);CHKERRQ(ierr);
207ac632422SStefano Zampini       }
20872b8c272SStefano Zampini       ierr = MatMult(deluxe_ctx->seq_mat[i],deluxe_ctx->seq_work1[i],deluxe_ctx->seq_work2[i]);CHKERRQ(ierr);
20916e386b8SStefano Zampini       if (deluxe_ctx->change) {
21072b8c272SStefano Zampini         if (deluxe_ctx->change_with_qr) {
21172b8c272SStefano Zampini           Mat change;
21272b8c272SStefano Zampini 
21372b8c272SStefano Zampini           ierr = KSPGetOperators(deluxe_ctx->change[i],&change,NULL);CHKERRQ(ierr);
21472b8c272SStefano Zampini           ierr = MatMult(change,deluxe_ctx->seq_work2[i],deluxe_ctx->seq_work1[i]);CHKERRQ(ierr);
21572b8c272SStefano Zampini         } else {
21672b8c272SStefano Zampini           ierr = KSPSolveTranspose(deluxe_ctx->change[i],deluxe_ctx->seq_work2[i],deluxe_ctx->seq_work1[i]);CHKERRQ(ierr);
21716e386b8SStefano Zampini         }
21872b8c272SStefano Zampini         ierr = VecScatterBegin(deluxe_ctx->seq_scctx[i],deluxe_ctx->seq_work1[i],y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
21972b8c272SStefano Zampini         ierr = VecScatterEnd(deluxe_ctx->seq_scctx[i],deluxe_ctx->seq_work1[i],y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
22072b8c272SStefano Zampini       } else {
22172b8c272SStefano Zampini         ierr = VecScatterBegin(deluxe_ctx->seq_scctx[i],deluxe_ctx->seq_work2[i],y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
22272b8c272SStefano Zampini         ierr = VecScatterEnd(deluxe_ctx->seq_scctx[i],deluxe_ctx->seq_work2[i],y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
22372b8c272SStefano Zampini       }
22472b8c272SStefano Zampini     }
225674ae819SStefano Zampini   }
226674ae819SStefano Zampini   PetscFunctionReturn(0);
227674ae819SStefano Zampini }
228674ae819SStefano Zampini 
229674ae819SStefano Zampini PetscErrorCode PCBDDCScalingRestriction(PC pc, Vec global_vector, Vec local_interface_vector)
230674ae819SStefano Zampini {
231674ae819SStefano Zampini   PC_BDDC        *pcbddc=(PC_BDDC*)pc->data;
232674ae819SStefano Zampini   PetscErrorCode ierr;
233674ae819SStefano Zampini 
234674ae819SStefano Zampini   PetscFunctionBegin;
235674ae819SStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
236674ae819SStefano Zampini   PetscValidHeaderSpecific(global_vector,VEC_CLASSID,2);
237674ae819SStefano Zampini   PetscValidHeaderSpecific(local_interface_vector,VEC_CLASSID,3);
238*2c71b3e2SJacob Faibussowitsch   PetscCheckFalse(local_interface_vector == pcbddc->work_scaling,PETSC_COMM_SELF,PETSC_ERR_SUP,"Local vector cannot be pcbddc->work_scaling!");
239163d334eSBarry Smith   ierr = PetscUseMethod(pc,"PCBDDCScalingRestriction_C",(PC,Vec,Vec),(pc,global_vector,local_interface_vector));CHKERRQ(ierr);
240674ae819SStefano Zampini   PetscFunctionReturn(0);
241674ae819SStefano Zampini }
242674ae819SStefano Zampini 
243674ae819SStefano Zampini PetscErrorCode PCBDDCScalingSetUp(PC pc)
244674ae819SStefano Zampini {
245674ae819SStefano Zampini   PC_IS*         pcis=(PC_IS*)pc->data;
246674ae819SStefano Zampini   PC_BDDC*       pcbddc=(PC_BDDC*)pc->data;
247674ae819SStefano Zampini   PetscErrorCode ierr;
248674ae819SStefano Zampini 
249674ae819SStefano Zampini   PetscFunctionBegin;
250674ae819SStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
25143371fb9SStefano Zampini   ierr = PetscLogEventBegin(PC_BDDC_Scaling[pcbddc->current_level],pc,0,0,0);CHKERRQ(ierr);
252674ae819SStefano Zampini   /* create work vector for the operator */
25334a97f8cSStefano Zampini   ierr = VecDestroy(&pcbddc->work_scaling);CHKERRQ(ierr);
254674ae819SStefano Zampini   ierr = VecDuplicate(pcis->vec1_B,&pcbddc->work_scaling);CHKERRQ(ierr);
25534a97f8cSStefano Zampini   /* always rebuild pcis->D */
25628d874f6SStefano Zampini   if (pcis->use_stiffness_scaling) {
25715f235b8SStefano Zampini     PetscScalar *a;
25815f235b8SStefano Zampini     PetscInt    i,n;
25915f235b8SStefano Zampini 
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);
2638d25b668SStefano Zampini     ierr = VecAbs(pcis->D);CHKERRQ(ierr);
26415f235b8SStefano Zampini     ierr = VecGetLocalSize(pcis->D,&n);CHKERRQ(ierr);
26515f235b8SStefano Zampini     ierr = VecGetArray(pcis->D,&a);CHKERRQ(ierr);
26615f235b8SStefano Zampini     for (i=0;i<n;i++) if (PetscAbsScalar(a[i])<PETSC_SMALL) a[i] = 1.0;
26715f235b8SStefano Zampini     ierr = VecRestoreArray(pcis->D,&a);CHKERRQ(ierr);
268674ae819SStefano Zampini   }
269674ae819SStefano Zampini   ierr = VecSet(pcis->vec1_global,0.0);CHKERRQ(ierr);
2708d25b668SStefano Zampini   ierr = VecScatterBegin(pcis->global_to_B,pcis->D,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
2718d25b668SStefano Zampini   ierr = VecScatterEnd(pcis->global_to_B,pcis->D,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
272674ae819SStefano Zampini   ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_global,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
273674ae819SStefano Zampini   ierr = VecScatterEnd(pcis->global_to_B,pcis->vec1_global,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
274674ae819SStefano Zampini   ierr = VecPointwiseDivide(pcis->D,pcis->D,pcis->vec1_B);CHKERRQ(ierr);
275674ae819SStefano Zampini   /* now setup */
276681e7c04SStefano Zampini   if (pcbddc->use_deluxe_scaling) {
27734a97f8cSStefano Zampini     if (!pcbddc->deluxe_ctx) {
27834a97f8cSStefano Zampini       ierr = PCBDDCScalingCreate_Deluxe(pc);CHKERRQ(ierr);
27934a97f8cSStefano Zampini     }
28034a97f8cSStefano Zampini     ierr = PCBDDCScalingSetUp_Deluxe(pc);CHKERRQ(ierr);
281674ae819SStefano Zampini     ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCScalingRestriction_C",PCBDDCScalingRestriction_Deluxe);CHKERRQ(ierr);
282674ae819SStefano Zampini     ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCScalingExtension_C",PCBDDCScalingExtension_Deluxe);CHKERRQ(ierr);
283674ae819SStefano Zampini   } else {
284674ae819SStefano Zampini     ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCScalingRestriction_C",PCBDDCScalingRestriction_Basic);CHKERRQ(ierr);
285674ae819SStefano Zampini     ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCScalingExtension_C",PCBDDCScalingExtension_Basic);CHKERRQ(ierr);
286674ae819SStefano Zampini   }
2878ead10e4SStefano Zampini   ierr = PetscLogEventEnd(PC_BDDC_Scaling[pcbddc->current_level],pc,0,0,0);CHKERRQ(ierr);
28834a97f8cSStefano Zampini 
289674ae819SStefano Zampini   /* test */
290674ae819SStefano Zampini   if (pcbddc->dbg_flag) {
29172b8c272SStefano Zampini     Mat         B0_B = NULL;
29272b8c272SStefano Zampini     Vec         B0_Bv = NULL, B0_Bv2 = NULL;
29334a97f8cSStefano Zampini     Vec         vec2_global;
294674ae819SStefano Zampini     PetscViewer viewer = pcbddc->dbg_viewer;
29534a97f8cSStefano Zampini     PetscReal   error;
296674ae819SStefano Zampini 
297674ae819SStefano Zampini     /* extension -> from local to parallel */
29834a97f8cSStefano Zampini     ierr = VecSet(pcis->vec1_global,0.0);CHKERRQ(ierr);
29934a97f8cSStefano Zampini     ierr = VecSetRandom(pcis->vec1_B,NULL);CHKERRQ(ierr);
30034a97f8cSStefano Zampini     ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_B,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
30134a97f8cSStefano Zampini     ierr = VecScatterEnd(pcis->global_to_B,pcis->vec1_B,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
30234a97f8cSStefano Zampini     ierr = VecDuplicate(pcis->vec1_global,&vec2_global);CHKERRQ(ierr);
30334a97f8cSStefano Zampini     ierr = VecCopy(pcis->vec1_global,vec2_global);CHKERRQ(ierr);
304674ae819SStefano Zampini     ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_global,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
305674ae819SStefano Zampini     ierr = VecScatterEnd(pcis->global_to_B,pcis->vec1_global,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
30672b8c272SStefano Zampini     if (pcbddc->benign_n) {
30772b8c272SStefano Zampini       IS is_dummy;
30872b8c272SStefano Zampini 
30972b8c272SStefano Zampini       ierr = ISCreateStride(PETSC_COMM_SELF,pcbddc->benign_n,0,1,&is_dummy);CHKERRQ(ierr);
3107dae84e0SHong Zhang       ierr = MatCreateSubMatrix(pcbddc->benign_B0,is_dummy,pcis->is_B_local,MAT_INITIAL_MATRIX,&B0_B);CHKERRQ(ierr);
31172b8c272SStefano Zampini       ierr = ISDestroy(&is_dummy);CHKERRQ(ierr);
31272b8c272SStefano Zampini       ierr = MatCreateVecs(B0_B,NULL,&B0_Bv);CHKERRQ(ierr);
31372b8c272SStefano Zampini       ierr = VecDuplicate(B0_Bv,&B0_Bv2);CHKERRQ(ierr);
31472b8c272SStefano Zampini       ierr = MatMult(B0_B,pcis->vec1_B,B0_Bv);CHKERRQ(ierr);
31572b8c272SStefano Zampini     }
316674ae819SStefano Zampini     ierr = PCBDDCScalingExtension(pc,pcis->vec1_B,pcis->vec1_global);CHKERRQ(ierr);
31772b8c272SStefano Zampini     if (pcbddc->benign_saddle_point) {
31872b8c272SStefano Zampini       PetscReal errorl = 0.;
31972b8c272SStefano Zampini       ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_global,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
32072b8c272SStefano Zampini       ierr = VecScatterEnd(pcis->global_to_B,pcis->vec1_global,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
32172b8c272SStefano Zampini       if (pcbddc->benign_n) {
32272b8c272SStefano Zampini         ierr = MatMult(B0_B,pcis->vec1_B,B0_Bv2);CHKERRQ(ierr);
32372b8c272SStefano Zampini         ierr = VecAXPY(B0_Bv,-1.0,B0_Bv2);CHKERRQ(ierr);
32472b8c272SStefano Zampini         ierr = VecNorm(B0_Bv,NORM_INFINITY,&errorl);CHKERRQ(ierr);
32572b8c272SStefano Zampini       }
326ffc4695bSBarry Smith       ierr = MPI_Allreduce(&errorl,&error,1,MPIU_REAL,MPI_SUM,PetscObjectComm((PetscObject)pc));CHKERRMPI(ierr);
32772b8c272SStefano Zampini       ierr = PetscViewerASCIIPrintf(viewer,"Error benign extension %1.14e\n",error);CHKERRQ(ierr);
32872b8c272SStefano Zampini     }
32934a97f8cSStefano Zampini     ierr = VecAXPY(pcis->vec1_global,-1.0,vec2_global);CHKERRQ(ierr);
33034a97f8cSStefano Zampini     ierr = VecNorm(pcis->vec1_global,NORM_INFINITY,&error);CHKERRQ(ierr);
331674ae819SStefano Zampini     ierr = PetscViewerASCIIPrintf(viewer,"Error scaling extension %1.14e\n",error);CHKERRQ(ierr);
33234a97f8cSStefano Zampini     ierr = VecDestroy(&vec2_global);CHKERRQ(ierr);
33334a97f8cSStefano Zampini 
334674ae819SStefano Zampini     /* restriction -> from parallel to local */
335674ae819SStefano Zampini     ierr = VecSet(pcis->vec1_global,0.0);CHKERRQ(ierr);
33634a97f8cSStefano Zampini     ierr = VecSetRandom(pcis->vec1_B,NULL);CHKERRQ(ierr);
337674ae819SStefano Zampini     ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_B,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
338674ae819SStefano Zampini     ierr = VecScatterEnd(pcis->global_to_B,pcis->vec1_B,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
33934a97f8cSStefano Zampini     ierr = PCBDDCScalingRestriction(pc,pcis->vec1_global,pcis->vec1_B);CHKERRQ(ierr);
34034a97f8cSStefano Zampini     ierr = VecScale(pcis->vec1_B,-1.0);CHKERRQ(ierr);
34134a97f8cSStefano Zampini     ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_B,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
34234a97f8cSStefano Zampini     ierr = VecScatterEnd(pcis->global_to_B,pcis->vec1_B,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
34334a97f8cSStefano Zampini     ierr = VecNorm(pcis->vec1_global,NORM_INFINITY,&error);CHKERRQ(ierr);
34434a97f8cSStefano Zampini     ierr = PetscViewerASCIIPrintf(viewer,"Error scaling restriction %1.14e\n",error);CHKERRQ(ierr);
34572b8c272SStefano Zampini     ierr = MatDestroy(&B0_B);CHKERRQ(ierr);
34672b8c272SStefano Zampini     ierr = VecDestroy(&B0_Bv);CHKERRQ(ierr);
34772b8c272SStefano Zampini     ierr = VecDestroy(&B0_Bv2);CHKERRQ(ierr);
348674ae819SStefano Zampini   }
349674ae819SStefano Zampini   PetscFunctionReturn(0);
350674ae819SStefano Zampini }
351674ae819SStefano Zampini 
352674ae819SStefano Zampini PetscErrorCode PCBDDCScalingDestroy(PC pc)
353674ae819SStefano Zampini {
354674ae819SStefano Zampini   PC_BDDC*       pcbddc=(PC_BDDC*)pc->data;
355674ae819SStefano Zampini   PetscErrorCode ierr;
356674ae819SStefano Zampini 
357674ae819SStefano Zampini   PetscFunctionBegin;
35834a97f8cSStefano Zampini   if (pcbddc->deluxe_ctx) {
35934a97f8cSStefano Zampini     ierr = PCBDDCScalingDestroy_Deluxe(pc);CHKERRQ(ierr);
360674ae819SStefano Zampini   }
361674ae819SStefano Zampini   ierr = VecDestroy(&pcbddc->work_scaling);CHKERRQ(ierr);
362674ae819SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCScalingRestriction_C",NULL);CHKERRQ(ierr);
363674ae819SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCScalingExtension_C",NULL);CHKERRQ(ierr);
364674ae819SStefano Zampini   PetscFunctionReturn(0);
365674ae819SStefano Zampini }
366674ae819SStefano Zampini 
36734a97f8cSStefano Zampini static PetscErrorCode PCBDDCScalingCreate_Deluxe(PC pc)
36834a97f8cSStefano Zampini {
36934a97f8cSStefano Zampini   PC_BDDC*            pcbddc=(PC_BDDC*)pc->data;
37034a97f8cSStefano Zampini   PCBDDCDeluxeScaling deluxe_ctx;
37134a97f8cSStefano Zampini   PetscErrorCode      ierr;
37234a97f8cSStefano Zampini 
37334a97f8cSStefano Zampini   PetscFunctionBegin;
37434a97f8cSStefano Zampini   ierr = PetscNew(&deluxe_ctx);CHKERRQ(ierr);
37534a97f8cSStefano Zampini   pcbddc->deluxe_ctx = deluxe_ctx;
37634a97f8cSStefano Zampini   PetscFunctionReturn(0);
37734a97f8cSStefano Zampini }
37834a97f8cSStefano Zampini 
37934a97f8cSStefano Zampini static PetscErrorCode PCBDDCScalingDestroy_Deluxe(PC pc)
38034a97f8cSStefano Zampini {
38134a97f8cSStefano Zampini   PC_BDDC*            pcbddc=(PC_BDDC*)pc->data;
38234a97f8cSStefano Zampini   PetscErrorCode      ierr;
38334a97f8cSStefano Zampini 
38434a97f8cSStefano Zampini   PetscFunctionBegin;
38534a97f8cSStefano Zampini   ierr = PCBDDCScalingReset_Deluxe_Solvers(pcbddc->deluxe_ctx);CHKERRQ(ierr);
38634a97f8cSStefano Zampini   ierr = PetscFree(pcbddc->deluxe_ctx);CHKERRQ(ierr);
38734a97f8cSStefano Zampini   PetscFunctionReturn(0);
38834a97f8cSStefano Zampini }
38934a97f8cSStefano Zampini 
39034a97f8cSStefano Zampini static PetscErrorCode PCBDDCScalingReset_Deluxe_Solvers(PCBDDCDeluxeScaling deluxe_ctx)
39134a97f8cSStefano Zampini {
39272b8c272SStefano Zampini   PetscInt       i;
39334a97f8cSStefano Zampini   PetscErrorCode ierr;
39434a97f8cSStefano Zampini 
39534a97f8cSStefano Zampini   PetscFunctionBegin;
39634a97f8cSStefano Zampini   ierr = PetscFree(deluxe_ctx->idx_simple_B);CHKERRQ(ierr);
39734a97f8cSStefano Zampini   deluxe_ctx->n_simple = 0;
39872b8c272SStefano Zampini   for (i=0;i<deluxe_ctx->seq_n;i++) {
39972b8c272SStefano Zampini     ierr = VecScatterDestroy(&deluxe_ctx->seq_scctx[i]);CHKERRQ(ierr);
40072b8c272SStefano Zampini     ierr = VecDestroy(&deluxe_ctx->seq_work1[i]);CHKERRQ(ierr);
40172b8c272SStefano Zampini     ierr = VecDestroy(&deluxe_ctx->seq_work2[i]);CHKERRQ(ierr);
40272b8c272SStefano Zampini     ierr = MatDestroy(&deluxe_ctx->seq_mat[i]);CHKERRQ(ierr);
40372b8c272SStefano Zampini     ierr = MatDestroy(&deluxe_ctx->seq_mat_inv_sum[i]);CHKERRQ(ierr);
40472b8c272SStefano Zampini   }
40572b8c272SStefano 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);
40672b8c272SStefano Zampini   ierr = PetscFree(deluxe_ctx->workspace);CHKERRQ(ierr);
40772b8c272SStefano Zampini   deluxe_ctx->seq_n = 0;
40834a97f8cSStefano Zampini   PetscFunctionReturn(0);
40934a97f8cSStefano Zampini }
41034a97f8cSStefano Zampini 
41134a97f8cSStefano Zampini static PetscErrorCode PCBDDCScalingSetUp_Deluxe(PC pc)
41234a97f8cSStefano Zampini {
41334a97f8cSStefano Zampini   PC_IS               *pcis=(PC_IS*)pc->data;
41434a97f8cSStefano Zampini   PC_BDDC             *pcbddc=(PC_BDDC*)pc->data;
41534a97f8cSStefano Zampini   PCBDDCDeluxeScaling deluxe_ctx=pcbddc->deluxe_ctx;
416b96c3477SStefano Zampini   PCBDDCSubSchurs     sub_schurs=pcbddc->sub_schurs;
41734a97f8cSStefano Zampini   PetscErrorCode      ierr;
41834a97f8cSStefano Zampini 
41934a97f8cSStefano Zampini   PetscFunctionBegin;
42072b8c272SStefano Zampini   /* reset data structures if the topology has changed */
42172b8c272SStefano Zampini   if (pcbddc->recompute_topography) {
42234a97f8cSStefano Zampini     ierr = PCBDDCScalingReset_Deluxe_Solvers(deluxe_ctx);CHKERRQ(ierr);
42372b8c272SStefano Zampini   }
42434a97f8cSStefano Zampini 
425b1b3d7a2SStefano Zampini   /* Compute data structures to solve sequential problems */
4265a95e1ceSStefano Zampini   ierr = PCBDDCScalingSetUp_Deluxe_Private(pc);CHKERRQ(ierr);
4275db18549SStefano Zampini 
428b1b3d7a2SStefano Zampini   /* diagonal scaling on interface dofs not contained in cc */
429d62866d3SStefano Zampini   if (sub_schurs->is_vertices || sub_schurs->is_dir) {
4301b968477SStefano Zampini     PetscInt n_com,n_dir;
4311b968477SStefano Zampini     n_com = 0;
432d62866d3SStefano Zampini     if (sub_schurs->is_vertices) {
433d62866d3SStefano Zampini       ierr = ISGetLocalSize(sub_schurs->is_vertices,&n_com);CHKERRQ(ierr);
4341b968477SStefano Zampini     }
4351b968477SStefano Zampini     n_dir = 0;
436d62866d3SStefano Zampini     if (sub_schurs->is_dir) {
437d62866d3SStefano Zampini       ierr = ISGetLocalSize(sub_schurs->is_dir,&n_dir);CHKERRQ(ierr);
4381b968477SStefano Zampini     }
43972b8c272SStefano Zampini     if (!deluxe_ctx->n_simple) {
4401b968477SStefano Zampini       deluxe_ctx->n_simple = n_dir + n_com;
4411b968477SStefano Zampini       ierr = PetscMalloc1(deluxe_ctx->n_simple,&deluxe_ctx->idx_simple_B);CHKERRQ(ierr);
442d62866d3SStefano Zampini       if (sub_schurs->is_vertices) {
4439bb4a8caSStefano Zampini         PetscInt       nmap;
4449bb4a8caSStefano Zampini         const PetscInt *idxs;
4451b968477SStefano Zampini 
446d62866d3SStefano Zampini         ierr = ISGetIndices(sub_schurs->is_vertices,&idxs);CHKERRQ(ierr);
4471b968477SStefano Zampini         ierr = ISGlobalToLocalMappingApply(pcis->BtoNmap,IS_GTOLM_DROP,n_com,idxs,&nmap,deluxe_ctx->idx_simple_B);CHKERRQ(ierr);
448*2c71b3e2SJacob Faibussowitsch         PetscCheckFalse(nmap != n_com,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Error when mapping simply scaled dofs (is_vertices)! %D != %D",nmap,n_com);
449d62866d3SStefano Zampini         ierr = ISRestoreIndices(sub_schurs->is_vertices,&idxs);CHKERRQ(ierr);
4501b968477SStefano Zampini       }
451d62866d3SStefano Zampini       if (sub_schurs->is_dir) {
4521b968477SStefano Zampini         PetscInt       nmap;
4531b968477SStefano Zampini         const PetscInt *idxs;
4541b968477SStefano Zampini 
455d62866d3SStefano Zampini         ierr = ISGetIndices(sub_schurs->is_dir,&idxs);CHKERRQ(ierr);
4561b968477SStefano Zampini         ierr = ISGlobalToLocalMappingApply(pcis->BtoNmap,IS_GTOLM_DROP,n_dir,idxs,&nmap,deluxe_ctx->idx_simple_B+n_com);CHKERRQ(ierr);
457*2c71b3e2SJacob Faibussowitsch         PetscCheckFalse(nmap != n_dir,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Error when mapping simply scaled dofs (sub_schurs->is_dir)! %D != %D",nmap,n_dir);
458d62866d3SStefano Zampini         ierr = ISRestoreIndices(sub_schurs->is_dir,&idxs);CHKERRQ(ierr);
4591b968477SStefano Zampini       }
4601b968477SStefano Zampini       ierr = PetscSortInt(deluxe_ctx->n_simple,deluxe_ctx->idx_simple_B);CHKERRQ(ierr);
4619bb4a8caSStefano Zampini     } else {
462*2c71b3e2SJacob Faibussowitsch       PetscCheckFalse(deluxe_ctx->n_simple != n_dir + n_com,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);
46372b8c272SStefano Zampini     }
46472b8c272SStefano Zampini   } else {
465b1b3d7a2SStefano Zampini     deluxe_ctx->n_simple = 0;
4660a545947SLisandro Dalcin     deluxe_ctx->idx_simple_B = NULL;
467b1b3d7a2SStefano Zampini   }
46834a97f8cSStefano Zampini   PetscFunctionReturn(0);
46934a97f8cSStefano Zampini }
47034a97f8cSStefano Zampini 
4715a95e1ceSStefano Zampini static PetscErrorCode PCBDDCScalingSetUp_Deluxe_Private(PC pc)
4725db18549SStefano Zampini {
4735db18549SStefano Zampini   PC_BDDC                *pcbddc=(PC_BDDC*)pc->data;
4745db18549SStefano Zampini   PCBDDCDeluxeScaling    deluxe_ctx=pcbddc->deluxe_ctx;
475b96c3477SStefano Zampini   PCBDDCSubSchurs        sub_schurs = pcbddc->sub_schurs;
47672b8c272SStefano Zampini   PetscScalar            *matdata,*matdata2;
47772b8c272SStefano Zampini   PetscInt               i,max_subset_size,cum,cum2;
47872b8c272SStefano Zampini   const PetscInt         *idxs;
47972b8c272SStefano Zampini   PetscBool              newsetup = PETSC_FALSE;
4805db18549SStefano Zampini   PetscErrorCode         ierr;
4815db18549SStefano Zampini 
4825db18549SStefano Zampini   PetscFunctionBegin;
483*2c71b3e2SJacob Faibussowitsch   PetscCheckFalse(!sub_schurs,PetscObjectComm((PetscObject)pc),PETSC_ERR_PLIB,"Missing PCBDDCSubSchurs");
48486bfa4cfSStefano Zampini   if (!sub_schurs->n_subs) PetscFunctionReturn(0);
4859221af80SStefano Zampini 
48672b8c272SStefano Zampini   /* Allocate arrays for subproblems */
48772b8c272SStefano Zampini   if (!deluxe_ctx->seq_n) {
48872b8c272SStefano Zampini     deluxe_ctx->seq_n = sub_schurs->n_subs;
48972b8c272SStefano 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);
49072b8c272SStefano Zampini     newsetup = PETSC_TRUE;
491*2c71b3e2SJacob Faibussowitsch   } else PetscCheckFalse(deluxe_ctx->seq_n != sub_schurs->n_subs,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);
4926080607fSStefano Zampini 
4932f41f7d1SStefano Zampini   /* the change of basis is just a reference to sub_schurs->change (if any) */
4942f41f7d1SStefano Zampini   deluxe_ctx->change         = sub_schurs->change;
49572b8c272SStefano Zampini   deluxe_ctx->change_with_qr = sub_schurs->change_with_qr;
4965db18549SStefano Zampini 
49772b8c272SStefano Zampini   /* Create objects for deluxe */
49872b8c272SStefano Zampini   max_subset_size = 0;
49972b8c272SStefano Zampini   for (i=0;i<sub_schurs->n_subs;i++) {
50072b8c272SStefano Zampini     PetscInt subset_size;
50172b8c272SStefano Zampini     ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr);
50272b8c272SStefano Zampini     max_subset_size = PetscMax(subset_size,max_subset_size);
50372b8c272SStefano Zampini   }
50472b8c272SStefano Zampini   if (newsetup) {
50572b8c272SStefano Zampini     ierr = PetscMalloc1(2*max_subset_size,&deluxe_ctx->workspace);CHKERRQ(ierr);
50672b8c272SStefano Zampini   }
50772b8c272SStefano Zampini   cum = cum2 = 0;
50872b8c272SStefano Zampini   ierr = ISGetIndices(sub_schurs->is_Ej_all,&idxs);CHKERRQ(ierr);
50972b8c272SStefano Zampini   ierr = MatSeqAIJGetArray(sub_schurs->S_Ej_all,&matdata);CHKERRQ(ierr);
51072b8c272SStefano Zampini   ierr = MatSeqAIJGetArray(sub_schurs->sum_S_Ej_all,&matdata2);CHKERRQ(ierr);
51172b8c272SStefano Zampini   for (i=0;i<deluxe_ctx->seq_n;i++) {
51272b8c272SStefano Zampini     PetscInt     subset_size;
513925dfd53SStefano Zampini 
51472b8c272SStefano Zampini     ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr);
51572b8c272SStefano Zampini     if (newsetup) {
51672b8c272SStefano Zampini       IS  sub;
51772b8c272SStefano Zampini       /* work vectors */
51872b8c272SStefano Zampini       ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,1,subset_size,deluxe_ctx->workspace,&deluxe_ctx->seq_work1[i]);CHKERRQ(ierr);
51972b8c272SStefano Zampini       ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,1,subset_size,deluxe_ctx->workspace+subset_size,&deluxe_ctx->seq_work2[i]);CHKERRQ(ierr);
52072b8c272SStefano Zampini 
52172b8c272SStefano Zampini       /* scatters */
52272b8c272SStefano Zampini       ierr = ISCreateGeneral(PETSC_COMM_SELF,subset_size,idxs+cum,PETSC_COPY_VALUES,&sub);CHKERRQ(ierr);
5239448b7f1SJunchao Zhang       ierr = VecScatterCreate(pcbddc->work_scaling,sub,deluxe_ctx->seq_work1[i],NULL,&deluxe_ctx->seq_scctx[i]);CHKERRQ(ierr);
52472b8c272SStefano Zampini       ierr = ISDestroy(&sub);CHKERRQ(ierr);
525839e9adbSstefano_zampini     }
52672b8c272SStefano Zampini 
52772b8c272SStefano Zampini     /* S_E_j */
528839e9adbSstefano_zampini     ierr = MatDestroy(&deluxe_ctx->seq_mat[i]);CHKERRQ(ierr);
52972b8c272SStefano Zampini     ierr = MatCreateSeqDense(PETSC_COMM_SELF,subset_size,subset_size,matdata+cum2,&deluxe_ctx->seq_mat[i]);CHKERRQ(ierr);
530839e9adbSstefano_zampini 
53172b8c272SStefano Zampini     /* \sum_k S^k_E_j */
53272b8c272SStefano Zampini     ierr = MatDestroy(&deluxe_ctx->seq_mat_inv_sum[i]);CHKERRQ(ierr);
53372b8c272SStefano Zampini     ierr = MatCreateSeqDense(PETSC_COMM_SELF,subset_size,subset_size,matdata2+cum2,&deluxe_ctx->seq_mat_inv_sum[i]);CHKERRQ(ierr);
5343829e7e6SStefano Zampini     ierr = MatSetOption(deluxe_ctx->seq_mat_inv_sum[i],MAT_SPD,sub_schurs->is_posdef);CHKERRQ(ierr);
5353829e7e6SStefano Zampini     ierr = MatSetOption(deluxe_ctx->seq_mat_inv_sum[i],MAT_HERMITIAN,sub_schurs->is_hermitian);CHKERRQ(ierr);
5363829e7e6SStefano Zampini     if (sub_schurs->is_hermitian) {
53772b8c272SStefano Zampini       ierr = MatCholeskyFactor(deluxe_ctx->seq_mat_inv_sum[i],NULL,NULL);CHKERRQ(ierr);
538d62866d3SStefano Zampini     } else {
53972b8c272SStefano Zampini       ierr = MatLUFactor(deluxe_ctx->seq_mat_inv_sum[i],NULL,NULL,NULL);CHKERRQ(ierr);
540d62866d3SStefano Zampini     }
541839e9adbSstefano_zampini     if (pcbddc->deluxe_singlemat) {
542839e9adbSstefano_zampini       Mat X,Y;
5433829e7e6SStefano Zampini       if (!sub_schurs->is_hermitian) {
544839e9adbSstefano_zampini         ierr = MatTranspose(deluxe_ctx->seq_mat[i],MAT_INITIAL_MATRIX,&X);CHKERRQ(ierr);
545839e9adbSstefano_zampini       } else {
546839e9adbSstefano_zampini         ierr = PetscObjectReference((PetscObject)deluxe_ctx->seq_mat[i]);CHKERRQ(ierr);
547839e9adbSstefano_zampini         X    = deluxe_ctx->seq_mat[i];
548839e9adbSstefano_zampini       }
549839e9adbSstefano_zampini       ierr = MatDuplicate(X,MAT_DO_NOT_COPY_VALUES,&Y);CHKERRQ(ierr);
5503829e7e6SStefano Zampini       if (!sub_schurs->is_hermitian) {
551839e9adbSstefano_zampini         ierr = PCBDDCMatTransposeMatSolve_SeqDense(deluxe_ctx->seq_mat_inv_sum[i],X,Y);CHKERRQ(ierr);
552839e9adbSstefano_zampini       } else {
553839e9adbSstefano_zampini         ierr = MatMatSolve(deluxe_ctx->seq_mat_inv_sum[i],X,Y);CHKERRQ(ierr);
554839e9adbSstefano_zampini       }
555729c86d3Sstefano_zampini 
556839e9adbSstefano_zampini       ierr = MatDestroy(&deluxe_ctx->seq_mat_inv_sum[i]);CHKERRQ(ierr);
557839e9adbSstefano_zampini       ierr = MatDestroy(&deluxe_ctx->seq_mat[i]);CHKERRQ(ierr);
558839e9adbSstefano_zampini       ierr = MatDestroy(&X);CHKERRQ(ierr);
5592f41f7d1SStefano Zampini       if (deluxe_ctx->change) {
5602f41f7d1SStefano Zampini         Mat C,CY;
561*2c71b3e2SJacob Faibussowitsch         PetscCheckFalse(!deluxe_ctx->change_with_qr,PETSC_COMM_SELF,PETSC_ERR_SUP,"Only QR based change of basis");
5622f41f7d1SStefano Zampini         ierr = KSPGetOperators(deluxe_ctx->change[i],&C,NULL);CHKERRQ(ierr);
5632f41f7d1SStefano Zampini         ierr = MatMatMult(C,Y,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&CY);CHKERRQ(ierr);
5642f41f7d1SStefano Zampini         ierr = MatMatTransposeMult(CY,C,MAT_REUSE_MATRIX,PETSC_DEFAULT,&Y);CHKERRQ(ierr);
5652f41f7d1SStefano Zampini         ierr = MatDestroy(&CY);CHKERRQ(ierr);
5664417c5e8SHong Zhang         ierr = MatProductClear(Y);CHKERRQ(ierr); /* clear internal matproduct structure of Y since CY is destroyed */
5672f41f7d1SStefano Zampini       }
56803dfb2d7SStefano Zampini       ierr = MatTranspose(Y,MAT_INPLACE_MATRIX,&Y);CHKERRQ(ierr);
569839e9adbSstefano_zampini       deluxe_ctx->seq_mat[i] = Y;
57072b8c272SStefano Zampini     }
57172b8c272SStefano Zampini     cum += subset_size;
57272b8c272SStefano Zampini     cum2 += subset_size*subset_size;
57372b8c272SStefano Zampini   }
57472b8c272SStefano Zampini   ierr = ISRestoreIndices(sub_schurs->is_Ej_all,&idxs);CHKERRQ(ierr);
57572b8c272SStefano Zampini   ierr = MatSeqAIJRestoreArray(sub_schurs->S_Ej_all,&matdata);CHKERRQ(ierr);
57672b8c272SStefano Zampini   ierr = MatSeqAIJRestoreArray(sub_schurs->sum_S_Ej_all,&matdata2);CHKERRQ(ierr);
5772f41f7d1SStefano Zampini   if (pcbddc->deluxe_singlemat) {
5782f41f7d1SStefano Zampini     deluxe_ctx->change         = NULL;
5792f41f7d1SStefano Zampini     deluxe_ctx->change_with_qr = PETSC_FALSE;
5802f41f7d1SStefano Zampini   }
5815db18549SStefano Zampini 
58272b8c272SStefano Zampini   if (deluxe_ctx->change && !deluxe_ctx->change_with_qr) {
58372b8c272SStefano Zampini     for (i=0;i<deluxe_ctx->seq_n;i++) {
58472b8c272SStefano Zampini       if (newsetup) {
58572b8c272SStefano Zampini         PC pc;
58672b8c272SStefano Zampini 
58772b8c272SStefano Zampini         ierr = KSPGetPC(deluxe_ctx->change[i],&pc);CHKERRQ(ierr);
58872b8c272SStefano Zampini         ierr = PCSetType(pc,PCLU);CHKERRQ(ierr);
58972b8c272SStefano Zampini         ierr = KSPSetFromOptions(deluxe_ctx->change[i]);CHKERRQ(ierr);
59072b8c272SStefano Zampini       }
59172b8c272SStefano Zampini       ierr = KSPSetUp(deluxe_ctx->change[i]);CHKERRQ(ierr);
59272b8c272SStefano Zampini     }
59316e386b8SStefano Zampini   }
5945db18549SStefano Zampini   PetscFunctionReturn(0);
5955db18549SStefano Zampini }
596