xref: /petsc/src/ksp/pc/impls/bddc/bddcscalingbasic.c (revision 925dfd53207892095bbd12b3cecf5c93aa26dcc1)
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);
585a95e1ceSStefano Zampini     ierr = MatMultTranspose(deluxe_ctx->seq_mat,deluxe_ctx->seq_work1,deluxe_ctx->seq_work2);CHKERRQ(ierr);
59ac632422SStefano Zampini     if (deluxe_ctx->seq_ksp) {
60ac632422SStefano Zampini       ierr = KSPSolveTranspose(deluxe_ctx->seq_ksp,deluxe_ctx->seq_work2,deluxe_ctx->seq_work2);CHKERRQ(ierr);
61ac632422SStefano Zampini     }
625a95e1ceSStefano Zampini     ierr = VecScatterBegin(deluxe_ctx->seq_scctx,deluxe_ctx->seq_work2,pcbddc->work_scaling,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
635a95e1ceSStefano Zampini     ierr = VecScatterEnd(deluxe_ctx->seq_scctx,deluxe_ctx->seq_work2,pcbddc->work_scaling,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
64674ae819SStefano Zampini   }
65674ae819SStefano Zampini   /* put local boundary part in global vector */
66674ae819SStefano Zampini   ierr = VecScatterBegin(pcis->global_to_B,pcbddc->work_scaling,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
67674ae819SStefano Zampini   ierr = VecScatterEnd(pcis->global_to_B,pcbddc->work_scaling,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
68674ae819SStefano Zampini   PetscFunctionReturn(0);
69674ae819SStefano Zampini }
70674ae819SStefano Zampini 
71674ae819SStefano Zampini #undef __FUNCT__
72674ae819SStefano Zampini #define __FUNCT__ "PCBDDCScalingExtension"
73674ae819SStefano Zampini PetscErrorCode PCBDDCScalingExtension(PC pc, Vec local_interface_vector, Vec global_vector)
74674ae819SStefano Zampini {
75674ae819SStefano Zampini   PC_BDDC *pcbddc=(PC_BDDC*)pc->data;
76674ae819SStefano Zampini   PetscErrorCode ierr;
77674ae819SStefano Zampini 
78674ae819SStefano Zampini   PetscFunctionBegin;
79674ae819SStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
80674ae819SStefano Zampini   PetscValidHeaderSpecific(local_interface_vector,VEC_CLASSID,2);
81674ae819SStefano Zampini   PetscValidHeaderSpecific(global_vector,VEC_CLASSID,3);
82674ae819SStefano Zampini   if (local_interface_vector == pcbddc->work_scaling) {
83a07ea27aSStefano Zampini     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Local vector cannot be pcbddc->work_scaling!\n");
84674ae819SStefano Zampini   }
85674ae819SStefano Zampini   ierr = PetscTryMethod(pc,"PCBDDCScalingExtension_C",(PC,Vec,Vec),(pc,local_interface_vector,global_vector));CHKERRQ(ierr);
86674ae819SStefano Zampini   PetscFunctionReturn(0);
87674ae819SStefano Zampini }
88674ae819SStefano Zampini 
89674ae819SStefano Zampini #undef __FUNCT__
90674ae819SStefano Zampini #define __FUNCT__ "PCBDDCScalingRestriction_Basic"
91674ae819SStefano Zampini static PetscErrorCode PCBDDCScalingRestriction_Basic(PC pc, Vec global_vector, Vec local_interface_vector)
92674ae819SStefano Zampini {
93674ae819SStefano Zampini   PetscErrorCode ierr;
94674ae819SStefano Zampini   PC_IS* pcis = (PC_IS*)pc->data;
95674ae819SStefano Zampini 
96674ae819SStefano Zampini   PetscFunctionBegin;
97674ae819SStefano Zampini   ierr = VecScatterBegin(pcis->global_to_B,global_vector,local_interface_vector,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
98674ae819SStefano Zampini   ierr = VecScatterEnd(pcis->global_to_B,global_vector,local_interface_vector,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
99674ae819SStefano Zampini   /* Apply partition of unity */
100674ae819SStefano Zampini   ierr = VecPointwiseMult(local_interface_vector,pcis->D,local_interface_vector);CHKERRQ(ierr);
101674ae819SStefano Zampini   PetscFunctionReturn(0);
102674ae819SStefano Zampini }
103674ae819SStefano Zampini 
104674ae819SStefano Zampini #undef __FUNCT__
105674ae819SStefano Zampini #define __FUNCT__ "PCBDDCScalingRestriction_Deluxe"
106674ae819SStefano Zampini static PetscErrorCode PCBDDCScalingRestriction_Deluxe(PC pc, Vec x, Vec y)
107674ae819SStefano Zampini {
108674ae819SStefano Zampini   PC_IS*              pcis=(PC_IS*)pc->data;
109674ae819SStefano Zampini   PC_BDDC*            pcbddc=(PC_BDDC*)pc->data;
110674ae819SStefano Zampini   PCBDDCDeluxeScaling deluxe_ctx = pcbddc->deluxe_ctx;
111674ae819SStefano Zampini   PetscErrorCode      ierr;
112674ae819SStefano Zampini 
113674ae819SStefano Zampini   PetscFunctionBegin;
114674ae819SStefano Zampini   /* get local boundary part of global vector */
115674ae819SStefano Zampini   ierr = VecScatterBegin(pcis->global_to_B,x,y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
116674ae819SStefano Zampini   ierr = VecScatterEnd(pcis->global_to_B,x,y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1175a95e1ceSStefano Zampini   if (deluxe_ctx->n_simple) { /* scale deluxe vertices using diagonal scaling */
1185a95e1ceSStefano Zampini     PetscInt          i;
1192b095fd8SStefano Zampini     PetscScalar       *array_y;
1202b095fd8SStefano Zampini     const PetscScalar *array_D;
121674ae819SStefano Zampini     ierr = VecGetArray(y,&array_y);CHKERRQ(ierr);
1222b095fd8SStefano Zampini     ierr = VecGetArrayRead(pcis->D,&array_D);CHKERRQ(ierr);
123674ae819SStefano Zampini     for (i=0;i<deluxe_ctx->n_simple;i++) {
124674ae819SStefano Zampini       array_y[deluxe_ctx->idx_simple_B[i]] *= array_D[deluxe_ctx->idx_simple_B[i]];
125674ae819SStefano Zampini     }
1262b095fd8SStefano Zampini     ierr = VecRestoreArrayRead(pcis->D,&array_D);CHKERRQ(ierr);
127674ae819SStefano Zampini     ierr = VecRestoreArray(y,&array_y);CHKERRQ(ierr);
12834a97f8cSStefano Zampini   }
12934a97f8cSStefano Zampini   /* sequential part : all problems and Schur applications collapsed into a single matrix vector multiplication and ksp solution */
1305a95e1ceSStefano Zampini   if (deluxe_ctx->seq_mat) {
131674ae819SStefano Zampini     ierr = VecScatterBegin(deluxe_ctx->seq_scctx,y,deluxe_ctx->seq_work1,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
132674ae819SStefano Zampini     ierr = VecScatterEnd(deluxe_ctx->seq_scctx,y,deluxe_ctx->seq_work1,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
133ac632422SStefano Zampini     if (deluxe_ctx->seq_ksp) {
134ac632422SStefano Zampini       ierr = KSPSolve(deluxe_ctx->seq_ksp,deluxe_ctx->seq_work1,deluxe_ctx->seq_work1);CHKERRQ(ierr);
135ac632422SStefano Zampini     }
1365a95e1ceSStefano Zampini     ierr = MatMult(deluxe_ctx->seq_mat,deluxe_ctx->seq_work1,deluxe_ctx->seq_work2);CHKERRQ(ierr);
1375a95e1ceSStefano Zampini     ierr = VecScatterBegin(deluxe_ctx->seq_scctx,deluxe_ctx->seq_work2,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1385a95e1ceSStefano Zampini     ierr = VecScatterEnd(deluxe_ctx->seq_scctx,deluxe_ctx->seq_work2,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
139674ae819SStefano Zampini   }
140674ae819SStefano Zampini   PetscFunctionReturn(0);
141674ae819SStefano Zampini }
142674ae819SStefano Zampini 
143674ae819SStefano Zampini #undef __FUNCT__
144674ae819SStefano Zampini #define __FUNCT__ "PCBDDCScalingRestriction"
145674ae819SStefano Zampini PetscErrorCode PCBDDCScalingRestriction(PC pc, Vec global_vector, Vec local_interface_vector)
146674ae819SStefano Zampini {
147674ae819SStefano Zampini   PC_BDDC        *pcbddc=(PC_BDDC*)pc->data;
148674ae819SStefano Zampini   PetscErrorCode ierr;
149674ae819SStefano Zampini 
150674ae819SStefano Zampini   PetscFunctionBegin;
151674ae819SStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
152674ae819SStefano Zampini   PetscValidHeaderSpecific(global_vector,VEC_CLASSID,2);
153674ae819SStefano Zampini   PetscValidHeaderSpecific(local_interface_vector,VEC_CLASSID,3);
154674ae819SStefano Zampini   if (local_interface_vector == pcbddc->work_scaling) {
155a07ea27aSStefano Zampini     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Local vector cannot be pcbddc->work_scaling!\n");
156674ae819SStefano Zampini   }
157674ae819SStefano Zampini   ierr = PetscTryMethod(pc,"PCBDDCScalingRestriction_C",(PC,Vec,Vec),(pc,global_vector,local_interface_vector));CHKERRQ(ierr);
158674ae819SStefano Zampini   PetscFunctionReturn(0);
159674ae819SStefano Zampini }
160674ae819SStefano Zampini 
161674ae819SStefano Zampini #undef __FUNCT__
162674ae819SStefano Zampini #define __FUNCT__ "PCBDDCScalingSetUp"
163674ae819SStefano Zampini PetscErrorCode PCBDDCScalingSetUp(PC pc)
164674ae819SStefano Zampini {
165674ae819SStefano Zampini   PC_IS*         pcis=(PC_IS*)pc->data;
166674ae819SStefano Zampini   PC_BDDC*       pcbddc=(PC_BDDC*)pc->data;
167674ae819SStefano Zampini   PetscErrorCode ierr;
168674ae819SStefano Zampini 
169674ae819SStefano Zampini   PetscFunctionBegin;
170674ae819SStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
171674ae819SStefano Zampini   /* create work vector for the operator */
17234a97f8cSStefano Zampini   ierr = VecDestroy(&pcbddc->work_scaling);CHKERRQ(ierr);
173674ae819SStefano Zampini   ierr = VecDuplicate(pcis->vec1_B,&pcbddc->work_scaling);CHKERRQ(ierr);
17434a97f8cSStefano Zampini   /* always rebuild pcis->D */
17528d874f6SStefano Zampini   if (pcis->use_stiffness_scaling) {
176674ae819SStefano Zampini     ierr = MatGetDiagonal(pcbddc->local_mat,pcis->vec1_N);CHKERRQ(ierr);
177674ae819SStefano Zampini     ierr = VecScatterBegin(pcis->N_to_B,pcis->vec1_N,pcis->D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
178674ae819SStefano Zampini     ierr = VecScatterEnd(pcis->N_to_B,pcis->vec1_N,pcis->D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
179674ae819SStefano Zampini   }
180674ae819SStefano Zampini   ierr = VecCopy(pcis->D,pcis->vec1_B);CHKERRQ(ierr);
181674ae819SStefano Zampini   ierr = VecSet(pcis->vec1_global,0.0);CHKERRQ(ierr);
182674ae819SStefano Zampini   ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_B,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
183674ae819SStefano Zampini   ierr = VecScatterEnd(pcis->global_to_B,pcis->vec1_B,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
184674ae819SStefano Zampini   ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_global,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
185674ae819SStefano Zampini   ierr = VecScatterEnd(pcis->global_to_B,pcis->vec1_global,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
186674ae819SStefano Zampini   ierr = VecPointwiseDivide(pcis->D,pcis->D,pcis->vec1_B);CHKERRQ(ierr);
187674ae819SStefano Zampini   /* now setup */
188681e7c04SStefano Zampini   if (pcbddc->use_deluxe_scaling) {
18934a97f8cSStefano Zampini     if (!pcbddc->deluxe_ctx) {
19034a97f8cSStefano Zampini       ierr = PCBDDCScalingCreate_Deluxe(pc);CHKERRQ(ierr);
19134a97f8cSStefano Zampini     }
19234a97f8cSStefano Zampini     ierr = PCBDDCScalingSetUp_Deluxe(pc);CHKERRQ(ierr);
193674ae819SStefano Zampini     ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCScalingRestriction_C",PCBDDCScalingRestriction_Deluxe);CHKERRQ(ierr);
194674ae819SStefano Zampini     ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCScalingExtension_C",PCBDDCScalingExtension_Deluxe);CHKERRQ(ierr);
195674ae819SStefano Zampini   } else {
196674ae819SStefano Zampini     ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCScalingRestriction_C",PCBDDCScalingRestriction_Basic);CHKERRQ(ierr);
197674ae819SStefano Zampini     ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCScalingExtension_C",PCBDDCScalingExtension_Basic);CHKERRQ(ierr);
198674ae819SStefano Zampini   }
19934a97f8cSStefano Zampini 
200674ae819SStefano Zampini   /* test */
201674ae819SStefano Zampini   if (pcbddc->dbg_flag) {
20234a97f8cSStefano Zampini     Vec         vec2_global;
203674ae819SStefano Zampini     PetscViewer viewer=pcbddc->dbg_viewer;
20434a97f8cSStefano Zampini     PetscReal   error;
205674ae819SStefano Zampini 
206674ae819SStefano Zampini     /* extension -> from local to parallel */
20734a97f8cSStefano Zampini     ierr = VecSet(pcis->vec1_global,0.0);CHKERRQ(ierr);
20834a97f8cSStefano Zampini     ierr = VecSetRandom(pcis->vec1_B,NULL);CHKERRQ(ierr);
20934a97f8cSStefano Zampini     ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_B,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
21034a97f8cSStefano Zampini     ierr = VecScatterEnd(pcis->global_to_B,pcis->vec1_B,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
21134a97f8cSStefano Zampini     ierr = VecDuplicate(pcis->vec1_global,&vec2_global);CHKERRQ(ierr);
21234a97f8cSStefano Zampini     ierr = VecCopy(pcis->vec1_global,vec2_global);CHKERRQ(ierr);
21334a97f8cSStefano Zampini 
214674ae819SStefano Zampini     ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_global,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
215674ae819SStefano Zampini     ierr = VecScatterEnd(pcis->global_to_B,pcis->vec1_global,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
216674ae819SStefano Zampini     ierr = PCBDDCScalingExtension(pc,pcis->vec1_B,pcis->vec1_global);CHKERRQ(ierr);
21734a97f8cSStefano Zampini     ierr = VecAXPY(pcis->vec1_global,-1.0,vec2_global);CHKERRQ(ierr);
21834a97f8cSStefano Zampini     ierr = VecNorm(pcis->vec1_global,NORM_INFINITY,&error);CHKERRQ(ierr);
219674ae819SStefano Zampini     ierr = PetscViewerASCIIPrintf(viewer,"Error scaling extension %1.14e\n",error);CHKERRQ(ierr);
22034a97f8cSStefano Zampini     if (error>1.e-8 && pcbddc->dbg_flag>1) {
221674ae819SStefano Zampini       ierr = VecView(pcis->vec1_global,viewer);CHKERRQ(ierr);
222674ae819SStefano Zampini     }
22334a97f8cSStefano Zampini     ierr = VecDestroy(&vec2_global);CHKERRQ(ierr);
22434a97f8cSStefano Zampini 
225674ae819SStefano Zampini     /* restriction -> from parallel to local */
226674ae819SStefano Zampini     ierr = VecSet(pcis->vec1_global,0.0);CHKERRQ(ierr);
22734a97f8cSStefano Zampini     ierr = VecSetRandom(pcis->vec1_B,NULL);CHKERRQ(ierr);
228674ae819SStefano Zampini     ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_B,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
229674ae819SStefano Zampini     ierr = VecScatterEnd(pcis->global_to_B,pcis->vec1_B,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
23034a97f8cSStefano Zampini 
23134a97f8cSStefano Zampini     ierr = PCBDDCScalingRestriction(pc,pcis->vec1_global,pcis->vec1_B);CHKERRQ(ierr);
23234a97f8cSStefano Zampini     ierr = VecScale(pcis->vec1_B,-1.0);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 = VecNorm(pcis->vec1_global,NORM_INFINITY,&error);CHKERRQ(ierr);
23634a97f8cSStefano Zampini     ierr = PetscViewerASCIIPrintf(viewer,"Error scaling restriction %1.14e\n",error);CHKERRQ(ierr);
23734a97f8cSStefano Zampini     if (error>1.e-8 && pcbddc->dbg_flag>1) {
238674ae819SStefano Zampini       ierr = VecView(pcis->vec1_global,viewer);CHKERRQ(ierr);
239674ae819SStefano Zampini     }
240674ae819SStefano Zampini   }
241674ae819SStefano Zampini   PetscFunctionReturn(0);
242674ae819SStefano Zampini }
243674ae819SStefano Zampini 
244674ae819SStefano Zampini #undef __FUNCT__
245674ae819SStefano Zampini #define __FUNCT__ "PCBDDCScalingDestroy"
246674ae819SStefano Zampini PetscErrorCode PCBDDCScalingDestroy(PC pc)
247674ae819SStefano Zampini {
248674ae819SStefano Zampini   PC_BDDC* pcbddc=(PC_BDDC*)pc->data;
249674ae819SStefano Zampini   PetscErrorCode ierr;
250674ae819SStefano Zampini 
251674ae819SStefano Zampini   PetscFunctionBegin;
25234a97f8cSStefano Zampini   if (pcbddc->deluxe_ctx) {
25334a97f8cSStefano Zampini     ierr = PCBDDCScalingDestroy_Deluxe(pc);CHKERRQ(ierr);
254674ae819SStefano Zampini   }
255674ae819SStefano Zampini   ierr = VecDestroy(&pcbddc->work_scaling);CHKERRQ(ierr);
256674ae819SStefano Zampini   /* remove functions */
257674ae819SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCScalingRestriction_C",NULL);CHKERRQ(ierr);
258674ae819SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCScalingExtension_C",NULL);CHKERRQ(ierr);
259674ae819SStefano Zampini   PetscFunctionReturn(0);
260674ae819SStefano Zampini }
261674ae819SStefano Zampini 
26234a97f8cSStefano Zampini #undef __FUNCT__
26334a97f8cSStefano Zampini #define __FUNCT__ "PCBDDCScalingCreate_Deluxe"
26434a97f8cSStefano Zampini static PetscErrorCode PCBDDCScalingCreate_Deluxe(PC pc)
26534a97f8cSStefano Zampini {
26634a97f8cSStefano Zampini   PC_BDDC*            pcbddc=(PC_BDDC*)pc->data;
26734a97f8cSStefano Zampini   PCBDDCDeluxeScaling deluxe_ctx;
26834a97f8cSStefano Zampini   PetscErrorCode      ierr;
26934a97f8cSStefano Zampini 
27034a97f8cSStefano Zampini   PetscFunctionBegin;
27134a97f8cSStefano Zampini   ierr = PetscNew(&deluxe_ctx);CHKERRQ(ierr);
27234a97f8cSStefano Zampini   pcbddc->deluxe_ctx = deluxe_ctx;
27334a97f8cSStefano Zampini   PetscFunctionReturn(0);
27434a97f8cSStefano Zampini }
27534a97f8cSStefano Zampini 
27634a97f8cSStefano Zampini #undef __FUNCT__
27734a97f8cSStefano Zampini #define __FUNCT__ "PCBDDCScalingDestroy_Deluxe"
27834a97f8cSStefano Zampini static PetscErrorCode PCBDDCScalingDestroy_Deluxe(PC pc)
27934a97f8cSStefano Zampini {
28034a97f8cSStefano Zampini   PC_BDDC*            pcbddc=(PC_BDDC*)pc->data;
28134a97f8cSStefano Zampini   PetscErrorCode      ierr;
28234a97f8cSStefano Zampini 
28334a97f8cSStefano Zampini   PetscFunctionBegin;
28434a97f8cSStefano Zampini   ierr = PCBDDCScalingReset_Deluxe_Solvers(pcbddc->deluxe_ctx);CHKERRQ(ierr);
28534a97f8cSStefano Zampini   ierr = PetscFree(pcbddc->deluxe_ctx);CHKERRQ(ierr);
28634a97f8cSStefano Zampini   PetscFunctionReturn(0);
28734a97f8cSStefano Zampini }
28834a97f8cSStefano Zampini 
28934a97f8cSStefano Zampini #undef __FUNCT__
29034a97f8cSStefano Zampini #define __FUNCT__ "PCBDDCScalingReset_Deluxe_Solvers"
29134a97f8cSStefano Zampini static PetscErrorCode PCBDDCScalingReset_Deluxe_Solvers(PCBDDCDeluxeScaling deluxe_ctx)
29234a97f8cSStefano Zampini {
29334a97f8cSStefano Zampini   PetscErrorCode      ierr;
29434a97f8cSStefano Zampini 
29534a97f8cSStefano Zampini   PetscFunctionBegin;
29634a97f8cSStefano Zampini   ierr = PetscFree(deluxe_ctx->idx_simple_B);CHKERRQ(ierr);
29734a97f8cSStefano Zampini   deluxe_ctx->n_simple = 0;
29834a97f8cSStefano Zampini   ierr = VecScatterDestroy(&deluxe_ctx->seq_scctx);CHKERRQ(ierr);
29934a97f8cSStefano Zampini   ierr = VecDestroy(&deluxe_ctx->seq_work1);CHKERRQ(ierr);
30034a97f8cSStefano Zampini   ierr = VecDestroy(&deluxe_ctx->seq_work2);CHKERRQ(ierr);
3015a95e1ceSStefano Zampini   ierr = MatDestroy(&deluxe_ctx->seq_mat);CHKERRQ(ierr);
302ac632422SStefano Zampini   ierr = KSPDestroy(&deluxe_ctx->seq_ksp);CHKERRQ(ierr);
30334a97f8cSStefano Zampini   PetscFunctionReturn(0);
30434a97f8cSStefano Zampini }
30534a97f8cSStefano Zampini 
30634a97f8cSStefano Zampini #undef __FUNCT__
30734a97f8cSStefano Zampini #define __FUNCT__ "PCBDDCScalingSetUp_Deluxe"
30834a97f8cSStefano Zampini static PetscErrorCode PCBDDCScalingSetUp_Deluxe(PC pc)
30934a97f8cSStefano Zampini {
31034a97f8cSStefano Zampini   PC_IS               *pcis=(PC_IS*)pc->data;
31134a97f8cSStefano Zampini   PC_BDDC             *pcbddc=(PC_BDDC*)pc->data;
31234a97f8cSStefano Zampini   PCBDDCDeluxeScaling deluxe_ctx=pcbddc->deluxe_ctx;
313b96c3477SStefano Zampini   PCBDDCSubSchurs     sub_schurs=pcbddc->sub_schurs;
31434a97f8cSStefano Zampini   PetscErrorCode      ierr;
31534a97f8cSStefano Zampini 
31634a97f8cSStefano Zampini   PetscFunctionBegin;
317b96c3477SStefano Zampini   /* (TODO: reuse) throw away the solvers */
31834a97f8cSStefano Zampini   ierr = PCBDDCScalingReset_Deluxe_Solvers(deluxe_ctx);CHKERRQ(ierr);
31934a97f8cSStefano Zampini 
320b1b3d7a2SStefano Zampini   /* Compute data structures to solve sequential problems */
3215a95e1ceSStefano Zampini   ierr = PCBDDCScalingSetUp_Deluxe_Private(pc);CHKERRQ(ierr);
3225db18549SStefano Zampini 
323b1b3d7a2SStefano Zampini   /* diagonal scaling on interface dofs not contained in cc */
324d62866d3SStefano Zampini   if (sub_schurs->is_vertices || sub_schurs->is_dir) {
3251b968477SStefano Zampini     PetscInt n_com,n_dir;
3261b968477SStefano Zampini     n_com = 0;
327d62866d3SStefano Zampini     if (sub_schurs->is_vertices) {
328d62866d3SStefano Zampini       ierr = ISGetLocalSize(sub_schurs->is_vertices,&n_com);CHKERRQ(ierr);
3291b968477SStefano Zampini     }
3301b968477SStefano Zampini     n_dir = 0;
331d62866d3SStefano Zampini     if (sub_schurs->is_dir) {
332d62866d3SStefano Zampini       ierr = ISGetLocalSize(sub_schurs->is_dir,&n_dir);CHKERRQ(ierr);
3331b968477SStefano Zampini     }
3341b968477SStefano Zampini     deluxe_ctx->n_simple = n_dir + n_com;
3351b968477SStefano Zampini     ierr = PetscMalloc1(deluxe_ctx->n_simple,&deluxe_ctx->idx_simple_B);CHKERRQ(ierr);
336d62866d3SStefano Zampini     if (sub_schurs->is_vertices) {
3379bb4a8caSStefano Zampini       PetscInt       nmap;
3389bb4a8caSStefano Zampini       const PetscInt *idxs;
3391b968477SStefano Zampini 
340d62866d3SStefano Zampini       ierr = ISGetIndices(sub_schurs->is_vertices,&idxs);CHKERRQ(ierr);
3411b968477SStefano Zampini       ierr = ISGlobalToLocalMappingApply(pcis->BtoNmap,IS_GTOLM_DROP,n_com,idxs,&nmap,deluxe_ctx->idx_simple_B);CHKERRQ(ierr);
3421b968477SStefano Zampini       if (nmap != n_com) {
343d62866d3SStefano Zampini         SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Error when mapping simply scaled dofs (is_vertices)! %d != %d",nmap,n_com);
3449bb4a8caSStefano Zampini       }
345d62866d3SStefano Zampini       ierr = ISRestoreIndices(sub_schurs->is_vertices,&idxs);CHKERRQ(ierr);
3461b968477SStefano Zampini     }
347d62866d3SStefano Zampini     if (sub_schurs->is_dir) {
3481b968477SStefano Zampini       PetscInt       nmap;
3491b968477SStefano Zampini       const PetscInt *idxs;
3501b968477SStefano Zampini 
351d62866d3SStefano Zampini       ierr = ISGetIndices(sub_schurs->is_dir,&idxs);CHKERRQ(ierr);
3521b968477SStefano Zampini       ierr = ISGlobalToLocalMappingApply(pcis->BtoNmap,IS_GTOLM_DROP,n_dir,idxs,&nmap,deluxe_ctx->idx_simple_B+n_com);CHKERRQ(ierr);
353d62866d3SStefano Zampini       if (nmap != n_dir) {
354d62866d3SStefano Zampini         SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Error when mapping simply scaled dofs (sub_schurs->is_dir)! %d != %d",nmap,n_dir);
355d62866d3SStefano Zampini       }
356d62866d3SStefano Zampini       ierr = ISRestoreIndices(sub_schurs->is_dir,&idxs);CHKERRQ(ierr);
3571b968477SStefano Zampini     }
3581b968477SStefano Zampini     ierr = PetscSortInt(deluxe_ctx->n_simple,deluxe_ctx->idx_simple_B);CHKERRQ(ierr);
3599bb4a8caSStefano Zampini   } else {
360b1b3d7a2SStefano Zampini     deluxe_ctx->n_simple = 0;
3619bb4a8caSStefano Zampini     deluxe_ctx->idx_simple_B = 0;
362b1b3d7a2SStefano Zampini   }
36334a97f8cSStefano Zampini   PetscFunctionReturn(0);
36434a97f8cSStefano Zampini }
36534a97f8cSStefano Zampini 
36634a97f8cSStefano Zampini #undef __FUNCT__
3675a95e1ceSStefano Zampini #define __FUNCT__ "PCBDDCScalingSetUp_Deluxe_Private"
3685a95e1ceSStefano Zampini static PetscErrorCode PCBDDCScalingSetUp_Deluxe_Private(PC pc)
3695db18549SStefano Zampini {
3705db18549SStefano Zampini   PC_BDDC                *pcbddc=(PC_BDDC*)pc->data;
3715db18549SStefano Zampini   PCBDDCDeluxeScaling    deluxe_ctx=pcbddc->deluxe_ctx;
372b96c3477SStefano Zampini   PCBDDCSubSchurs        sub_schurs = pcbddc->sub_schurs;
3735db18549SStefano Zampini   PetscErrorCode         ierr;
3745db18549SStefano Zampini 
3755db18549SStefano Zampini   PetscFunctionBegin;
3765a95e1ceSStefano Zampini   if (!sub_schurs->n_subs) {
3779221af80SStefano Zampini     PetscFunctionReturn(0);
3789221af80SStefano Zampini   }
3799221af80SStefano Zampini 
3805db18549SStefano Zampini   /* Create work vectors for sequential part of deluxe */
381ac632422SStefano Zampini   ierr = MatCreateVecs(sub_schurs->S_Ej_all,&deluxe_ctx->seq_work1,&deluxe_ctx->seq_work2);CHKERRQ(ierr);
3825db18549SStefano Zampini 
3835db18549SStefano Zampini   /* Compute deluxe sequential scatter */
384*925dfd53SStefano Zampini   if (sub_schurs->reuse_mumps) {
385d62866d3SStefano Zampini     PCBDDCReuseMumps reuse_mumps = sub_schurs->reuse_mumps;
386*925dfd53SStefano Zampini 
387*925dfd53SStefano Zampini     if (!reuse_mumps->has_vertices && !sub_schurs->is_dir) {
388d62866d3SStefano Zampini       ierr = PetscObjectReference((PetscObject)reuse_mumps->correction_scatter_B);CHKERRQ(ierr);
389d62866d3SStefano Zampini       deluxe_ctx->seq_scctx = reuse_mumps->correction_scatter_B;
390d62866d3SStefano Zampini     } else {
3915db18549SStefano Zampini       ierr = VecScatterCreate(pcbddc->work_scaling,sub_schurs->is_Ej_all,deluxe_ctx->seq_work1,NULL,&deluxe_ctx->seq_scctx);CHKERRQ(ierr);
392d62866d3SStefano Zampini     }
393*925dfd53SStefano Zampini   } else {
394*925dfd53SStefano Zampini     ierr = VecScatterCreate(pcbddc->work_scaling,sub_schurs->is_Ej_all,deluxe_ctx->seq_work1,NULL,&deluxe_ctx->seq_scctx);CHKERRQ(ierr);
395*925dfd53SStefano Zampini   }
3965db18549SStefano Zampini 
3975a95e1ceSStefano Zampini   /* Create Mat object for deluxe scaling */
3985a95e1ceSStefano Zampini   ierr = PetscObjectReference((PetscObject)sub_schurs->S_Ej_all);CHKERRQ(ierr);
3995a95e1ceSStefano Zampini   deluxe_ctx->seq_mat = sub_schurs->S_Ej_all;
400ac632422SStefano Zampini   if (sub_schurs->sum_S_Ej_all) { /* if this matrix is present, then we need to create the KSP object to invert it */
401ac632422SStefano Zampini     PC               pc_temp;
402ac632422SStefano Zampini     MatSolverPackage solver=NULL;
403ac632422SStefano Zampini     char             ksp_prefix[256];
404ac632422SStefano Zampini     size_t           len;
405ac632422SStefano Zampini 
406ac632422SStefano Zampini     ierr = KSPCreate(PETSC_COMM_SELF,&deluxe_ctx->seq_ksp);CHKERRQ(ierr);
407ac632422SStefano Zampini     ierr = KSPSetOperators(deluxe_ctx->seq_ksp,sub_schurs->sum_S_Ej_all,sub_schurs->sum_S_Ej_all);CHKERRQ(ierr);
408ac632422SStefano Zampini     ierr = KSPSetType(deluxe_ctx->seq_ksp,KSPPREONLY);CHKERRQ(ierr);
409ac632422SStefano Zampini     ierr = KSPGetPC(deluxe_ctx->seq_ksp,&pc_temp);CHKERRQ(ierr);
410ac632422SStefano Zampini     ierr = PCSetType(pc_temp,PCLU);CHKERRQ(ierr);
411ac632422SStefano Zampini     ierr = KSPGetPC(pcbddc->ksp_D,&pc_temp);CHKERRQ(ierr);
412ac632422SStefano Zampini     ierr = PCFactorGetMatSolverPackage(pc_temp,(const MatSolverPackage*)&solver);CHKERRQ(ierr);
413ac632422SStefano Zampini     if (solver) {
414ac632422SStefano Zampini       PC     new_pc;
415ac632422SStefano Zampini       PCType type;
416ac632422SStefano Zampini 
417ac632422SStefano Zampini       ierr = PCGetType(pc_temp,&type);CHKERRQ(ierr);
418ac632422SStefano Zampini       ierr = KSPGetPC(deluxe_ctx->seq_ksp,&new_pc);CHKERRQ(ierr);
419ac632422SStefano Zampini       ierr = PCSetType(new_pc,type);CHKERRQ(ierr);
420ac632422SStefano Zampini       ierr = PCFactorSetMatSolverPackage(new_pc,solver);CHKERRQ(ierr);
421ac632422SStefano Zampini     }
422ac632422SStefano Zampini     ierr = PetscStrlen(((PetscObject)(pcbddc->ksp_D))->prefix,&len);CHKERRQ(ierr);
423ac632422SStefano Zampini     len -= 10; /* remove "dirichlet_" */
424ac632422SStefano Zampini     ierr = PetscStrncpy(ksp_prefix,((PetscObject)(pcbddc->ksp_D))->prefix,len+1);CHKERRQ(ierr);
425ac632422SStefano Zampini     ierr = PetscStrcat(ksp_prefix,"deluxe_");CHKERRQ(ierr);
426ac632422SStefano Zampini     ierr = KSPSetOptionsPrefix(deluxe_ctx->seq_ksp,ksp_prefix);CHKERRQ(ierr);
427ac632422SStefano Zampini     ierr = KSPSetFromOptions(deluxe_ctx->seq_ksp);CHKERRQ(ierr);
428ac632422SStefano Zampini     ierr = KSPSetUp(deluxe_ctx->seq_ksp);CHKERRQ(ierr);
429ac632422SStefano Zampini   }
4305db18549SStefano Zampini   PetscFunctionReturn(0);
4315db18549SStefano Zampini }
432