xref: /petsc/src/ksp/pc/impls/bddc/bddcscalingbasic.c (revision af25d912c737f21fca88746c26cd8b5ca8321f21)
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>
4674ae819SStefano Zampini 
534a97f8cSStefano Zampini /* prototypes for deluxe functions */
634a97f8cSStefano Zampini static PetscErrorCode PCBDDCScalingCreate_Deluxe(PC);
734a97f8cSStefano Zampini static PetscErrorCode PCBDDCScalingDestroy_Deluxe(PC);
834a97f8cSStefano Zampini static PetscErrorCode PCBDDCScalingSetUp_Deluxe(PC);
95a95e1ceSStefano Zampini static PetscErrorCode PCBDDCScalingSetUp_Deluxe_Private(PC);
1034a97f8cSStefano Zampini static PetscErrorCode PCBDDCScalingReset_Deluxe_Solvers(PCBDDCDeluxeScaling);
11674ae819SStefano Zampini 
12674ae819SStefano Zampini #undef __FUNCT__
13674ae819SStefano Zampini #define __FUNCT__ "PCBDDCScalingExtension_Basic"
14674ae819SStefano Zampini static PetscErrorCode PCBDDCScalingExtension_Basic(PC pc, Vec local_interface_vector, Vec global_vector)
15674ae819SStefano Zampini {
16674ae819SStefano Zampini   PC_IS*         pcis = (PC_IS*)pc->data;
17674ae819SStefano Zampini   PC_BDDC*       pcbddc = (PC_BDDC*)pc->data;
18674ae819SStefano Zampini   PetscErrorCode ierr;
19674ae819SStefano Zampini 
20674ae819SStefano Zampini   PetscFunctionBegin;
21674ae819SStefano Zampini   /* Apply partition of unity */
22674ae819SStefano Zampini   ierr = VecPointwiseMult(pcbddc->work_scaling,pcis->D,local_interface_vector);CHKERRQ(ierr);
23674ae819SStefano Zampini   ierr = VecSet(global_vector,0.0);CHKERRQ(ierr);
24674ae819SStefano Zampini   ierr = VecScatterBegin(pcis->global_to_B,pcbddc->work_scaling,global_vector,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
25674ae819SStefano Zampini   ierr = VecScatterEnd(pcis->global_to_B,pcbddc->work_scaling,global_vector,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
26674ae819SStefano Zampini   PetscFunctionReturn(0);
27674ae819SStefano Zampini }
28674ae819SStefano Zampini 
29674ae819SStefano Zampini #undef __FUNCT__
30674ae819SStefano Zampini #define __FUNCT__ "PCBDDCScalingExtension_Deluxe"
31674ae819SStefano Zampini static PetscErrorCode PCBDDCScalingExtension_Deluxe(PC pc, Vec x, Vec y)
32674ae819SStefano Zampini {
33674ae819SStefano Zampini   PC_IS*              pcis=(PC_IS*)pc->data;
34674ae819SStefano Zampini   PC_BDDC*            pcbddc=(PC_BDDC*)pc->data;
35674ae819SStefano Zampini   PCBDDCDeluxeScaling deluxe_ctx = pcbddc->deluxe_ctx;
36674ae819SStefano Zampini   PetscErrorCode      ierr;
37674ae819SStefano Zampini 
38674ae819SStefano Zampini   PetscFunctionBegin;
395a95e1ceSStefano Zampini   ierr = VecSet(pcbddc->work_scaling,0.0);CHKERRQ(ierr);
405a95e1ceSStefano Zampini   ierr = VecSet(y,0.0);CHKERRQ(ierr);
415a95e1ceSStefano Zampini   if (deluxe_ctx->n_simple) { /* scale deluxe vertices using diagonal scaling */
425a95e1ceSStefano Zampini     PetscInt          i;
432b095fd8SStefano Zampini     const PetscScalar *array_x,*array_D;
442b095fd8SStefano Zampini     PetscScalar       *array;
452b095fd8SStefano Zampini     ierr = VecGetArrayRead(x,&array_x);CHKERRQ(ierr);
462b095fd8SStefano Zampini     ierr = VecGetArrayRead(pcis->D,&array_D);CHKERRQ(ierr);
47674ae819SStefano Zampini     ierr = VecGetArray(pcbddc->work_scaling,&array);CHKERRQ(ierr);
48674ae819SStefano Zampini     for (i=0;i<deluxe_ctx->n_simple;i++) {
49674ae819SStefano Zampini       array[deluxe_ctx->idx_simple_B[i]] = array_x[deluxe_ctx->idx_simple_B[i]]*array_D[deluxe_ctx->idx_simple_B[i]];
50674ae819SStefano Zampini     }
51674ae819SStefano Zampini     ierr = VecRestoreArray(pcbddc->work_scaling,&array);CHKERRQ(ierr);
522b095fd8SStefano Zampini     ierr = VecRestoreArrayRead(pcis->D,&array_D);CHKERRQ(ierr);
532b095fd8SStefano Zampini     ierr = VecRestoreArrayRead(x,&array_x);CHKERRQ(ierr);
5434a97f8cSStefano Zampini   }
55ac632422SStefano Zampini   /* sequential part : all problems and Schur applications collapsed into a single matrix vector multiplication or a matvec and a solve */
565a95e1ceSStefano Zampini   if (deluxe_ctx->seq_mat) {
5772b8c272SStefano Zampini     PetscInt i;
5872b8c272SStefano Zampini     for (i=0;i<deluxe_ctx->seq_n;i++) {
5916e386b8SStefano Zampini       if (deluxe_ctx->change) {
6072b8c272SStefano Zampini         ierr = VecScatterBegin(deluxe_ctx->seq_scctx[i],x,deluxe_ctx->seq_work2[i],INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
6172b8c272SStefano Zampini         ierr = VecScatterEnd(deluxe_ctx->seq_scctx[i],x,deluxe_ctx->seq_work2[i],INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
6272b8c272SStefano Zampini         if (deluxe_ctx->change_with_qr) {
6372b8c272SStefano Zampini           Mat change;
6472b8c272SStefano Zampini 
6572b8c272SStefano Zampini           ierr = KSPGetOperators(deluxe_ctx->change[i],&change,NULL);CHKERRQ(ierr);
6672b8c272SStefano Zampini           ierr = MatMultTranspose(change,deluxe_ctx->seq_work2[i],deluxe_ctx->seq_work1[i]);CHKERRQ(ierr);
6772b8c272SStefano Zampini         } else {
6872b8c272SStefano Zampini           ierr = KSPSolve(deluxe_ctx->change[i],deluxe_ctx->seq_work2[i],deluxe_ctx->seq_work1[i]);CHKERRQ(ierr);
6916e386b8SStefano Zampini         }
7072b8c272SStefano Zampini       } else {
7172b8c272SStefano Zampini         ierr = VecScatterBegin(deluxe_ctx->seq_scctx[i],x,deluxe_ctx->seq_work1[i],INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
7272b8c272SStefano Zampini         ierr = VecScatterEnd(deluxe_ctx->seq_scctx[i],x,deluxe_ctx->seq_work1[i],INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
7372b8c272SStefano Zampini       }
7472b8c272SStefano Zampini       ierr = MatMultTranspose(deluxe_ctx->seq_mat[i],deluxe_ctx->seq_work1[i],deluxe_ctx->seq_work2[i]);CHKERRQ(ierr);
7572b8c272SStefano Zampini       if (deluxe_ctx->seq_mat_inv_sum[i]) {
7672b8c272SStefano Zampini         PetscScalar *x;
7772b8c272SStefano Zampini 
7872b8c272SStefano Zampini         ierr = VecGetArray(deluxe_ctx->seq_work2[i],&x);CHKERRQ(ierr);
7972b8c272SStefano Zampini         ierr = VecPlaceArray(deluxe_ctx->seq_work1[i],x);CHKERRQ(ierr);
8072b8c272SStefano Zampini         ierr = VecRestoreArray(deluxe_ctx->seq_work2[i],&x);CHKERRQ(ierr);
8172b8c272SStefano Zampini         ierr = MatSolveTranspose(deluxe_ctx->seq_mat_inv_sum[i],deluxe_ctx->seq_work1[i],deluxe_ctx->seq_work2[i]);CHKERRQ(ierr);
8272b8c272SStefano Zampini         ierr = VecResetArray(deluxe_ctx->seq_work1[i]);CHKERRQ(ierr);
83ac632422SStefano Zampini       }
8416e386b8SStefano Zampini       if (deluxe_ctx->change) {
8516e386b8SStefano Zampini         Mat change;
8616e386b8SStefano Zampini 
8772b8c272SStefano Zampini         ierr = KSPGetOperators(deluxe_ctx->change[i],&change,NULL);CHKERRQ(ierr);
8872b8c272SStefano Zampini         ierr = MatMult(change,deluxe_ctx->seq_work2[i],deluxe_ctx->seq_work1[i]);CHKERRQ(ierr);
8972b8c272SStefano Zampini         ierr = VecScatterBegin(deluxe_ctx->seq_scctx[i],deluxe_ctx->seq_work1[i],pcbddc->work_scaling,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
9072b8c272SStefano Zampini         ierr = VecScatterEnd(deluxe_ctx->seq_scctx[i],deluxe_ctx->seq_work1[i],pcbddc->work_scaling,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
9116e386b8SStefano Zampini       } else {
9272b8c272SStefano Zampini         ierr = VecScatterBegin(deluxe_ctx->seq_scctx[i],deluxe_ctx->seq_work2[i],pcbddc->work_scaling,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
9372b8c272SStefano Zampini         ierr = VecScatterEnd(deluxe_ctx->seq_scctx[i],deluxe_ctx->seq_work2[i],pcbddc->work_scaling,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
9472b8c272SStefano Zampini       }
95674ae819SStefano Zampini     }
9616e386b8SStefano Zampini   }
97674ae819SStefano Zampini   /* put local boundary part in global vector */
98674ae819SStefano Zampini   ierr = VecScatterBegin(pcis->global_to_B,pcbddc->work_scaling,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
99674ae819SStefano Zampini   ierr = VecScatterEnd(pcis->global_to_B,pcbddc->work_scaling,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
100674ae819SStefano Zampini   PetscFunctionReturn(0);
101674ae819SStefano Zampini }
102674ae819SStefano Zampini 
103674ae819SStefano Zampini #undef __FUNCT__
104674ae819SStefano Zampini #define __FUNCT__ "PCBDDCScalingExtension"
105674ae819SStefano Zampini PetscErrorCode PCBDDCScalingExtension(PC pc, Vec local_interface_vector, Vec global_vector)
106674ae819SStefano Zampini {
107674ae819SStefano Zampini   PC_BDDC *pcbddc=(PC_BDDC*)pc->data;
108674ae819SStefano Zampini   PetscErrorCode ierr;
109674ae819SStefano Zampini 
110674ae819SStefano Zampini   PetscFunctionBegin;
111674ae819SStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
112674ae819SStefano Zampini   PetscValidHeaderSpecific(local_interface_vector,VEC_CLASSID,2);
113674ae819SStefano Zampini   PetscValidHeaderSpecific(global_vector,VEC_CLASSID,3);
1146c4ed002SBarry Smith   if (local_interface_vector == pcbddc->work_scaling) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Local vector cannot be pcbddc->work_scaling!\n");
115674ae819SStefano Zampini   ierr = PetscTryMethod(pc,"PCBDDCScalingExtension_C",(PC,Vec,Vec),(pc,local_interface_vector,global_vector));CHKERRQ(ierr);
116674ae819SStefano Zampini   PetscFunctionReturn(0);
117674ae819SStefano Zampini }
118674ae819SStefano Zampini 
119674ae819SStefano Zampini #undef __FUNCT__
120674ae819SStefano Zampini #define __FUNCT__ "PCBDDCScalingRestriction_Basic"
121674ae819SStefano Zampini static PetscErrorCode PCBDDCScalingRestriction_Basic(PC pc, Vec global_vector, Vec local_interface_vector)
122674ae819SStefano Zampini {
123674ae819SStefano Zampini   PetscErrorCode ierr;
124674ae819SStefano Zampini   PC_IS* pcis = (PC_IS*)pc->data;
125674ae819SStefano Zampini 
126674ae819SStefano Zampini   PetscFunctionBegin;
127674ae819SStefano Zampini   ierr = VecScatterBegin(pcis->global_to_B,global_vector,local_interface_vector,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
128674ae819SStefano Zampini   ierr = VecScatterEnd(pcis->global_to_B,global_vector,local_interface_vector,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
129674ae819SStefano Zampini   /* Apply partition of unity */
130674ae819SStefano Zampini   ierr = VecPointwiseMult(local_interface_vector,pcis->D,local_interface_vector);CHKERRQ(ierr);
131674ae819SStefano Zampini   PetscFunctionReturn(0);
132674ae819SStefano Zampini }
133674ae819SStefano Zampini 
134674ae819SStefano Zampini #undef __FUNCT__
135674ae819SStefano Zampini #define __FUNCT__ "PCBDDCScalingRestriction_Deluxe"
136674ae819SStefano Zampini static PetscErrorCode PCBDDCScalingRestriction_Deluxe(PC pc, Vec x, Vec y)
137674ae819SStefano Zampini {
138674ae819SStefano Zampini   PC_IS*              pcis=(PC_IS*)pc->data;
139674ae819SStefano Zampini   PC_BDDC*            pcbddc=(PC_BDDC*)pc->data;
140674ae819SStefano Zampini   PCBDDCDeluxeScaling deluxe_ctx = pcbddc->deluxe_ctx;
141674ae819SStefano Zampini   PetscErrorCode      ierr;
142674ae819SStefano Zampini 
143674ae819SStefano Zampini   PetscFunctionBegin;
144674ae819SStefano Zampini   /* get local boundary part of global vector */
145674ae819SStefano Zampini   ierr = VecScatterBegin(pcis->global_to_B,x,y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
146674ae819SStefano Zampini   ierr = VecScatterEnd(pcis->global_to_B,x,y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1475a95e1ceSStefano Zampini   if (deluxe_ctx->n_simple) { /* scale deluxe vertices using diagonal scaling */
1485a95e1ceSStefano Zampini     PetscInt          i;
1492b095fd8SStefano Zampini     PetscScalar       *array_y;
1502b095fd8SStefano Zampini     const PetscScalar *array_D;
151674ae819SStefano Zampini     ierr = VecGetArray(y,&array_y);CHKERRQ(ierr);
1522b095fd8SStefano Zampini     ierr = VecGetArrayRead(pcis->D,&array_D);CHKERRQ(ierr);
153674ae819SStefano Zampini     for (i=0;i<deluxe_ctx->n_simple;i++) {
154674ae819SStefano Zampini       array_y[deluxe_ctx->idx_simple_B[i]] *= array_D[deluxe_ctx->idx_simple_B[i]];
155674ae819SStefano Zampini     }
1562b095fd8SStefano Zampini     ierr = VecRestoreArrayRead(pcis->D,&array_D);CHKERRQ(ierr);
157674ae819SStefano Zampini     ierr = VecRestoreArray(y,&array_y);CHKERRQ(ierr);
15834a97f8cSStefano Zampini   }
15934a97f8cSStefano Zampini   /* sequential part : all problems and Schur applications collapsed into a single matrix vector multiplication and ksp solution */
1605a95e1ceSStefano Zampini   if (deluxe_ctx->seq_mat) {
16172b8c272SStefano Zampini     PetscInt i;
16272b8c272SStefano Zampini     for (i=0;i<deluxe_ctx->seq_n;i++) {
16316e386b8SStefano Zampini       if (deluxe_ctx->change) {
16416e386b8SStefano Zampini         Mat change;
16516e386b8SStefano Zampini 
16672b8c272SStefano Zampini         ierr = VecScatterBegin(deluxe_ctx->seq_scctx[i],y,deluxe_ctx->seq_work2[i],INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
16772b8c272SStefano Zampini         ierr = VecScatterEnd(deluxe_ctx->seq_scctx[i],y,deluxe_ctx->seq_work2[i],INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
16872b8c272SStefano Zampini         ierr = KSPGetOperators(deluxe_ctx->change[i],&change,NULL);CHKERRQ(ierr);
16972b8c272SStefano Zampini         ierr = MatMultTranspose(change,deluxe_ctx->seq_work2[i],deluxe_ctx->seq_work1[i]);CHKERRQ(ierr);
17016e386b8SStefano Zampini       } else {
17172b8c272SStefano Zampini         ierr = VecScatterBegin(deluxe_ctx->seq_scctx[i],y,deluxe_ctx->seq_work1[i],INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
17272b8c272SStefano Zampini         ierr = VecScatterEnd(deluxe_ctx->seq_scctx[i],y,deluxe_ctx->seq_work1[i],INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
17316e386b8SStefano Zampini       }
17472b8c272SStefano Zampini       if (deluxe_ctx->seq_mat_inv_sum[i]) {
17572b8c272SStefano Zampini         PetscScalar *x;
17672b8c272SStefano Zampini 
17772b8c272SStefano Zampini         ierr = VecGetArray(deluxe_ctx->seq_work1[i],&x);CHKERRQ(ierr);
17872b8c272SStefano Zampini         ierr = VecPlaceArray(deluxe_ctx->seq_work2[i],x);CHKERRQ(ierr);
17972b8c272SStefano Zampini         ierr = VecRestoreArray(deluxe_ctx->seq_work1[i],&x);CHKERRQ(ierr);
18072b8c272SStefano Zampini         ierr = MatSolve(deluxe_ctx->seq_mat_inv_sum[i],deluxe_ctx->seq_work1[i],deluxe_ctx->seq_work2[i]);CHKERRQ(ierr);
18172b8c272SStefano Zampini         ierr = VecResetArray(deluxe_ctx->seq_work2[i]);CHKERRQ(ierr);
182ac632422SStefano Zampini       }
18372b8c272SStefano Zampini       ierr = MatMult(deluxe_ctx->seq_mat[i],deluxe_ctx->seq_work1[i],deluxe_ctx->seq_work2[i]);CHKERRQ(ierr);
18416e386b8SStefano Zampini       if (deluxe_ctx->change) {
18572b8c272SStefano Zampini         if (deluxe_ctx->change_with_qr) {
18672b8c272SStefano Zampini           Mat change;
18772b8c272SStefano Zampini 
18872b8c272SStefano Zampini           ierr = KSPGetOperators(deluxe_ctx->change[i],&change,NULL);CHKERRQ(ierr);
18972b8c272SStefano Zampini           ierr = MatMult(change,deluxe_ctx->seq_work2[i],deluxe_ctx->seq_work1[i]);CHKERRQ(ierr);
19072b8c272SStefano Zampini         } else {
19172b8c272SStefano Zampini           ierr = KSPSolveTranspose(deluxe_ctx->change[i],deluxe_ctx->seq_work2[i],deluxe_ctx->seq_work1[i]);CHKERRQ(ierr);
19216e386b8SStefano Zampini         }
19372b8c272SStefano Zampini         ierr = VecScatterBegin(deluxe_ctx->seq_scctx[i],deluxe_ctx->seq_work1[i],y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
19472b8c272SStefano Zampini         ierr = VecScatterEnd(deluxe_ctx->seq_scctx[i],deluxe_ctx->seq_work1[i],y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
19572b8c272SStefano Zampini       } else {
19672b8c272SStefano Zampini         ierr = VecScatterBegin(deluxe_ctx->seq_scctx[i],deluxe_ctx->seq_work2[i],y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
19772b8c272SStefano Zampini         ierr = VecScatterEnd(deluxe_ctx->seq_scctx[i],deluxe_ctx->seq_work2[i],y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
19872b8c272SStefano Zampini       }
19972b8c272SStefano Zampini     }
200674ae819SStefano Zampini   }
201674ae819SStefano Zampini   PetscFunctionReturn(0);
202674ae819SStefano Zampini }
203674ae819SStefano Zampini 
204674ae819SStefano Zampini #undef __FUNCT__
205674ae819SStefano Zampini #define __FUNCT__ "PCBDDCScalingRestriction"
206674ae819SStefano Zampini PetscErrorCode PCBDDCScalingRestriction(PC pc, Vec global_vector, Vec local_interface_vector)
207674ae819SStefano Zampini {
208674ae819SStefano Zampini   PC_BDDC        *pcbddc=(PC_BDDC*)pc->data;
209674ae819SStefano Zampini   PetscErrorCode ierr;
210674ae819SStefano Zampini 
211674ae819SStefano Zampini   PetscFunctionBegin;
212674ae819SStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
213674ae819SStefano Zampini   PetscValidHeaderSpecific(global_vector,VEC_CLASSID,2);
214674ae819SStefano Zampini   PetscValidHeaderSpecific(local_interface_vector,VEC_CLASSID,3);
2156c4ed002SBarry Smith   if (local_interface_vector == pcbddc->work_scaling) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Local vector cannot be pcbddc->work_scaling!\n");
216674ae819SStefano Zampini   ierr = PetscTryMethod(pc,"PCBDDCScalingRestriction_C",(PC,Vec,Vec),(pc,global_vector,local_interface_vector));CHKERRQ(ierr);
217674ae819SStefano Zampini   PetscFunctionReturn(0);
218674ae819SStefano Zampini }
219674ae819SStefano Zampini 
220674ae819SStefano Zampini #undef __FUNCT__
221674ae819SStefano Zampini #define __FUNCT__ "PCBDDCScalingSetUp"
222674ae819SStefano Zampini PetscErrorCode PCBDDCScalingSetUp(PC pc)
223674ae819SStefano Zampini {
224674ae819SStefano Zampini   PC_IS*         pcis=(PC_IS*)pc->data;
225674ae819SStefano Zampini   PC_BDDC*       pcbddc=(PC_BDDC*)pc->data;
226674ae819SStefano Zampini   PetscErrorCode ierr;
227674ae819SStefano Zampini 
228674ae819SStefano Zampini   PetscFunctionBegin;
229674ae819SStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
230674ae819SStefano Zampini   /* create work vector for the operator */
23134a97f8cSStefano Zampini   ierr = VecDestroy(&pcbddc->work_scaling);CHKERRQ(ierr);
232674ae819SStefano Zampini   ierr = VecDuplicate(pcis->vec1_B,&pcbddc->work_scaling);CHKERRQ(ierr);
23334a97f8cSStefano Zampini   /* always rebuild pcis->D */
23428d874f6SStefano Zampini   if (pcis->use_stiffness_scaling) {
235674ae819SStefano Zampini     ierr = MatGetDiagonal(pcbddc->local_mat,pcis->vec1_N);CHKERRQ(ierr);
236674ae819SStefano Zampini     ierr = VecScatterBegin(pcis->N_to_B,pcis->vec1_N,pcis->D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
237674ae819SStefano Zampini     ierr = VecScatterEnd(pcis->N_to_B,pcis->vec1_N,pcis->D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
238674ae819SStefano Zampini   }
239674ae819SStefano Zampini   ierr = VecCopy(pcis->D,pcis->vec1_B);CHKERRQ(ierr);
240674ae819SStefano Zampini   ierr = VecSet(pcis->vec1_global,0.0);CHKERRQ(ierr);
241674ae819SStefano Zampini   ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_B,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
242674ae819SStefano Zampini   ierr = VecScatterEnd(pcis->global_to_B,pcis->vec1_B,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
243674ae819SStefano Zampini   ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_global,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
244674ae819SStefano Zampini   ierr = VecScatterEnd(pcis->global_to_B,pcis->vec1_global,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
245674ae819SStefano Zampini   ierr = VecPointwiseDivide(pcis->D,pcis->D,pcis->vec1_B);CHKERRQ(ierr);
246674ae819SStefano Zampini   /* now setup */
247681e7c04SStefano Zampini   if (pcbddc->use_deluxe_scaling) {
24834a97f8cSStefano Zampini     if (!pcbddc->deluxe_ctx) {
24934a97f8cSStefano Zampini       ierr = PCBDDCScalingCreate_Deluxe(pc);CHKERRQ(ierr);
25034a97f8cSStefano Zampini     }
25134a97f8cSStefano Zampini     ierr = PCBDDCScalingSetUp_Deluxe(pc);CHKERRQ(ierr);
252674ae819SStefano Zampini     ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCScalingRestriction_C",PCBDDCScalingRestriction_Deluxe);CHKERRQ(ierr);
253674ae819SStefano Zampini     ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCScalingExtension_C",PCBDDCScalingExtension_Deluxe);CHKERRQ(ierr);
254674ae819SStefano Zampini   } else {
255674ae819SStefano Zampini     ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCScalingRestriction_C",PCBDDCScalingRestriction_Basic);CHKERRQ(ierr);
256674ae819SStefano Zampini     ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCScalingExtension_C",PCBDDCScalingExtension_Basic);CHKERRQ(ierr);
257674ae819SStefano Zampini   }
25834a97f8cSStefano Zampini 
259674ae819SStefano Zampini   /* test */
260674ae819SStefano Zampini   if (pcbddc->dbg_flag) {
26172b8c272SStefano Zampini     Mat         B0_B = NULL;
26272b8c272SStefano Zampini     Vec         B0_Bv = NULL, B0_Bv2 = NULL;
26334a97f8cSStefano Zampini     Vec         vec2_global;
264674ae819SStefano Zampini     PetscViewer viewer = pcbddc->dbg_viewer;
26534a97f8cSStefano Zampini     PetscReal   error;
266674ae819SStefano Zampini 
267674ae819SStefano Zampini     /* extension -> from local to parallel */
26834a97f8cSStefano Zampini     ierr = VecSet(pcis->vec1_global,0.0);CHKERRQ(ierr);
26934a97f8cSStefano Zampini     ierr = VecSetRandom(pcis->vec1_B,NULL);CHKERRQ(ierr);
27034a97f8cSStefano Zampini     ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_B,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
27134a97f8cSStefano Zampini     ierr = VecScatterEnd(pcis->global_to_B,pcis->vec1_B,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
27234a97f8cSStefano Zampini     ierr = VecDuplicate(pcis->vec1_global,&vec2_global);CHKERRQ(ierr);
27334a97f8cSStefano Zampini     ierr = VecCopy(pcis->vec1_global,vec2_global);CHKERRQ(ierr);
274674ae819SStefano Zampini     ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_global,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
275674ae819SStefano Zampini     ierr = VecScatterEnd(pcis->global_to_B,pcis->vec1_global,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
27672b8c272SStefano Zampini     if (pcbddc->benign_n) {
27772b8c272SStefano Zampini       IS is_dummy;
27872b8c272SStefano Zampini 
27972b8c272SStefano Zampini       ierr = ISCreateStride(PETSC_COMM_SELF,pcbddc->benign_n,0,1,&is_dummy);CHKERRQ(ierr);
28072b8c272SStefano Zampini       ierr = MatGetSubMatrix(pcbddc->benign_B0,is_dummy,pcis->is_B_local,MAT_INITIAL_MATRIX,&B0_B);CHKERRQ(ierr);
28172b8c272SStefano Zampini       ierr = ISDestroy(&is_dummy);CHKERRQ(ierr);
28272b8c272SStefano Zampini       ierr = MatCreateVecs(B0_B,NULL,&B0_Bv);CHKERRQ(ierr);
28372b8c272SStefano Zampini       ierr = VecDuplicate(B0_Bv,&B0_Bv2);CHKERRQ(ierr);
28472b8c272SStefano Zampini       ierr = MatMult(B0_B,pcis->vec1_B,B0_Bv);CHKERRQ(ierr);
28572b8c272SStefano Zampini     }
286674ae819SStefano Zampini     ierr = PCBDDCScalingExtension(pc,pcis->vec1_B,pcis->vec1_global);CHKERRQ(ierr);
28772b8c272SStefano Zampini     if (pcbddc->benign_saddle_point) {
28872b8c272SStefano Zampini       PetscReal errorl = 0.;
28972b8c272SStefano Zampini       ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_global,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
29072b8c272SStefano Zampini       ierr = VecScatterEnd(pcis->global_to_B,pcis->vec1_global,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
29172b8c272SStefano Zampini       if (pcbddc->benign_n) {
29272b8c272SStefano Zampini         ierr = MatMult(B0_B,pcis->vec1_B,B0_Bv2);CHKERRQ(ierr);
29372b8c272SStefano Zampini         ierr = VecAXPY(B0_Bv,-1.0,B0_Bv2);CHKERRQ(ierr);
29472b8c272SStefano Zampini         ierr = VecNorm(B0_Bv,NORM_INFINITY,&errorl);CHKERRQ(ierr);
29572b8c272SStefano Zampini       }
29672b8c272SStefano Zampini       ierr = MPI_Allreduce(&errorl,&error,1,MPIU_REAL,MPI_SUM,PetscObjectComm((PetscObject)pc));CHKERRQ(ierr);
29772b8c272SStefano Zampini       ierr = PetscViewerASCIIPrintf(viewer,"Error benign extension %1.14e\n",error);CHKERRQ(ierr);
29872b8c272SStefano Zampini     }
29934a97f8cSStefano Zampini     ierr = VecAXPY(pcis->vec1_global,-1.0,vec2_global);CHKERRQ(ierr);
30034a97f8cSStefano Zampini     ierr = VecNorm(pcis->vec1_global,NORM_INFINITY,&error);CHKERRQ(ierr);
301674ae819SStefano Zampini     ierr = PetscViewerASCIIPrintf(viewer,"Error scaling extension %1.14e\n",error);CHKERRQ(ierr);
30234a97f8cSStefano Zampini     ierr = VecDestroy(&vec2_global);CHKERRQ(ierr);
30334a97f8cSStefano Zampini 
304674ae819SStefano Zampini     /* restriction -> from parallel to local */
305674ae819SStefano Zampini     ierr = VecSet(pcis->vec1_global,0.0);CHKERRQ(ierr);
30634a97f8cSStefano Zampini     ierr = VecSetRandom(pcis->vec1_B,NULL);CHKERRQ(ierr);
307674ae819SStefano Zampini     ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_B,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
308674ae819SStefano Zampini     ierr = VecScatterEnd(pcis->global_to_B,pcis->vec1_B,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
30934a97f8cSStefano Zampini     ierr = PCBDDCScalingRestriction(pc,pcis->vec1_global,pcis->vec1_B);CHKERRQ(ierr);
31034a97f8cSStefano Zampini     ierr = VecScale(pcis->vec1_B,-1.0);CHKERRQ(ierr);
31134a97f8cSStefano Zampini     ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_B,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
31234a97f8cSStefano Zampini     ierr = VecScatterEnd(pcis->global_to_B,pcis->vec1_B,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
31334a97f8cSStefano Zampini     ierr = VecNorm(pcis->vec1_global,NORM_INFINITY,&error);CHKERRQ(ierr);
31434a97f8cSStefano Zampini     ierr = PetscViewerASCIIPrintf(viewer,"Error scaling restriction %1.14e\n",error);CHKERRQ(ierr);
31572b8c272SStefano Zampini     ierr = MatDestroy(&B0_B);CHKERRQ(ierr);
31672b8c272SStefano Zampini     ierr = VecDestroy(&B0_Bv);CHKERRQ(ierr);
31772b8c272SStefano Zampini     ierr = VecDestroy(&B0_Bv2);CHKERRQ(ierr);
318674ae819SStefano Zampini   }
319674ae819SStefano Zampini   PetscFunctionReturn(0);
320674ae819SStefano Zampini }
321674ae819SStefano Zampini 
322674ae819SStefano Zampini #undef __FUNCT__
323674ae819SStefano Zampini #define __FUNCT__ "PCBDDCScalingDestroy"
324674ae819SStefano Zampini PetscErrorCode PCBDDCScalingDestroy(PC pc)
325674ae819SStefano Zampini {
326674ae819SStefano Zampini   PC_BDDC* pcbddc=(PC_BDDC*)pc->data;
327674ae819SStefano Zampini   PetscErrorCode ierr;
328674ae819SStefano Zampini 
329674ae819SStefano Zampini   PetscFunctionBegin;
33034a97f8cSStefano Zampini   if (pcbddc->deluxe_ctx) {
33134a97f8cSStefano Zampini     ierr = PCBDDCScalingDestroy_Deluxe(pc);CHKERRQ(ierr);
332674ae819SStefano Zampini   }
333674ae819SStefano Zampini   ierr = VecDestroy(&pcbddc->work_scaling);CHKERRQ(ierr);
334674ae819SStefano Zampini   /* remove functions */
335674ae819SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCScalingRestriction_C",NULL);CHKERRQ(ierr);
336674ae819SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCScalingExtension_C",NULL);CHKERRQ(ierr);
337674ae819SStefano Zampini   PetscFunctionReturn(0);
338674ae819SStefano Zampini }
339674ae819SStefano Zampini 
34034a97f8cSStefano Zampini #undef __FUNCT__
34134a97f8cSStefano Zampini #define __FUNCT__ "PCBDDCScalingCreate_Deluxe"
34234a97f8cSStefano Zampini static PetscErrorCode PCBDDCScalingCreate_Deluxe(PC pc)
34334a97f8cSStefano Zampini {
34434a97f8cSStefano Zampini   PC_BDDC*            pcbddc=(PC_BDDC*)pc->data;
34534a97f8cSStefano Zampini   PCBDDCDeluxeScaling deluxe_ctx;
34634a97f8cSStefano Zampini   PetscErrorCode      ierr;
34734a97f8cSStefano Zampini 
34834a97f8cSStefano Zampini   PetscFunctionBegin;
34934a97f8cSStefano Zampini   ierr = PetscNew(&deluxe_ctx);CHKERRQ(ierr);
35034a97f8cSStefano Zampini   pcbddc->deluxe_ctx = deluxe_ctx;
35134a97f8cSStefano Zampini   PetscFunctionReturn(0);
35234a97f8cSStefano Zampini }
35334a97f8cSStefano Zampini 
35434a97f8cSStefano Zampini #undef __FUNCT__
35534a97f8cSStefano Zampini #define __FUNCT__ "PCBDDCScalingDestroy_Deluxe"
35634a97f8cSStefano Zampini static PetscErrorCode PCBDDCScalingDestroy_Deluxe(PC pc)
35734a97f8cSStefano Zampini {
35834a97f8cSStefano Zampini   PC_BDDC*            pcbddc=(PC_BDDC*)pc->data;
35934a97f8cSStefano Zampini   PetscErrorCode      ierr;
36034a97f8cSStefano Zampini 
36134a97f8cSStefano Zampini   PetscFunctionBegin;
36234a97f8cSStefano Zampini   ierr = PCBDDCScalingReset_Deluxe_Solvers(pcbddc->deluxe_ctx);CHKERRQ(ierr);
36334a97f8cSStefano Zampini   ierr = PetscFree(pcbddc->deluxe_ctx);CHKERRQ(ierr);
36434a97f8cSStefano Zampini   PetscFunctionReturn(0);
36534a97f8cSStefano Zampini }
36634a97f8cSStefano Zampini 
36734a97f8cSStefano Zampini #undef __FUNCT__
36834a97f8cSStefano Zampini #define __FUNCT__ "PCBDDCScalingReset_Deluxe_Solvers"
36934a97f8cSStefano Zampini static PetscErrorCode PCBDDCScalingReset_Deluxe_Solvers(PCBDDCDeluxeScaling deluxe_ctx)
37034a97f8cSStefano Zampini {
37172b8c272SStefano Zampini   PetscInt       i;
37234a97f8cSStefano Zampini   PetscErrorCode ierr;
37334a97f8cSStefano Zampini 
37434a97f8cSStefano Zampini   PetscFunctionBegin;
37534a97f8cSStefano Zampini   ierr = PetscFree(deluxe_ctx->idx_simple_B);CHKERRQ(ierr);
37634a97f8cSStefano Zampini   deluxe_ctx->n_simple = 0;
37772b8c272SStefano Zampini   for (i=0;i<deluxe_ctx->seq_n;i++) {
37872b8c272SStefano Zampini     ierr = VecScatterDestroy(&deluxe_ctx->seq_scctx[i]);CHKERRQ(ierr);
37972b8c272SStefano Zampini     ierr = VecDestroy(&deluxe_ctx->seq_work1[i]);CHKERRQ(ierr);
38072b8c272SStefano Zampini     ierr = VecDestroy(&deluxe_ctx->seq_work2[i]);CHKERRQ(ierr);
38172b8c272SStefano Zampini     ierr = MatDestroy(&deluxe_ctx->seq_mat[i]);CHKERRQ(ierr);
38272b8c272SStefano Zampini     ierr = MatDestroy(&deluxe_ctx->seq_mat_inv_sum[i]);CHKERRQ(ierr);
38372b8c272SStefano Zampini   }
38472b8c272SStefano 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);
38572b8c272SStefano Zampini   ierr = PetscFree(deluxe_ctx->workspace);CHKERRQ(ierr);
38672b8c272SStefano Zampini   deluxe_ctx->seq_n = 0;
38734a97f8cSStefano Zampini   PetscFunctionReturn(0);
38834a97f8cSStefano Zampini }
38934a97f8cSStefano Zampini 
39034a97f8cSStefano Zampini #undef __FUNCT__
39134a97f8cSStefano Zampini #define __FUNCT__ "PCBDDCScalingSetUp_Deluxe"
39234a97f8cSStefano Zampini static PetscErrorCode PCBDDCScalingSetUp_Deluxe(PC pc)
39334a97f8cSStefano Zampini {
39434a97f8cSStefano Zampini   PC_IS               *pcis=(PC_IS*)pc->data;
39534a97f8cSStefano Zampini   PC_BDDC             *pcbddc=(PC_BDDC*)pc->data;
39634a97f8cSStefano Zampini   PCBDDCDeluxeScaling deluxe_ctx=pcbddc->deluxe_ctx;
397b96c3477SStefano Zampini   PCBDDCSubSchurs     sub_schurs=pcbddc->sub_schurs;
39834a97f8cSStefano Zampini   PetscErrorCode      ierr;
39934a97f8cSStefano Zampini 
40034a97f8cSStefano Zampini   PetscFunctionBegin;
40172b8c272SStefano Zampini   /* reset data structures if the topology has changed */
40272b8c272SStefano Zampini   if (pcbddc->recompute_topography) {
40334a97f8cSStefano Zampini     ierr = PCBDDCScalingReset_Deluxe_Solvers(deluxe_ctx);CHKERRQ(ierr);
40472b8c272SStefano Zampini   }
40534a97f8cSStefano Zampini 
406b1b3d7a2SStefano Zampini   /* Compute data structures to solve sequential problems */
4075a95e1ceSStefano Zampini   ierr = PCBDDCScalingSetUp_Deluxe_Private(pc);CHKERRQ(ierr);
4085db18549SStefano Zampini 
409b1b3d7a2SStefano Zampini   /* diagonal scaling on interface dofs not contained in cc */
410d62866d3SStefano Zampini   if (sub_schurs->is_vertices || sub_schurs->is_dir) {
4111b968477SStefano Zampini     PetscInt n_com,n_dir;
4121b968477SStefano Zampini     n_com = 0;
413d62866d3SStefano Zampini     if (sub_schurs->is_vertices) {
414d62866d3SStefano Zampini       ierr = ISGetLocalSize(sub_schurs->is_vertices,&n_com);CHKERRQ(ierr);
4151b968477SStefano Zampini     }
4161b968477SStefano Zampini     n_dir = 0;
417d62866d3SStefano Zampini     if (sub_schurs->is_dir) {
418d62866d3SStefano Zampini       ierr = ISGetLocalSize(sub_schurs->is_dir,&n_dir);CHKERRQ(ierr);
4191b968477SStefano Zampini     }
42072b8c272SStefano Zampini     if (!deluxe_ctx->n_simple) {
4211b968477SStefano Zampini       deluxe_ctx->n_simple = n_dir + n_com;
4221b968477SStefano Zampini       ierr = PetscMalloc1(deluxe_ctx->n_simple,&deluxe_ctx->idx_simple_B);CHKERRQ(ierr);
423d62866d3SStefano Zampini       if (sub_schurs->is_vertices) {
4249bb4a8caSStefano Zampini         PetscInt       nmap;
4259bb4a8caSStefano Zampini         const PetscInt *idxs;
4261b968477SStefano Zampini 
427d62866d3SStefano Zampini         ierr = ISGetIndices(sub_schurs->is_vertices,&idxs);CHKERRQ(ierr);
4281b968477SStefano Zampini         ierr = ISGlobalToLocalMappingApply(pcis->BtoNmap,IS_GTOLM_DROP,n_com,idxs,&nmap,deluxe_ctx->idx_simple_B);CHKERRQ(ierr);
429*af25d912SStefano 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);
430d62866d3SStefano Zampini         ierr = ISRestoreIndices(sub_schurs->is_vertices,&idxs);CHKERRQ(ierr);
4311b968477SStefano Zampini       }
432d62866d3SStefano Zampini       if (sub_schurs->is_dir) {
4331b968477SStefano Zampini         PetscInt       nmap;
4341b968477SStefano Zampini         const PetscInt *idxs;
4351b968477SStefano Zampini 
436d62866d3SStefano Zampini         ierr = ISGetIndices(sub_schurs->is_dir,&idxs);CHKERRQ(ierr);
4371b968477SStefano Zampini         ierr = ISGlobalToLocalMappingApply(pcis->BtoNmap,IS_GTOLM_DROP,n_dir,idxs,&nmap,deluxe_ctx->idx_simple_B+n_com);CHKERRQ(ierr);
438*af25d912SStefano 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);
439d62866d3SStefano Zampini         ierr = ISRestoreIndices(sub_schurs->is_dir,&idxs);CHKERRQ(ierr);
4401b968477SStefano Zampini       }
4411b968477SStefano Zampini       ierr = PetscSortInt(deluxe_ctx->n_simple,deluxe_ctx->idx_simple_B);CHKERRQ(ierr);
4429bb4a8caSStefano Zampini     } else {
443*af25d912SStefano 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);
44472b8c272SStefano Zampini     }
44572b8c272SStefano Zampini   } else {
446b1b3d7a2SStefano Zampini     deluxe_ctx->n_simple = 0;
4479bb4a8caSStefano Zampini     deluxe_ctx->idx_simple_B = 0;
448b1b3d7a2SStefano Zampini   }
44934a97f8cSStefano Zampini   PetscFunctionReturn(0);
45034a97f8cSStefano Zampini }
45134a97f8cSStefano Zampini 
45234a97f8cSStefano Zampini #undef __FUNCT__
4535a95e1ceSStefano Zampini #define __FUNCT__ "PCBDDCScalingSetUp_Deluxe_Private"
4545a95e1ceSStefano Zampini static PetscErrorCode PCBDDCScalingSetUp_Deluxe_Private(PC pc)
4555db18549SStefano Zampini {
4565db18549SStefano Zampini   PC_BDDC                *pcbddc=(PC_BDDC*)pc->data;
4575db18549SStefano Zampini   PCBDDCDeluxeScaling    deluxe_ctx=pcbddc->deluxe_ctx;
458b96c3477SStefano Zampini   PCBDDCSubSchurs        sub_schurs = pcbddc->sub_schurs;
45972b8c272SStefano Zampini   PetscScalar            *matdata,*matdata2;
46072b8c272SStefano Zampini   PetscInt               i,max_subset_size,cum,cum2;
46172b8c272SStefano Zampini   const PetscInt         *idxs;
46272b8c272SStefano Zampini   PetscBool              newsetup = PETSC_FALSE;
4635db18549SStefano Zampini   PetscErrorCode         ierr;
4645db18549SStefano Zampini 
4655db18549SStefano Zampini   PetscFunctionBegin;
4665a95e1ceSStefano Zampini   if (!sub_schurs->n_subs) {
4679221af80SStefano Zampini     PetscFunctionReturn(0);
4689221af80SStefano Zampini   }
4699221af80SStefano Zampini 
47072b8c272SStefano Zampini   /* Allocate arrays for subproblems */
47172b8c272SStefano Zampini   if (!deluxe_ctx->seq_n) {
47272b8c272SStefano Zampini     deluxe_ctx->seq_n = sub_schurs->n_subs;
47372b8c272SStefano 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);
47472b8c272SStefano Zampini     newsetup = PETSC_TRUE;
47572b8c272SStefano Zampini   } else if (deluxe_ctx->seq_n != sub_schurs->n_subs) {
47672b8c272SStefano 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);
47772b8c272SStefano Zampini   }
47872b8c272SStefano Zampini   deluxe_ctx->change_with_qr = sub_schurs->change_with_qr;
4795db18549SStefano Zampini 
48072b8c272SStefano Zampini   /* Create objects for deluxe */
48172b8c272SStefano Zampini   max_subset_size = 0;
48272b8c272SStefano Zampini   for (i=0;i<sub_schurs->n_subs;i++) {
48372b8c272SStefano Zampini     PetscInt subset_size;
48472b8c272SStefano Zampini     ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr);
48572b8c272SStefano Zampini     max_subset_size = PetscMax(subset_size,max_subset_size);
48672b8c272SStefano Zampini   }
48772b8c272SStefano Zampini   if (newsetup) {
48872b8c272SStefano Zampini     ierr = PetscMalloc1(2*max_subset_size,&deluxe_ctx->workspace);CHKERRQ(ierr);
48972b8c272SStefano Zampini   }
49072b8c272SStefano Zampini   cum = cum2 = 0;
49172b8c272SStefano Zampini   ierr = ISGetIndices(sub_schurs->is_Ej_all,&idxs);CHKERRQ(ierr);
49272b8c272SStefano Zampini   ierr = MatSeqAIJGetArray(sub_schurs->S_Ej_all,&matdata);CHKERRQ(ierr);
49372b8c272SStefano Zampini   matdata2 = NULL;
49472b8c272SStefano Zampini   if (sub_schurs->sum_S_Ej_all) {
49572b8c272SStefano Zampini     ierr = MatSeqAIJGetArray(sub_schurs->sum_S_Ej_all,&matdata2);CHKERRQ(ierr);
49672b8c272SStefano Zampini   }
49772b8c272SStefano Zampini   for (i=0;i<deluxe_ctx->seq_n;i++) {
49872b8c272SStefano Zampini     PetscInt     subset_size;
499925dfd53SStefano Zampini 
50072b8c272SStefano Zampini     ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr);
50172b8c272SStefano Zampini     if (newsetup) {
50272b8c272SStefano Zampini       IS  sub;
50372b8c272SStefano Zampini       /* work vectors */
50472b8c272SStefano Zampini       ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,1,subset_size,deluxe_ctx->workspace,&deluxe_ctx->seq_work1[i]);CHKERRQ(ierr);
50572b8c272SStefano Zampini       ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,1,subset_size,deluxe_ctx->workspace+subset_size,&deluxe_ctx->seq_work2[i]);CHKERRQ(ierr);
50672b8c272SStefano Zampini 
50772b8c272SStefano Zampini       /* scatters */
50872b8c272SStefano Zampini       ierr = ISCreateGeneral(PETSC_COMM_SELF,subset_size,idxs+cum,PETSC_COPY_VALUES,&sub);CHKERRQ(ierr);
50972b8c272SStefano Zampini       ierr = VecScatterCreate(pcbddc->work_scaling,sub,deluxe_ctx->seq_work1[i],NULL,&deluxe_ctx->seq_scctx[i]);CHKERRQ(ierr);
51072b8c272SStefano Zampini       ierr = ISDestroy(&sub);CHKERRQ(ierr);
51172b8c272SStefano Zampini 
51272b8c272SStefano Zampini       /* S_E_j */
51372b8c272SStefano Zampini       ierr = MatCreateSeqDense(PETSC_COMM_SELF,subset_size,subset_size,matdata+cum2,&deluxe_ctx->seq_mat[i]);CHKERRQ(ierr);
51472b8c272SStefano Zampini     }
51572b8c272SStefano Zampini     /* \sum_k S^k_E_j */
51672b8c272SStefano Zampini     if (matdata2) {
51772b8c272SStefano Zampini       ierr = MatDestroy(&deluxe_ctx->seq_mat_inv_sum[i]);CHKERRQ(ierr);
51872b8c272SStefano Zampini       ierr = MatCreateSeqDense(PETSC_COMM_SELF,subset_size,subset_size,matdata2+cum2,&deluxe_ctx->seq_mat_inv_sum[i]);CHKERRQ(ierr);
51972b8c272SStefano Zampini       if (sub_schurs->is_hermitian && sub_schurs->is_posdef) {
52072b8c272SStefano Zampini         ierr = MatCholeskyFactor(deluxe_ctx->seq_mat_inv_sum[i],NULL,NULL);CHKERRQ(ierr);
521d62866d3SStefano Zampini       } else {
52272b8c272SStefano Zampini         ierr = MatLUFactor(deluxe_ctx->seq_mat_inv_sum[i],NULL,NULL,NULL);CHKERRQ(ierr);
523d62866d3SStefano Zampini       }
52472b8c272SStefano Zampini     }
52572b8c272SStefano Zampini     cum += subset_size;
52672b8c272SStefano Zampini     cum2 += subset_size*subset_size;
52772b8c272SStefano Zampini   }
52872b8c272SStefano Zampini   ierr = ISRestoreIndices(sub_schurs->is_Ej_all,&idxs);CHKERRQ(ierr);
52972b8c272SStefano Zampini   ierr = MatSeqAIJRestoreArray(sub_schurs->S_Ej_all,&matdata);CHKERRQ(ierr);
53072b8c272SStefano Zampini   if (sub_schurs->sum_S_Ej_all) {
53172b8c272SStefano Zampini     ierr = MatSeqAIJRestoreArray(sub_schurs->sum_S_Ej_all,&matdata2);CHKERRQ(ierr);
532925dfd53SStefano Zampini   }
5335db18549SStefano Zampini 
53472b8c272SStefano Zampini   /* the change of basis is just a reference to sub_schurs->change (if any) */
53516e386b8SStefano Zampini   deluxe_ctx->change = sub_schurs->change;
53672b8c272SStefano Zampini   if (deluxe_ctx->change && !deluxe_ctx->change_with_qr) {
53772b8c272SStefano Zampini     for (i=0;i<deluxe_ctx->seq_n;i++) {
53872b8c272SStefano Zampini       if (newsetup) {
53972b8c272SStefano Zampini         PC pc;
54072b8c272SStefano Zampini 
54172b8c272SStefano Zampini         ierr = KSPGetPC(deluxe_ctx->change[i],&pc);CHKERRQ(ierr);
54272b8c272SStefano Zampini         ierr = PCSetType(pc,PCLU);CHKERRQ(ierr);
54372b8c272SStefano Zampini         ierr = KSPSetFromOptions(deluxe_ctx->change[i]);CHKERRQ(ierr);
54472b8c272SStefano Zampini       }
54572b8c272SStefano Zampini       ierr = KSPSetUp(deluxe_ctx->change[i]);CHKERRQ(ierr);
54672b8c272SStefano Zampini     }
54716e386b8SStefano Zampini   }
5485db18549SStefano Zampini   PetscFunctionReturn(0);
5495db18549SStefano Zampini }
550