xref: /petsc/src/ksp/pc/impls/bddc/bddcscalingbasic.c (revision 16e386b8f3e90b5e6399ce06f317f0b6cb4eaacd)
1ab5c6b0cSJed Brown #include <../src/ksp/pc/impls/bddc/bddc.h>
2ab5c6b0cSJed Brown #include <../src/ksp/pc/impls/bddc/bddcprivate.h>
3674ae819SStefano Zampini 
434a97f8cSStefano Zampini /* prototypes for deluxe functions */
534a97f8cSStefano Zampini static PetscErrorCode PCBDDCScalingCreate_Deluxe(PC);
634a97f8cSStefano Zampini static PetscErrorCode PCBDDCScalingDestroy_Deluxe(PC);
734a97f8cSStefano Zampini static PetscErrorCode PCBDDCScalingSetUp_Deluxe(PC);
85a95e1ceSStefano Zampini static PetscErrorCode PCBDDCScalingSetUp_Deluxe_Private(PC);
934a97f8cSStefano Zampini static PetscErrorCode PCBDDCScalingReset_Deluxe_Solvers(PCBDDCDeluxeScaling);
10674ae819SStefano Zampini 
11674ae819SStefano Zampini #undef __FUNCT__
12674ae819SStefano Zampini #define __FUNCT__ "PCBDDCScalingExtension_Basic"
13674ae819SStefano Zampini static PetscErrorCode PCBDDCScalingExtension_Basic(PC pc, Vec local_interface_vector, Vec global_vector)
14674ae819SStefano Zampini {
15674ae819SStefano Zampini   PC_IS*         pcis = (PC_IS*)pc->data;
16674ae819SStefano Zampini   PC_BDDC*       pcbddc = (PC_BDDC*)pc->data;
17674ae819SStefano Zampini   PetscErrorCode ierr;
18674ae819SStefano Zampini 
19674ae819SStefano Zampini   PetscFunctionBegin;
20674ae819SStefano Zampini   /* Apply partition of unity */
21674ae819SStefano Zampini   ierr = VecPointwiseMult(pcbddc->work_scaling,pcis->D,local_interface_vector);CHKERRQ(ierr);
22674ae819SStefano Zampini   ierr = VecSet(global_vector,0.0);CHKERRQ(ierr);
23674ae819SStefano Zampini   ierr = VecScatterBegin(pcis->global_to_B,pcbddc->work_scaling,global_vector,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
24674ae819SStefano Zampini   ierr = VecScatterEnd(pcis->global_to_B,pcbddc->work_scaling,global_vector,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
25674ae819SStefano Zampini   PetscFunctionReturn(0);
26674ae819SStefano Zampini }
27674ae819SStefano Zampini 
28674ae819SStefano Zampini #undef __FUNCT__
29674ae819SStefano Zampini #define __FUNCT__ "PCBDDCScalingExtension_Deluxe"
30674ae819SStefano Zampini static PetscErrorCode PCBDDCScalingExtension_Deluxe(PC pc, Vec x, Vec y)
31674ae819SStefano Zampini {
32674ae819SStefano Zampini   PC_IS*              pcis=(PC_IS*)pc->data;
33674ae819SStefano Zampini   PC_BDDC*            pcbddc=(PC_BDDC*)pc->data;
34674ae819SStefano Zampini   PCBDDCDeluxeScaling deluxe_ctx = pcbddc->deluxe_ctx;
35674ae819SStefano Zampini   PetscErrorCode      ierr;
36674ae819SStefano Zampini 
37674ae819SStefano Zampini   PetscFunctionBegin;
385a95e1ceSStefano Zampini   ierr = VecSet(pcbddc->work_scaling,0.0);CHKERRQ(ierr);
395a95e1ceSStefano Zampini   ierr = VecSet(y,0.0);CHKERRQ(ierr);
405a95e1ceSStefano Zampini   if (deluxe_ctx->n_simple) { /* scale deluxe vertices using diagonal scaling */
415a95e1ceSStefano Zampini     PetscInt          i;
422b095fd8SStefano Zampini     const PetscScalar *array_x,*array_D;
432b095fd8SStefano Zampini     PetscScalar       *array;
442b095fd8SStefano Zampini     ierr = VecGetArrayRead(x,&array_x);CHKERRQ(ierr);
452b095fd8SStefano Zampini     ierr = VecGetArrayRead(pcis->D,&array_D);CHKERRQ(ierr);
46674ae819SStefano Zampini     ierr = VecGetArray(pcbddc->work_scaling,&array);CHKERRQ(ierr);
47674ae819SStefano Zampini     for (i=0;i<deluxe_ctx->n_simple;i++) {
48674ae819SStefano Zampini       array[deluxe_ctx->idx_simple_B[i]] = array_x[deluxe_ctx->idx_simple_B[i]]*array_D[deluxe_ctx->idx_simple_B[i]];
49674ae819SStefano Zampini     }
50674ae819SStefano Zampini     ierr = VecRestoreArray(pcbddc->work_scaling,&array);CHKERRQ(ierr);
512b095fd8SStefano Zampini     ierr = VecRestoreArrayRead(pcis->D,&array_D);CHKERRQ(ierr);
522b095fd8SStefano Zampini     ierr = VecRestoreArrayRead(x,&array_x);CHKERRQ(ierr);
5334a97f8cSStefano Zampini   }
54ac632422SStefano Zampini   /* sequential part : all problems and Schur applications collapsed into a single matrix vector multiplication or a matvec and a solve */
555a95e1ceSStefano Zampini   if (deluxe_ctx->seq_mat) {
56674ae819SStefano Zampini     ierr = VecScatterBegin(deluxe_ctx->seq_scctx,x,deluxe_ctx->seq_work1,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
57674ae819SStefano Zampini     ierr = VecScatterEnd(deluxe_ctx->seq_scctx,x,deluxe_ctx->seq_work1,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
58*16e386b8SStefano Zampini     if (deluxe_ctx->change) {
59*16e386b8SStefano Zampini       ierr = KSPSolve(deluxe_ctx->change,deluxe_ctx->seq_work1,deluxe_ctx->seq_work1);CHKERRQ(ierr);
60*16e386b8SStefano Zampini     }
615a95e1ceSStefano Zampini     ierr = MatMultTranspose(deluxe_ctx->seq_mat,deluxe_ctx->seq_work1,deluxe_ctx->seq_work2);CHKERRQ(ierr);
62ac632422SStefano Zampini     if (deluxe_ctx->seq_ksp) {
63ac632422SStefano Zampini       ierr = KSPSolveTranspose(deluxe_ctx->seq_ksp,deluxe_ctx->seq_work2,deluxe_ctx->seq_work2);CHKERRQ(ierr);
64ac632422SStefano Zampini     }
65*16e386b8SStefano Zampini     if (deluxe_ctx->change) {
66*16e386b8SStefano Zampini       Mat change;
67*16e386b8SStefano Zampini 
68*16e386b8SStefano Zampini       ierr = KSPGetOperators(deluxe_ctx->change,&change,NULL);CHKERRQ(ierr);
69*16e386b8SStefano Zampini       ierr = MatMult(change,deluxe_ctx->seq_work2,deluxe_ctx->seq_work1);CHKERRQ(ierr);
70*16e386b8SStefano Zampini       ierr = VecScatterBegin(deluxe_ctx->seq_scctx,deluxe_ctx->seq_work1,pcbddc->work_scaling,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
71*16e386b8SStefano Zampini       ierr = VecScatterEnd(deluxe_ctx->seq_scctx,deluxe_ctx->seq_work1,pcbddc->work_scaling,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
72*16e386b8SStefano Zampini     } else {
735a95e1ceSStefano Zampini       ierr = VecScatterBegin(deluxe_ctx->seq_scctx,deluxe_ctx->seq_work2,pcbddc->work_scaling,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
745a95e1ceSStefano Zampini       ierr = VecScatterEnd(deluxe_ctx->seq_scctx,deluxe_ctx->seq_work2,pcbddc->work_scaling,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
75674ae819SStefano Zampini     }
76*16e386b8SStefano Zampini   }
77674ae819SStefano Zampini   /* put local boundary part in global vector */
78674ae819SStefano Zampini   ierr = VecScatterBegin(pcis->global_to_B,pcbddc->work_scaling,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
79674ae819SStefano Zampini   ierr = VecScatterEnd(pcis->global_to_B,pcbddc->work_scaling,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
80674ae819SStefano Zampini   PetscFunctionReturn(0);
81674ae819SStefano Zampini }
82674ae819SStefano Zampini 
83674ae819SStefano Zampini #undef __FUNCT__
84674ae819SStefano Zampini #define __FUNCT__ "PCBDDCScalingExtension"
85674ae819SStefano Zampini PetscErrorCode PCBDDCScalingExtension(PC pc, Vec local_interface_vector, Vec global_vector)
86674ae819SStefano Zampini {
87674ae819SStefano Zampini   PC_BDDC *pcbddc=(PC_BDDC*)pc->data;
88674ae819SStefano Zampini   PetscErrorCode ierr;
89674ae819SStefano Zampini 
90674ae819SStefano Zampini   PetscFunctionBegin;
91674ae819SStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
92674ae819SStefano Zampini   PetscValidHeaderSpecific(local_interface_vector,VEC_CLASSID,2);
93674ae819SStefano Zampini   PetscValidHeaderSpecific(global_vector,VEC_CLASSID,3);
94674ae819SStefano Zampini   if (local_interface_vector == pcbddc->work_scaling) {
95a07ea27aSStefano Zampini     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Local vector cannot be pcbddc->work_scaling!\n");
96674ae819SStefano Zampini   }
97674ae819SStefano Zampini   ierr = PetscTryMethod(pc,"PCBDDCScalingExtension_C",(PC,Vec,Vec),(pc,local_interface_vector,global_vector));CHKERRQ(ierr);
98674ae819SStefano Zampini   PetscFunctionReturn(0);
99674ae819SStefano Zampini }
100674ae819SStefano Zampini 
101674ae819SStefano Zampini #undef __FUNCT__
102674ae819SStefano Zampini #define __FUNCT__ "PCBDDCScalingRestriction_Basic"
103674ae819SStefano Zampini static PetscErrorCode PCBDDCScalingRestriction_Basic(PC pc, Vec global_vector, Vec local_interface_vector)
104674ae819SStefano Zampini {
105674ae819SStefano Zampini   PetscErrorCode ierr;
106674ae819SStefano Zampini   PC_IS* pcis = (PC_IS*)pc->data;
107674ae819SStefano Zampini 
108674ae819SStefano Zampini   PetscFunctionBegin;
109674ae819SStefano Zampini   ierr = VecScatterBegin(pcis->global_to_B,global_vector,local_interface_vector,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
110674ae819SStefano Zampini   ierr = VecScatterEnd(pcis->global_to_B,global_vector,local_interface_vector,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
111674ae819SStefano Zampini   /* Apply partition of unity */
112674ae819SStefano Zampini   ierr = VecPointwiseMult(local_interface_vector,pcis->D,local_interface_vector);CHKERRQ(ierr);
113674ae819SStefano Zampini   PetscFunctionReturn(0);
114674ae819SStefano Zampini }
115674ae819SStefano Zampini 
116674ae819SStefano Zampini #undef __FUNCT__
117674ae819SStefano Zampini #define __FUNCT__ "PCBDDCScalingRestriction_Deluxe"
118674ae819SStefano Zampini static PetscErrorCode PCBDDCScalingRestriction_Deluxe(PC pc, Vec x, Vec y)
119674ae819SStefano Zampini {
120674ae819SStefano Zampini   PC_IS*              pcis=(PC_IS*)pc->data;
121674ae819SStefano Zampini   PC_BDDC*            pcbddc=(PC_BDDC*)pc->data;
122674ae819SStefano Zampini   PCBDDCDeluxeScaling deluxe_ctx = pcbddc->deluxe_ctx;
123674ae819SStefano Zampini   PetscErrorCode      ierr;
124674ae819SStefano Zampini 
125674ae819SStefano Zampini   PetscFunctionBegin;
126674ae819SStefano Zampini   /* get local boundary part of global vector */
127674ae819SStefano Zampini   ierr = VecScatterBegin(pcis->global_to_B,x,y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
128674ae819SStefano Zampini   ierr = VecScatterEnd(pcis->global_to_B,x,y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1295a95e1ceSStefano Zampini   if (deluxe_ctx->n_simple) { /* scale deluxe vertices using diagonal scaling */
1305a95e1ceSStefano Zampini     PetscInt          i;
1312b095fd8SStefano Zampini     PetscScalar       *array_y;
1322b095fd8SStefano Zampini     const PetscScalar *array_D;
133674ae819SStefano Zampini     ierr = VecGetArray(y,&array_y);CHKERRQ(ierr);
1342b095fd8SStefano Zampini     ierr = VecGetArrayRead(pcis->D,&array_D);CHKERRQ(ierr);
135674ae819SStefano Zampini     for (i=0;i<deluxe_ctx->n_simple;i++) {
136674ae819SStefano Zampini       array_y[deluxe_ctx->idx_simple_B[i]] *= array_D[deluxe_ctx->idx_simple_B[i]];
137674ae819SStefano Zampini     }
1382b095fd8SStefano Zampini     ierr = VecRestoreArrayRead(pcis->D,&array_D);CHKERRQ(ierr);
139674ae819SStefano Zampini     ierr = VecRestoreArray(y,&array_y);CHKERRQ(ierr);
14034a97f8cSStefano Zampini   }
14134a97f8cSStefano Zampini   /* sequential part : all problems and Schur applications collapsed into a single matrix vector multiplication and ksp solution */
1425a95e1ceSStefano Zampini   if (deluxe_ctx->seq_mat) {
143*16e386b8SStefano Zampini     if (deluxe_ctx->change) {
144*16e386b8SStefano Zampini       Mat change;
145*16e386b8SStefano Zampini 
146*16e386b8SStefano Zampini       ierr = VecScatterBegin(deluxe_ctx->seq_scctx,y,deluxe_ctx->seq_work2,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
147*16e386b8SStefano Zampini       ierr = VecScatterEnd(deluxe_ctx->seq_scctx,y,deluxe_ctx->seq_work2,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
148*16e386b8SStefano Zampini       ierr = KSPGetOperators(deluxe_ctx->change,&change,NULL);CHKERRQ(ierr);
149*16e386b8SStefano Zampini       ierr = MatMultTranspose(change,deluxe_ctx->seq_work2,deluxe_ctx->seq_work1);CHKERRQ(ierr);
150*16e386b8SStefano Zampini     } else {
151674ae819SStefano Zampini       ierr = VecScatterBegin(deluxe_ctx->seq_scctx,y,deluxe_ctx->seq_work1,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
152674ae819SStefano Zampini       ierr = VecScatterEnd(deluxe_ctx->seq_scctx,y,deluxe_ctx->seq_work1,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
153*16e386b8SStefano Zampini     }
154ac632422SStefano Zampini     if (deluxe_ctx->seq_ksp) {
155ac632422SStefano Zampini       ierr = KSPSolve(deluxe_ctx->seq_ksp,deluxe_ctx->seq_work1,deluxe_ctx->seq_work1);CHKERRQ(ierr);
156ac632422SStefano Zampini     }
1575a95e1ceSStefano Zampini     ierr = MatMult(deluxe_ctx->seq_mat,deluxe_ctx->seq_work1,deluxe_ctx->seq_work2);CHKERRQ(ierr);
158*16e386b8SStefano Zampini     if (deluxe_ctx->change) {
159*16e386b8SStefano Zampini       ierr = KSPSolveTranspose(deluxe_ctx->change,deluxe_ctx->seq_work2,deluxe_ctx->seq_work2);CHKERRQ(ierr);
160*16e386b8SStefano Zampini     }
1615a95e1ceSStefano Zampini     ierr = VecScatterBegin(deluxe_ctx->seq_scctx,deluxe_ctx->seq_work2,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1625a95e1ceSStefano Zampini     ierr = VecScatterEnd(deluxe_ctx->seq_scctx,deluxe_ctx->seq_work2,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
163674ae819SStefano Zampini   }
164674ae819SStefano Zampini   PetscFunctionReturn(0);
165674ae819SStefano Zampini }
166674ae819SStefano Zampini 
167674ae819SStefano Zampini #undef __FUNCT__
168674ae819SStefano Zampini #define __FUNCT__ "PCBDDCScalingRestriction"
169674ae819SStefano Zampini PetscErrorCode PCBDDCScalingRestriction(PC pc, Vec global_vector, Vec local_interface_vector)
170674ae819SStefano Zampini {
171674ae819SStefano Zampini   PC_BDDC        *pcbddc=(PC_BDDC*)pc->data;
172674ae819SStefano Zampini   PetscErrorCode ierr;
173674ae819SStefano Zampini 
174674ae819SStefano Zampini   PetscFunctionBegin;
175674ae819SStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
176674ae819SStefano Zampini   PetscValidHeaderSpecific(global_vector,VEC_CLASSID,2);
177674ae819SStefano Zampini   PetscValidHeaderSpecific(local_interface_vector,VEC_CLASSID,3);
178674ae819SStefano Zampini   if (local_interface_vector == pcbddc->work_scaling) {
179a07ea27aSStefano Zampini     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Local vector cannot be pcbddc->work_scaling!\n");
180674ae819SStefano Zampini   }
181674ae819SStefano Zampini   ierr = PetscTryMethod(pc,"PCBDDCScalingRestriction_C",(PC,Vec,Vec),(pc,global_vector,local_interface_vector));CHKERRQ(ierr);
182674ae819SStefano Zampini   PetscFunctionReturn(0);
183674ae819SStefano Zampini }
184674ae819SStefano Zampini 
185674ae819SStefano Zampini #undef __FUNCT__
186674ae819SStefano Zampini #define __FUNCT__ "PCBDDCScalingSetUp"
187674ae819SStefano Zampini PetscErrorCode PCBDDCScalingSetUp(PC pc)
188674ae819SStefano Zampini {
189674ae819SStefano Zampini   PC_IS*         pcis=(PC_IS*)pc->data;
190674ae819SStefano Zampini   PC_BDDC*       pcbddc=(PC_BDDC*)pc->data;
191674ae819SStefano Zampini   PetscErrorCode ierr;
192674ae819SStefano Zampini 
193674ae819SStefano Zampini   PetscFunctionBegin;
194674ae819SStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
195674ae819SStefano Zampini   /* create work vector for the operator */
19634a97f8cSStefano Zampini   ierr = VecDestroy(&pcbddc->work_scaling);CHKERRQ(ierr);
197674ae819SStefano Zampini   ierr = VecDuplicate(pcis->vec1_B,&pcbddc->work_scaling);CHKERRQ(ierr);
19834a97f8cSStefano Zampini   /* always rebuild pcis->D */
19928d874f6SStefano Zampini   if (pcis->use_stiffness_scaling) {
200674ae819SStefano Zampini     ierr = MatGetDiagonal(pcbddc->local_mat,pcis->vec1_N);CHKERRQ(ierr);
201674ae819SStefano Zampini     ierr = VecScatterBegin(pcis->N_to_B,pcis->vec1_N,pcis->D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
202674ae819SStefano Zampini     ierr = VecScatterEnd(pcis->N_to_B,pcis->vec1_N,pcis->D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
203674ae819SStefano Zampini   }
204674ae819SStefano Zampini   ierr = VecCopy(pcis->D,pcis->vec1_B);CHKERRQ(ierr);
205674ae819SStefano Zampini   ierr = VecSet(pcis->vec1_global,0.0);CHKERRQ(ierr);
206674ae819SStefano Zampini   ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_B,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
207674ae819SStefano Zampini   ierr = VecScatterEnd(pcis->global_to_B,pcis->vec1_B,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
208674ae819SStefano Zampini   ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_global,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
209674ae819SStefano Zampini   ierr = VecScatterEnd(pcis->global_to_B,pcis->vec1_global,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
210674ae819SStefano Zampini   ierr = VecPointwiseDivide(pcis->D,pcis->D,pcis->vec1_B);CHKERRQ(ierr);
211674ae819SStefano Zampini   /* now setup */
212681e7c04SStefano Zampini   if (pcbddc->use_deluxe_scaling) {
21334a97f8cSStefano Zampini     if (!pcbddc->deluxe_ctx) {
21434a97f8cSStefano Zampini       ierr = PCBDDCScalingCreate_Deluxe(pc);CHKERRQ(ierr);
21534a97f8cSStefano Zampini     }
21634a97f8cSStefano Zampini     ierr = PCBDDCScalingSetUp_Deluxe(pc);CHKERRQ(ierr);
217674ae819SStefano Zampini     ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCScalingRestriction_C",PCBDDCScalingRestriction_Deluxe);CHKERRQ(ierr);
218674ae819SStefano Zampini     ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCScalingExtension_C",PCBDDCScalingExtension_Deluxe);CHKERRQ(ierr);
219674ae819SStefano Zampini   } else {
220674ae819SStefano Zampini     ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCScalingRestriction_C",PCBDDCScalingRestriction_Basic);CHKERRQ(ierr);
221674ae819SStefano Zampini     ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCScalingExtension_C",PCBDDCScalingExtension_Basic);CHKERRQ(ierr);
222674ae819SStefano Zampini   }
22334a97f8cSStefano Zampini 
224674ae819SStefano Zampini   /* test */
225674ae819SStefano Zampini   if (pcbddc->dbg_flag) {
22634a97f8cSStefano Zampini     Vec         vec2_global;
227674ae819SStefano Zampini     PetscViewer viewer=pcbddc->dbg_viewer;
22834a97f8cSStefano Zampini     PetscReal   error;
229674ae819SStefano Zampini 
230674ae819SStefano Zampini     /* extension -> from local to parallel */
23134a97f8cSStefano Zampini     ierr = VecSet(pcis->vec1_global,0.0);CHKERRQ(ierr);
23234a97f8cSStefano Zampini     ierr = VecSetRandom(pcis->vec1_B,NULL);CHKERRQ(ierr);
23334a97f8cSStefano Zampini     ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_B,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
23434a97f8cSStefano Zampini     ierr = VecScatterEnd(pcis->global_to_B,pcis->vec1_B,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
23534a97f8cSStefano Zampini     ierr = VecDuplicate(pcis->vec1_global,&vec2_global);CHKERRQ(ierr);
23634a97f8cSStefano Zampini     ierr = VecCopy(pcis->vec1_global,vec2_global);CHKERRQ(ierr);
23734a97f8cSStefano Zampini 
238674ae819SStefano Zampini     ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_global,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
239674ae819SStefano Zampini     ierr = VecScatterEnd(pcis->global_to_B,pcis->vec1_global,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
240674ae819SStefano Zampini     ierr = PCBDDCScalingExtension(pc,pcis->vec1_B,pcis->vec1_global);CHKERRQ(ierr);
24134a97f8cSStefano Zampini     ierr = VecAXPY(pcis->vec1_global,-1.0,vec2_global);CHKERRQ(ierr);
24234a97f8cSStefano Zampini     ierr = VecNorm(pcis->vec1_global,NORM_INFINITY,&error);CHKERRQ(ierr);
243674ae819SStefano Zampini     ierr = PetscViewerASCIIPrintf(viewer,"Error scaling extension %1.14e\n",error);CHKERRQ(ierr);
24434a97f8cSStefano Zampini     if (error>1.e-8 && pcbddc->dbg_flag>1) {
245674ae819SStefano Zampini       ierr = VecView(pcis->vec1_global,viewer);CHKERRQ(ierr);
246674ae819SStefano Zampini     }
24734a97f8cSStefano Zampini     ierr = VecDestroy(&vec2_global);CHKERRQ(ierr);
24834a97f8cSStefano Zampini 
249674ae819SStefano Zampini     /* restriction -> from parallel to local */
250674ae819SStefano Zampini     ierr = VecSet(pcis->vec1_global,0.0);CHKERRQ(ierr);
25134a97f8cSStefano Zampini     ierr = VecSetRandom(pcis->vec1_B,NULL);CHKERRQ(ierr);
252674ae819SStefano Zampini     ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_B,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
253674ae819SStefano Zampini     ierr = VecScatterEnd(pcis->global_to_B,pcis->vec1_B,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
25434a97f8cSStefano Zampini 
25534a97f8cSStefano Zampini     ierr = PCBDDCScalingRestriction(pc,pcis->vec1_global,pcis->vec1_B);CHKERRQ(ierr);
25634a97f8cSStefano Zampini     ierr = VecScale(pcis->vec1_B,-1.0);CHKERRQ(ierr);
25734a97f8cSStefano Zampini     ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_B,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
25834a97f8cSStefano Zampini     ierr = VecScatterEnd(pcis->global_to_B,pcis->vec1_B,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
25934a97f8cSStefano Zampini     ierr = VecNorm(pcis->vec1_global,NORM_INFINITY,&error);CHKERRQ(ierr);
26034a97f8cSStefano Zampini     ierr = PetscViewerASCIIPrintf(viewer,"Error scaling restriction %1.14e\n",error);CHKERRQ(ierr);
26134a97f8cSStefano Zampini     if (error>1.e-8 && pcbddc->dbg_flag>1) {
262674ae819SStefano Zampini       ierr = VecView(pcis->vec1_global,viewer);CHKERRQ(ierr);
263674ae819SStefano Zampini     }
264674ae819SStefano Zampini   }
265674ae819SStefano Zampini   PetscFunctionReturn(0);
266674ae819SStefano Zampini }
267674ae819SStefano Zampini 
268674ae819SStefano Zampini #undef __FUNCT__
269674ae819SStefano Zampini #define __FUNCT__ "PCBDDCScalingDestroy"
270674ae819SStefano Zampini PetscErrorCode PCBDDCScalingDestroy(PC pc)
271674ae819SStefano Zampini {
272674ae819SStefano Zampini   PC_BDDC* pcbddc=(PC_BDDC*)pc->data;
273674ae819SStefano Zampini   PetscErrorCode ierr;
274674ae819SStefano Zampini 
275674ae819SStefano Zampini   PetscFunctionBegin;
27634a97f8cSStefano Zampini   if (pcbddc->deluxe_ctx) {
27734a97f8cSStefano Zampini     ierr = PCBDDCScalingDestroy_Deluxe(pc);CHKERRQ(ierr);
278674ae819SStefano Zampini   }
279674ae819SStefano Zampini   ierr = VecDestroy(&pcbddc->work_scaling);CHKERRQ(ierr);
280674ae819SStefano Zampini   /* remove functions */
281674ae819SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCScalingRestriction_C",NULL);CHKERRQ(ierr);
282674ae819SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCScalingExtension_C",NULL);CHKERRQ(ierr);
283674ae819SStefano Zampini   PetscFunctionReturn(0);
284674ae819SStefano Zampini }
285674ae819SStefano Zampini 
28634a97f8cSStefano Zampini #undef __FUNCT__
28734a97f8cSStefano Zampini #define __FUNCT__ "PCBDDCScalingCreate_Deluxe"
28834a97f8cSStefano Zampini static PetscErrorCode PCBDDCScalingCreate_Deluxe(PC pc)
28934a97f8cSStefano Zampini {
29034a97f8cSStefano Zampini   PC_BDDC*            pcbddc=(PC_BDDC*)pc->data;
29134a97f8cSStefano Zampini   PCBDDCDeluxeScaling deluxe_ctx;
29234a97f8cSStefano Zampini   PetscErrorCode      ierr;
29334a97f8cSStefano Zampini 
29434a97f8cSStefano Zampini   PetscFunctionBegin;
29534a97f8cSStefano Zampini   ierr = PetscNew(&deluxe_ctx);CHKERRQ(ierr);
29634a97f8cSStefano Zampini   pcbddc->deluxe_ctx = deluxe_ctx;
29734a97f8cSStefano Zampini   PetscFunctionReturn(0);
29834a97f8cSStefano Zampini }
29934a97f8cSStefano Zampini 
30034a97f8cSStefano Zampini #undef __FUNCT__
30134a97f8cSStefano Zampini #define __FUNCT__ "PCBDDCScalingDestroy_Deluxe"
30234a97f8cSStefano Zampini static PetscErrorCode PCBDDCScalingDestroy_Deluxe(PC pc)
30334a97f8cSStefano Zampini {
30434a97f8cSStefano Zampini   PC_BDDC*            pcbddc=(PC_BDDC*)pc->data;
30534a97f8cSStefano Zampini   PetscErrorCode      ierr;
30634a97f8cSStefano Zampini 
30734a97f8cSStefano Zampini   PetscFunctionBegin;
30834a97f8cSStefano Zampini   ierr = PCBDDCScalingReset_Deluxe_Solvers(pcbddc->deluxe_ctx);CHKERRQ(ierr);
309*16e386b8SStefano Zampini   ierr = KSPDestroy(&pcbddc->deluxe_ctx->change);CHKERRQ(ierr);
31034a97f8cSStefano Zampini   ierr = PetscFree(pcbddc->deluxe_ctx);CHKERRQ(ierr);
31134a97f8cSStefano Zampini   PetscFunctionReturn(0);
31234a97f8cSStefano Zampini }
31334a97f8cSStefano Zampini 
31434a97f8cSStefano Zampini #undef __FUNCT__
31534a97f8cSStefano Zampini #define __FUNCT__ "PCBDDCScalingReset_Deluxe_Solvers"
31634a97f8cSStefano Zampini static PetscErrorCode PCBDDCScalingReset_Deluxe_Solvers(PCBDDCDeluxeScaling deluxe_ctx)
31734a97f8cSStefano Zampini {
31834a97f8cSStefano Zampini   PetscErrorCode      ierr;
31934a97f8cSStefano Zampini 
32034a97f8cSStefano Zampini   PetscFunctionBegin;
32134a97f8cSStefano Zampini   ierr = PetscFree(deluxe_ctx->idx_simple_B);CHKERRQ(ierr);
32234a97f8cSStefano Zampini   deluxe_ctx->n_simple = 0;
32334a97f8cSStefano Zampini   ierr = VecScatterDestroy(&deluxe_ctx->seq_scctx);CHKERRQ(ierr);
32434a97f8cSStefano Zampini   ierr = VecDestroy(&deluxe_ctx->seq_work1);CHKERRQ(ierr);
32534a97f8cSStefano Zampini   ierr = VecDestroy(&deluxe_ctx->seq_work2);CHKERRQ(ierr);
3265a95e1ceSStefano Zampini   ierr = MatDestroy(&deluxe_ctx->seq_mat);CHKERRQ(ierr);
327ac632422SStefano Zampini   ierr = KSPDestroy(&deluxe_ctx->seq_ksp);CHKERRQ(ierr);
32834a97f8cSStefano Zampini   PetscFunctionReturn(0);
32934a97f8cSStefano Zampini }
33034a97f8cSStefano Zampini 
33134a97f8cSStefano Zampini #undef __FUNCT__
33234a97f8cSStefano Zampini #define __FUNCT__ "PCBDDCScalingSetUp_Deluxe"
33334a97f8cSStefano Zampini static PetscErrorCode PCBDDCScalingSetUp_Deluxe(PC pc)
33434a97f8cSStefano Zampini {
33534a97f8cSStefano Zampini   PC_IS               *pcis=(PC_IS*)pc->data;
33634a97f8cSStefano Zampini   PC_BDDC             *pcbddc=(PC_BDDC*)pc->data;
33734a97f8cSStefano Zampini   PCBDDCDeluxeScaling deluxe_ctx=pcbddc->deluxe_ctx;
338b96c3477SStefano Zampini   PCBDDCSubSchurs     sub_schurs=pcbddc->sub_schurs;
33934a97f8cSStefano Zampini   PetscErrorCode      ierr;
34034a97f8cSStefano Zampini 
34134a97f8cSStefano Zampini   PetscFunctionBegin;
342b96c3477SStefano Zampini   /* (TODO: reuse) throw away the solvers */
34334a97f8cSStefano Zampini   ierr = PCBDDCScalingReset_Deluxe_Solvers(deluxe_ctx);CHKERRQ(ierr);
34434a97f8cSStefano Zampini 
345b1b3d7a2SStefano Zampini   /* Compute data structures to solve sequential problems */
3465a95e1ceSStefano Zampini   ierr = PCBDDCScalingSetUp_Deluxe_Private(pc);CHKERRQ(ierr);
3475db18549SStefano Zampini 
348b1b3d7a2SStefano Zampini   /* diagonal scaling on interface dofs not contained in cc */
349d62866d3SStefano Zampini   if (sub_schurs->is_vertices || sub_schurs->is_dir) {
3501b968477SStefano Zampini     PetscInt n_com,n_dir;
3511b968477SStefano Zampini     n_com = 0;
352d62866d3SStefano Zampini     if (sub_schurs->is_vertices) {
353d62866d3SStefano Zampini       ierr = ISGetLocalSize(sub_schurs->is_vertices,&n_com);CHKERRQ(ierr);
3541b968477SStefano Zampini     }
3551b968477SStefano Zampini     n_dir = 0;
356d62866d3SStefano Zampini     if (sub_schurs->is_dir) {
357d62866d3SStefano Zampini       ierr = ISGetLocalSize(sub_schurs->is_dir,&n_dir);CHKERRQ(ierr);
3581b968477SStefano Zampini     }
3591b968477SStefano Zampini     deluxe_ctx->n_simple = n_dir + n_com;
3601b968477SStefano Zampini     ierr = PetscMalloc1(deluxe_ctx->n_simple,&deluxe_ctx->idx_simple_B);CHKERRQ(ierr);
361d62866d3SStefano Zampini     if (sub_schurs->is_vertices) {
3629bb4a8caSStefano Zampini       PetscInt       nmap;
3639bb4a8caSStefano Zampini       const PetscInt *idxs;
3641b968477SStefano Zampini 
365d62866d3SStefano Zampini       ierr = ISGetIndices(sub_schurs->is_vertices,&idxs);CHKERRQ(ierr);
3661b968477SStefano Zampini       ierr = ISGlobalToLocalMappingApply(pcis->BtoNmap,IS_GTOLM_DROP,n_com,idxs,&nmap,deluxe_ctx->idx_simple_B);CHKERRQ(ierr);
3671b968477SStefano Zampini       if (nmap != n_com) {
368d62866d3SStefano Zampini         SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Error when mapping simply scaled dofs (is_vertices)! %d != %d",nmap,n_com);
3699bb4a8caSStefano Zampini       }
370d62866d3SStefano Zampini       ierr = ISRestoreIndices(sub_schurs->is_vertices,&idxs);CHKERRQ(ierr);
3711b968477SStefano Zampini     }
372d62866d3SStefano Zampini     if (sub_schurs->is_dir) {
3731b968477SStefano Zampini       PetscInt       nmap;
3741b968477SStefano Zampini       const PetscInt *idxs;
3751b968477SStefano Zampini 
376d62866d3SStefano Zampini       ierr = ISGetIndices(sub_schurs->is_dir,&idxs);CHKERRQ(ierr);
3771b968477SStefano Zampini       ierr = ISGlobalToLocalMappingApply(pcis->BtoNmap,IS_GTOLM_DROP,n_dir,idxs,&nmap,deluxe_ctx->idx_simple_B+n_com);CHKERRQ(ierr);
378d62866d3SStefano Zampini       if (nmap != n_dir) {
379d62866d3SStefano Zampini         SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Error when mapping simply scaled dofs (sub_schurs->is_dir)! %d != %d",nmap,n_dir);
380d62866d3SStefano Zampini       }
381d62866d3SStefano Zampini       ierr = ISRestoreIndices(sub_schurs->is_dir,&idxs);CHKERRQ(ierr);
3821b968477SStefano Zampini     }
3831b968477SStefano Zampini     ierr = PetscSortInt(deluxe_ctx->n_simple,deluxe_ctx->idx_simple_B);CHKERRQ(ierr);
3849bb4a8caSStefano Zampini   } else {
385b1b3d7a2SStefano Zampini     deluxe_ctx->n_simple = 0;
3869bb4a8caSStefano Zampini     deluxe_ctx->idx_simple_B = 0;
387b1b3d7a2SStefano Zampini   }
38834a97f8cSStefano Zampini   PetscFunctionReturn(0);
38934a97f8cSStefano Zampini }
39034a97f8cSStefano Zampini 
39134a97f8cSStefano Zampini #undef __FUNCT__
3925a95e1ceSStefano Zampini #define __FUNCT__ "PCBDDCScalingSetUp_Deluxe_Private"
3935a95e1ceSStefano Zampini static PetscErrorCode PCBDDCScalingSetUp_Deluxe_Private(PC pc)
3945db18549SStefano Zampini {
3955db18549SStefano Zampini   PC_BDDC                *pcbddc=(PC_BDDC*)pc->data;
3965db18549SStefano Zampini   PCBDDCDeluxeScaling    deluxe_ctx=pcbddc->deluxe_ctx;
397b96c3477SStefano Zampini   PCBDDCSubSchurs        sub_schurs = pcbddc->sub_schurs;
3985db18549SStefano Zampini   PetscErrorCode         ierr;
3995db18549SStefano Zampini 
4005db18549SStefano Zampini   PetscFunctionBegin;
4015a95e1ceSStefano Zampini   if (!sub_schurs->n_subs) {
4029221af80SStefano Zampini     PetscFunctionReturn(0);
4039221af80SStefano Zampini   }
4049221af80SStefano Zampini 
4055db18549SStefano Zampini   /* Create work vectors for sequential part of deluxe */
406ac632422SStefano Zampini   ierr = MatCreateVecs(sub_schurs->S_Ej_all,&deluxe_ctx->seq_work1,&deluxe_ctx->seq_work2);CHKERRQ(ierr);
4075db18549SStefano Zampini 
4085db18549SStefano Zampini   /* Compute deluxe sequential scatter */
409df4d28bfSStefano Zampini   if (sub_schurs->reuse_solver) {
410df4d28bfSStefano Zampini     PCBDDCReuseSolvers reuse_solver = sub_schurs->reuse_solver;
411925dfd53SStefano Zampini 
412df4d28bfSStefano Zampini     if (!reuse_solver->has_vertices && !sub_schurs->is_dir) {
413df4d28bfSStefano Zampini       ierr = PetscObjectReference((PetscObject)reuse_solver->correction_scatter_B);CHKERRQ(ierr);
414df4d28bfSStefano Zampini       deluxe_ctx->seq_scctx = reuse_solver->correction_scatter_B;
415d62866d3SStefano Zampini     } else {
4165db18549SStefano Zampini       ierr = VecScatterCreate(pcbddc->work_scaling,sub_schurs->is_Ej_all,deluxe_ctx->seq_work1,NULL,&deluxe_ctx->seq_scctx);CHKERRQ(ierr);
417d62866d3SStefano Zampini     }
418925dfd53SStefano Zampini   } else {
419925dfd53SStefano Zampini     ierr = VecScatterCreate(pcbddc->work_scaling,sub_schurs->is_Ej_all,deluxe_ctx->seq_work1,NULL,&deluxe_ctx->seq_scctx);CHKERRQ(ierr);
420925dfd53SStefano Zampini   }
4215db18549SStefano Zampini 
4225a95e1ceSStefano Zampini   /* Create Mat object for deluxe scaling */
4235a95e1ceSStefano Zampini   ierr = PetscObjectReference((PetscObject)sub_schurs->S_Ej_all);CHKERRQ(ierr);
4245a95e1ceSStefano Zampini   deluxe_ctx->seq_mat = sub_schurs->S_Ej_all;
425ac632422SStefano Zampini   if (sub_schurs->sum_S_Ej_all) { /* if this matrix is present, then we need to create the KSP object to invert it */
426ac632422SStefano Zampini     PC               pc_temp;
427ac632422SStefano Zampini     MatSolverPackage solver=NULL;
428ac632422SStefano Zampini     char             ksp_prefix[256];
429ac632422SStefano Zampini     size_t           len;
430ac632422SStefano Zampini 
431ac632422SStefano Zampini     ierr = KSPCreate(PETSC_COMM_SELF,&deluxe_ctx->seq_ksp);CHKERRQ(ierr);
432ac632422SStefano Zampini     ierr = KSPSetOperators(deluxe_ctx->seq_ksp,sub_schurs->sum_S_Ej_all,sub_schurs->sum_S_Ej_all);CHKERRQ(ierr);
433ac632422SStefano Zampini     ierr = KSPSetType(deluxe_ctx->seq_ksp,KSPPREONLY);CHKERRQ(ierr);
434ac632422SStefano Zampini     ierr = KSPGetPC(deluxe_ctx->seq_ksp,&pc_temp);CHKERRQ(ierr);
435ac632422SStefano Zampini     ierr = PCSetType(pc_temp,PCLU);CHKERRQ(ierr);
436ac632422SStefano Zampini     ierr = KSPGetPC(pcbddc->ksp_D,&pc_temp);CHKERRQ(ierr);
437ac632422SStefano Zampini     ierr = PCFactorGetMatSolverPackage(pc_temp,(const MatSolverPackage*)&solver);CHKERRQ(ierr);
438ac632422SStefano Zampini     if (solver) {
439ac632422SStefano Zampini       PC     new_pc;
440ac632422SStefano Zampini       PCType type;
441ac632422SStefano Zampini 
442ac632422SStefano Zampini       ierr = PCGetType(pc_temp,&type);CHKERRQ(ierr);
443ac632422SStefano Zampini       ierr = KSPGetPC(deluxe_ctx->seq_ksp,&new_pc);CHKERRQ(ierr);
444ac632422SStefano Zampini       ierr = PCSetType(new_pc,type);CHKERRQ(ierr);
445ac632422SStefano Zampini       ierr = PCFactorSetMatSolverPackage(new_pc,solver);CHKERRQ(ierr);
446ac632422SStefano Zampini     }
447ac632422SStefano Zampini     ierr = PetscStrlen(((PetscObject)(pcbddc->ksp_D))->prefix,&len);CHKERRQ(ierr);
448ac632422SStefano Zampini     len -= 10; /* remove "dirichlet_" */
449ac632422SStefano Zampini     ierr = PetscStrncpy(ksp_prefix,((PetscObject)(pcbddc->ksp_D))->prefix,len+1);CHKERRQ(ierr);
450ac632422SStefano Zampini     ierr = PetscStrcat(ksp_prefix,"deluxe_");CHKERRQ(ierr);
451ac632422SStefano Zampini     ierr = KSPSetOptionsPrefix(deluxe_ctx->seq_ksp,ksp_prefix);CHKERRQ(ierr);
452ac632422SStefano Zampini     ierr = KSPSetFromOptions(deluxe_ctx->seq_ksp);CHKERRQ(ierr);
453ac632422SStefano Zampini     ierr = KSPSetUp(deluxe_ctx->seq_ksp);CHKERRQ(ierr);
454ac632422SStefano Zampini   }
455*16e386b8SStefano Zampini   if (sub_schurs->change && !deluxe_ctx->change) {
456*16e386b8SStefano Zampini     PC               pc_temp;
457*16e386b8SStefano Zampini     char             ksp_prefix[256];
458*16e386b8SStefano Zampini     size_t           len;
459*16e386b8SStefano Zampini 
460*16e386b8SStefano Zampini     ierr = PetscObjectReference((PetscObject)sub_schurs->change);CHKERRQ(ierr);
461*16e386b8SStefano Zampini     deluxe_ctx->change = sub_schurs->change;
462*16e386b8SStefano Zampini     ierr = KSPSetType(deluxe_ctx->change,KSPPREONLY);CHKERRQ(ierr);
463*16e386b8SStefano Zampini     ierr = KSPGetPC(deluxe_ctx->change,&pc_temp);CHKERRQ(ierr);
464*16e386b8SStefano Zampini     ierr = PCSetType(pc_temp,PCLU);CHKERRQ(ierr);
465*16e386b8SStefano Zampini     ierr = PetscStrlen(((PetscObject)(pcbddc->ksp_D))->prefix,&len);CHKERRQ(ierr);
466*16e386b8SStefano Zampini     len -= 10; /* remove "dirichlet_" */
467*16e386b8SStefano Zampini     ierr = PetscStrncpy(ksp_prefix,((PetscObject)(pcbddc->ksp_D))->prefix,len+1);CHKERRQ(ierr);
468*16e386b8SStefano Zampini     ierr = PetscStrcat(ksp_prefix,"deluxe_change_");CHKERRQ(ierr);
469*16e386b8SStefano Zampini     ierr = KSPSetOptionsPrefix(deluxe_ctx->change,ksp_prefix);CHKERRQ(ierr);
470*16e386b8SStefano Zampini     ierr = KSPSetFromOptions(deluxe_ctx->change);CHKERRQ(ierr);
471*16e386b8SStefano Zampini     ierr = KSPSetUp(deluxe_ctx->change);CHKERRQ(ierr);
472*16e386b8SStefano Zampini   }
4735db18549SStefano Zampini   PetscFunctionReturn(0);
4745db18549SStefano Zampini }
475