xref: /petsc/src/ksp/pc/impls/bddc/bddcscalingbasic.c (revision 163d334ea1df412efbe58c2e8fce360a2afeaf87)
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);
826c4ed002SBarry Smith   if (local_interface_vector == pcbddc->work_scaling) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Local vector cannot be pcbddc->work_scaling!\n");
83*163d334eSBarry Smith   ierr = PetscUseMethod(pc,"PCBDDCScalingExtension_C",(PC,Vec,Vec),(pc,local_interface_vector,global_vector));CHKERRQ(ierr);
84674ae819SStefano Zampini   PetscFunctionReturn(0);
85674ae819SStefano Zampini }
86674ae819SStefano Zampini 
87674ae819SStefano Zampini #undef __FUNCT__
88674ae819SStefano Zampini #define __FUNCT__ "PCBDDCScalingRestriction_Basic"
89674ae819SStefano Zampini static PetscErrorCode PCBDDCScalingRestriction_Basic(PC pc, Vec global_vector, Vec local_interface_vector)
90674ae819SStefano Zampini {
91674ae819SStefano Zampini   PetscErrorCode ierr;
92674ae819SStefano Zampini   PC_IS          *pcis = (PC_IS*)pc->data;
93674ae819SStefano Zampini 
94674ae819SStefano Zampini   PetscFunctionBegin;
95674ae819SStefano Zampini   ierr = VecScatterBegin(pcis->global_to_B,global_vector,local_interface_vector,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
96674ae819SStefano Zampini   ierr = VecScatterEnd(pcis->global_to_B,global_vector,local_interface_vector,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
97674ae819SStefano Zampini   /* Apply partition of unity */
98674ae819SStefano Zampini   ierr = VecPointwiseMult(local_interface_vector,pcis->D,local_interface_vector);CHKERRQ(ierr);
99674ae819SStefano Zampini   PetscFunctionReturn(0);
100674ae819SStefano Zampini }
101674ae819SStefano Zampini 
102674ae819SStefano Zampini #undef __FUNCT__
103674ae819SStefano Zampini #define __FUNCT__ "PCBDDCScalingRestriction_Deluxe"
104674ae819SStefano Zampini static PetscErrorCode PCBDDCScalingRestriction_Deluxe(PC pc, Vec x, Vec y)
105674ae819SStefano Zampini {
106674ae819SStefano Zampini   PC_IS*              pcis=(PC_IS*)pc->data;
107674ae819SStefano Zampini   PC_BDDC*            pcbddc=(PC_BDDC*)pc->data;
108674ae819SStefano Zampini   PCBDDCDeluxeScaling deluxe_ctx = pcbddc->deluxe_ctx;
109674ae819SStefano Zampini   PetscErrorCode      ierr;
110674ae819SStefano Zampini 
111674ae819SStefano Zampini   PetscFunctionBegin;
112674ae819SStefano Zampini   /* get local boundary part of global vector */
113674ae819SStefano Zampini   ierr = VecScatterBegin(pcis->global_to_B,x,y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
114674ae819SStefano Zampini   ierr = VecScatterEnd(pcis->global_to_B,x,y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1155a95e1ceSStefano Zampini   if (deluxe_ctx->n_simple) { /* scale deluxe vertices using diagonal scaling */
1165a95e1ceSStefano Zampini     PetscInt          i;
1172b095fd8SStefano Zampini     PetscScalar       *array_y;
1182b095fd8SStefano Zampini     const PetscScalar *array_D;
119674ae819SStefano Zampini     ierr = VecGetArray(y,&array_y);CHKERRQ(ierr);
1202b095fd8SStefano Zampini     ierr = VecGetArrayRead(pcis->D,&array_D);CHKERRQ(ierr);
121674ae819SStefano Zampini     for (i=0;i<deluxe_ctx->n_simple;i++) {
122674ae819SStefano Zampini       array_y[deluxe_ctx->idx_simple_B[i]] *= array_D[deluxe_ctx->idx_simple_B[i]];
123674ae819SStefano Zampini     }
1242b095fd8SStefano Zampini     ierr = VecRestoreArrayRead(pcis->D,&array_D);CHKERRQ(ierr);
125674ae819SStefano Zampini     ierr = VecRestoreArray(y,&array_y);CHKERRQ(ierr);
12634a97f8cSStefano Zampini   }
12734a97f8cSStefano Zampini   /* sequential part : all problems and Schur applications collapsed into a single matrix vector multiplication and ksp solution */
1285a95e1ceSStefano Zampini   if (deluxe_ctx->seq_mat) {
129674ae819SStefano Zampini     ierr = VecScatterBegin(deluxe_ctx->seq_scctx,y,deluxe_ctx->seq_work1,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
130674ae819SStefano Zampini     ierr = VecScatterEnd(deluxe_ctx->seq_scctx,y,deluxe_ctx->seq_work1,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
131ac632422SStefano Zampini     if (deluxe_ctx->seq_ksp) {
132ac632422SStefano Zampini       ierr = KSPSolve(deluxe_ctx->seq_ksp,deluxe_ctx->seq_work1,deluxe_ctx->seq_work1);CHKERRQ(ierr);
133ac632422SStefano Zampini     }
1345a95e1ceSStefano Zampini     ierr = MatMult(deluxe_ctx->seq_mat,deluxe_ctx->seq_work1,deluxe_ctx->seq_work2);CHKERRQ(ierr);
1355a95e1ceSStefano Zampini     ierr = VecScatterBegin(deluxe_ctx->seq_scctx,deluxe_ctx->seq_work2,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1365a95e1ceSStefano Zampini     ierr = VecScatterEnd(deluxe_ctx->seq_scctx,deluxe_ctx->seq_work2,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
137674ae819SStefano Zampini   }
138674ae819SStefano Zampini   PetscFunctionReturn(0);
139674ae819SStefano Zampini }
140674ae819SStefano Zampini 
141674ae819SStefano Zampini #undef __FUNCT__
142674ae819SStefano Zampini #define __FUNCT__ "PCBDDCScalingRestriction"
143674ae819SStefano Zampini PetscErrorCode PCBDDCScalingRestriction(PC pc, Vec global_vector, Vec local_interface_vector)
144674ae819SStefano Zampini {
145674ae819SStefano Zampini   PC_BDDC        *pcbddc=(PC_BDDC*)pc->data;
146674ae819SStefano Zampini   PetscErrorCode ierr;
147674ae819SStefano Zampini 
148674ae819SStefano Zampini   PetscFunctionBegin;
149674ae819SStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
150674ae819SStefano Zampini   PetscValidHeaderSpecific(global_vector,VEC_CLASSID,2);
151674ae819SStefano Zampini   PetscValidHeaderSpecific(local_interface_vector,VEC_CLASSID,3);
1526c4ed002SBarry Smith   if (local_interface_vector == pcbddc->work_scaling) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Local vector cannot be pcbddc->work_scaling!\n");
153*163d334eSBarry Smith   ierr = PetscUseMethod(pc,"PCBDDCScalingRestriction_C",(PC,Vec,Vec),(pc,global_vector,local_interface_vector));CHKERRQ(ierr);
154674ae819SStefano Zampini   PetscFunctionReturn(0);
155674ae819SStefano Zampini }
156674ae819SStefano Zampini 
157674ae819SStefano Zampini #undef __FUNCT__
158674ae819SStefano Zampini #define __FUNCT__ "PCBDDCScalingSetUp"
159674ae819SStefano Zampini PetscErrorCode PCBDDCScalingSetUp(PC pc)
160674ae819SStefano Zampini {
161674ae819SStefano Zampini   PC_IS*         pcis=(PC_IS*)pc->data;
162674ae819SStefano Zampini   PC_BDDC*       pcbddc=(PC_BDDC*)pc->data;
163674ae819SStefano Zampini   PetscErrorCode ierr;
164674ae819SStefano Zampini 
165674ae819SStefano Zampini   PetscFunctionBegin;
166674ae819SStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
167674ae819SStefano Zampini   /* create work vector for the operator */
16834a97f8cSStefano Zampini   ierr = VecDestroy(&pcbddc->work_scaling);CHKERRQ(ierr);
169674ae819SStefano Zampini   ierr = VecDuplicate(pcis->vec1_B,&pcbddc->work_scaling);CHKERRQ(ierr);
17034a97f8cSStefano Zampini   /* always rebuild pcis->D */
17128d874f6SStefano Zampini   if (pcis->use_stiffness_scaling) {
172674ae819SStefano Zampini     ierr = MatGetDiagonal(pcbddc->local_mat,pcis->vec1_N);CHKERRQ(ierr);
173674ae819SStefano Zampini     ierr = VecScatterBegin(pcis->N_to_B,pcis->vec1_N,pcis->D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
174674ae819SStefano Zampini     ierr = VecScatterEnd(pcis->N_to_B,pcis->vec1_N,pcis->D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
175674ae819SStefano Zampini   }
176674ae819SStefano Zampini   ierr = VecCopy(pcis->D,pcis->vec1_B);CHKERRQ(ierr);
177674ae819SStefano Zampini   ierr = VecSet(pcis->vec1_global,0.0);CHKERRQ(ierr);
178674ae819SStefano Zampini   ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_B,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
179674ae819SStefano Zampini   ierr = VecScatterEnd(pcis->global_to_B,pcis->vec1_B,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
180674ae819SStefano Zampini   ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_global,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
181674ae819SStefano Zampini   ierr = VecScatterEnd(pcis->global_to_B,pcis->vec1_global,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
182674ae819SStefano Zampini   ierr = VecPointwiseDivide(pcis->D,pcis->D,pcis->vec1_B);CHKERRQ(ierr);
183674ae819SStefano Zampini   /* now setup */
184681e7c04SStefano Zampini   if (pcbddc->use_deluxe_scaling) {
18534a97f8cSStefano Zampini     if (!pcbddc->deluxe_ctx) {
18634a97f8cSStefano Zampini       ierr = PCBDDCScalingCreate_Deluxe(pc);CHKERRQ(ierr);
18734a97f8cSStefano Zampini     }
18834a97f8cSStefano Zampini     ierr = PCBDDCScalingSetUp_Deluxe(pc);CHKERRQ(ierr);
189674ae819SStefano Zampini     ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCScalingRestriction_C",PCBDDCScalingRestriction_Deluxe);CHKERRQ(ierr);
190674ae819SStefano Zampini     ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCScalingExtension_C",PCBDDCScalingExtension_Deluxe);CHKERRQ(ierr);
191674ae819SStefano Zampini   } else {
192674ae819SStefano Zampini     ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCScalingRestriction_C",PCBDDCScalingRestriction_Basic);CHKERRQ(ierr);
193674ae819SStefano Zampini     ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCScalingExtension_C",PCBDDCScalingExtension_Basic);CHKERRQ(ierr);
194674ae819SStefano Zampini   }
19534a97f8cSStefano Zampini 
196674ae819SStefano Zampini   /* test */
197674ae819SStefano Zampini   if (pcbddc->dbg_flag) {
19834a97f8cSStefano Zampini     Vec         vec2_global;
199674ae819SStefano Zampini     PetscViewer viewer=pcbddc->dbg_viewer;
20034a97f8cSStefano Zampini     PetscReal   error;
201674ae819SStefano Zampini 
202674ae819SStefano Zampini     /* extension -> from local to parallel */
20334a97f8cSStefano Zampini     ierr = VecSet(pcis->vec1_global,0.0);CHKERRQ(ierr);
20434a97f8cSStefano Zampini     ierr = VecSetRandom(pcis->vec1_B,NULL);CHKERRQ(ierr);
20534a97f8cSStefano Zampini     ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_B,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
20634a97f8cSStefano Zampini     ierr = VecScatterEnd(pcis->global_to_B,pcis->vec1_B,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
20734a97f8cSStefano Zampini     ierr = VecDuplicate(pcis->vec1_global,&vec2_global);CHKERRQ(ierr);
20834a97f8cSStefano Zampini     ierr = VecCopy(pcis->vec1_global,vec2_global);CHKERRQ(ierr);
20934a97f8cSStefano Zampini 
210674ae819SStefano Zampini     ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_global,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
211674ae819SStefano Zampini     ierr = VecScatterEnd(pcis->global_to_B,pcis->vec1_global,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
212674ae819SStefano Zampini     ierr = PCBDDCScalingExtension(pc,pcis->vec1_B,pcis->vec1_global);CHKERRQ(ierr);
21334a97f8cSStefano Zampini     ierr = VecAXPY(pcis->vec1_global,-1.0,vec2_global);CHKERRQ(ierr);
21434a97f8cSStefano Zampini     ierr = VecNorm(pcis->vec1_global,NORM_INFINITY,&error);CHKERRQ(ierr);
215674ae819SStefano Zampini     ierr = PetscViewerASCIIPrintf(viewer,"Error scaling extension %1.14e\n",error);CHKERRQ(ierr);
21634a97f8cSStefano Zampini     if (error>1.e-8 && pcbddc->dbg_flag>1) {
217674ae819SStefano Zampini       ierr = VecView(pcis->vec1_global,viewer);CHKERRQ(ierr);
218674ae819SStefano Zampini     }
21934a97f8cSStefano Zampini     ierr = VecDestroy(&vec2_global);CHKERRQ(ierr);
22034a97f8cSStefano Zampini 
221674ae819SStefano Zampini     /* restriction -> from parallel to local */
222674ae819SStefano Zampini     ierr = VecSet(pcis->vec1_global,0.0);CHKERRQ(ierr);
22334a97f8cSStefano Zampini     ierr = VecSetRandom(pcis->vec1_B,NULL);CHKERRQ(ierr);
224674ae819SStefano Zampini     ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_B,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
225674ae819SStefano Zampini     ierr = VecScatterEnd(pcis->global_to_B,pcis->vec1_B,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
22634a97f8cSStefano Zampini 
22734a97f8cSStefano Zampini     ierr = PCBDDCScalingRestriction(pc,pcis->vec1_global,pcis->vec1_B);CHKERRQ(ierr);
22834a97f8cSStefano Zampini     ierr = VecScale(pcis->vec1_B,-1.0);CHKERRQ(ierr);
22934a97f8cSStefano Zampini     ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_B,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
23034a97f8cSStefano Zampini     ierr = VecScatterEnd(pcis->global_to_B,pcis->vec1_B,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
23134a97f8cSStefano Zampini     ierr = VecNorm(pcis->vec1_global,NORM_INFINITY,&error);CHKERRQ(ierr);
23234a97f8cSStefano Zampini     ierr = PetscViewerASCIIPrintf(viewer,"Error scaling restriction %1.14e\n",error);CHKERRQ(ierr);
23334a97f8cSStefano Zampini     if (error>1.e-8 && pcbddc->dbg_flag>1) {
234674ae819SStefano Zampini       ierr = VecView(pcis->vec1_global,viewer);CHKERRQ(ierr);
235674ae819SStefano Zampini     }
236674ae819SStefano Zampini   }
237674ae819SStefano Zampini   PetscFunctionReturn(0);
238674ae819SStefano Zampini }
239674ae819SStefano Zampini 
240674ae819SStefano Zampini #undef __FUNCT__
241674ae819SStefano Zampini #define __FUNCT__ "PCBDDCScalingDestroy"
242674ae819SStefano Zampini PetscErrorCode PCBDDCScalingDestroy(PC pc)
243674ae819SStefano Zampini {
244674ae819SStefano Zampini   PC_BDDC* pcbddc=(PC_BDDC*)pc->data;
245674ae819SStefano Zampini   PetscErrorCode ierr;
246674ae819SStefano Zampini 
247674ae819SStefano Zampini   PetscFunctionBegin;
24834a97f8cSStefano Zampini   if (pcbddc->deluxe_ctx) {
24934a97f8cSStefano Zampini     ierr = PCBDDCScalingDestroy_Deluxe(pc);CHKERRQ(ierr);
250674ae819SStefano Zampini   }
251674ae819SStefano Zampini   ierr = VecDestroy(&pcbddc->work_scaling);CHKERRQ(ierr);
252674ae819SStefano Zampini   /* remove functions */
253674ae819SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCScalingRestriction_C",NULL);CHKERRQ(ierr);
254674ae819SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCScalingExtension_C",NULL);CHKERRQ(ierr);
255674ae819SStefano Zampini   PetscFunctionReturn(0);
256674ae819SStefano Zampini }
257674ae819SStefano Zampini 
25834a97f8cSStefano Zampini #undef __FUNCT__
25934a97f8cSStefano Zampini #define __FUNCT__ "PCBDDCScalingCreate_Deluxe"
26034a97f8cSStefano Zampini static PetscErrorCode PCBDDCScalingCreate_Deluxe(PC pc)
26134a97f8cSStefano Zampini {
26234a97f8cSStefano Zampini   PC_BDDC*            pcbddc=(PC_BDDC*)pc->data;
26334a97f8cSStefano Zampini   PCBDDCDeluxeScaling deluxe_ctx;
26434a97f8cSStefano Zampini   PetscErrorCode      ierr;
26534a97f8cSStefano Zampini 
26634a97f8cSStefano Zampini   PetscFunctionBegin;
26734a97f8cSStefano Zampini   ierr = PetscNew(&deluxe_ctx);CHKERRQ(ierr);
26834a97f8cSStefano Zampini   pcbddc->deluxe_ctx = deluxe_ctx;
26934a97f8cSStefano Zampini   PetscFunctionReturn(0);
27034a97f8cSStefano Zampini }
27134a97f8cSStefano Zampini 
27234a97f8cSStefano Zampini #undef __FUNCT__
27334a97f8cSStefano Zampini #define __FUNCT__ "PCBDDCScalingDestroy_Deluxe"
27434a97f8cSStefano Zampini static PetscErrorCode PCBDDCScalingDestroy_Deluxe(PC pc)
27534a97f8cSStefano Zampini {
27634a97f8cSStefano Zampini   PC_BDDC*            pcbddc=(PC_BDDC*)pc->data;
27734a97f8cSStefano Zampini   PetscErrorCode      ierr;
27834a97f8cSStefano Zampini 
27934a97f8cSStefano Zampini   PetscFunctionBegin;
28034a97f8cSStefano Zampini   ierr = PCBDDCScalingReset_Deluxe_Solvers(pcbddc->deluxe_ctx);CHKERRQ(ierr);
28134a97f8cSStefano Zampini   ierr = PetscFree(pcbddc->deluxe_ctx);CHKERRQ(ierr);
28234a97f8cSStefano Zampini   PetscFunctionReturn(0);
28334a97f8cSStefano Zampini }
28434a97f8cSStefano Zampini 
28534a97f8cSStefano Zampini #undef __FUNCT__
28634a97f8cSStefano Zampini #define __FUNCT__ "PCBDDCScalingReset_Deluxe_Solvers"
28734a97f8cSStefano Zampini static PetscErrorCode PCBDDCScalingReset_Deluxe_Solvers(PCBDDCDeluxeScaling deluxe_ctx)
28834a97f8cSStefano Zampini {
28934a97f8cSStefano Zampini   PetscErrorCode      ierr;
29034a97f8cSStefano Zampini 
29134a97f8cSStefano Zampini   PetscFunctionBegin;
29234a97f8cSStefano Zampini   ierr = PetscFree(deluxe_ctx->idx_simple_B);CHKERRQ(ierr);
29334a97f8cSStefano Zampini   deluxe_ctx->n_simple = 0;
29434a97f8cSStefano Zampini   ierr = VecScatterDestroy(&deluxe_ctx->seq_scctx);CHKERRQ(ierr);
29534a97f8cSStefano Zampini   ierr = VecDestroy(&deluxe_ctx->seq_work1);CHKERRQ(ierr);
29634a97f8cSStefano Zampini   ierr = VecDestroy(&deluxe_ctx->seq_work2);CHKERRQ(ierr);
2975a95e1ceSStefano Zampini   ierr = MatDestroy(&deluxe_ctx->seq_mat);CHKERRQ(ierr);
298ac632422SStefano Zampini   ierr = KSPDestroy(&deluxe_ctx->seq_ksp);CHKERRQ(ierr);
29934a97f8cSStefano Zampini   PetscFunctionReturn(0);
30034a97f8cSStefano Zampini }
30134a97f8cSStefano Zampini 
30234a97f8cSStefano Zampini #undef __FUNCT__
30334a97f8cSStefano Zampini #define __FUNCT__ "PCBDDCScalingSetUp_Deluxe"
30434a97f8cSStefano Zampini static PetscErrorCode PCBDDCScalingSetUp_Deluxe(PC pc)
30534a97f8cSStefano Zampini {
30634a97f8cSStefano Zampini   PC_IS               *pcis=(PC_IS*)pc->data;
30734a97f8cSStefano Zampini   PC_BDDC             *pcbddc=(PC_BDDC*)pc->data;
30834a97f8cSStefano Zampini   PCBDDCDeluxeScaling deluxe_ctx=pcbddc->deluxe_ctx;
309b96c3477SStefano Zampini   PCBDDCSubSchurs     sub_schurs=pcbddc->sub_schurs;
31034a97f8cSStefano Zampini   PetscErrorCode      ierr;
31134a97f8cSStefano Zampini 
31234a97f8cSStefano Zampini   PetscFunctionBegin;
313b96c3477SStefano Zampini   /* (TODO: reuse) throw away the solvers */
31434a97f8cSStefano Zampini   ierr = PCBDDCScalingReset_Deluxe_Solvers(deluxe_ctx);CHKERRQ(ierr);
31534a97f8cSStefano Zampini 
316b1b3d7a2SStefano Zampini   /* Compute data structures to solve sequential problems */
3175a95e1ceSStefano Zampini   ierr = PCBDDCScalingSetUp_Deluxe_Private(pc);CHKERRQ(ierr);
3185db18549SStefano Zampini 
319b1b3d7a2SStefano Zampini   /* diagonal scaling on interface dofs not contained in cc */
320d62866d3SStefano Zampini   if (sub_schurs->is_vertices || sub_schurs->is_dir) {
3211b968477SStefano Zampini     PetscInt n_com,n_dir;
3221b968477SStefano Zampini     n_com = 0;
323d62866d3SStefano Zampini     if (sub_schurs->is_vertices) {
324d62866d3SStefano Zampini       ierr = ISGetLocalSize(sub_schurs->is_vertices,&n_com);CHKERRQ(ierr);
3251b968477SStefano Zampini     }
3261b968477SStefano Zampini     n_dir = 0;
327d62866d3SStefano Zampini     if (sub_schurs->is_dir) {
328d62866d3SStefano Zampini       ierr = ISGetLocalSize(sub_schurs->is_dir,&n_dir);CHKERRQ(ierr);
3291b968477SStefano Zampini     }
3301b968477SStefano Zampini     deluxe_ctx->n_simple = n_dir + n_com;
3311b968477SStefano Zampini     ierr = PetscMalloc1(deluxe_ctx->n_simple,&deluxe_ctx->idx_simple_B);CHKERRQ(ierr);
332d62866d3SStefano Zampini     if (sub_schurs->is_vertices) {
3339bb4a8caSStefano Zampini       PetscInt       nmap;
3349bb4a8caSStefano Zampini       const PetscInt *idxs;
3351b968477SStefano Zampini 
336d62866d3SStefano Zampini       ierr = ISGetIndices(sub_schurs->is_vertices,&idxs);CHKERRQ(ierr);
3371b968477SStefano Zampini       ierr = ISGlobalToLocalMappingApply(pcis->BtoNmap,IS_GTOLM_DROP,n_com,idxs,&nmap,deluxe_ctx->idx_simple_B);CHKERRQ(ierr);
3386c4ed002SBarry Smith       if (nmap != n_com) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Error when mapping simply scaled dofs (is_vertices)! %D != %D",nmap,n_com);
339d62866d3SStefano Zampini       ierr = ISRestoreIndices(sub_schurs->is_vertices,&idxs);CHKERRQ(ierr);
3401b968477SStefano Zampini     }
341d62866d3SStefano Zampini     if (sub_schurs->is_dir) {
3421b968477SStefano Zampini       PetscInt       nmap;
3431b968477SStefano Zampini       const PetscInt *idxs;
3441b968477SStefano Zampini 
345d62866d3SStefano Zampini       ierr = ISGetIndices(sub_schurs->is_dir,&idxs);CHKERRQ(ierr);
3461b968477SStefano Zampini       ierr = ISGlobalToLocalMappingApply(pcis->BtoNmap,IS_GTOLM_DROP,n_dir,idxs,&nmap,deluxe_ctx->idx_simple_B+n_com);CHKERRQ(ierr);
3476c4ed002SBarry Smith       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);
348d62866d3SStefano Zampini       ierr = ISRestoreIndices(sub_schurs->is_dir,&idxs);CHKERRQ(ierr);
3491b968477SStefano Zampini     }
3501b968477SStefano Zampini     ierr = PetscSortInt(deluxe_ctx->n_simple,deluxe_ctx->idx_simple_B);CHKERRQ(ierr);
3519bb4a8caSStefano Zampini   } else {
352b1b3d7a2SStefano Zampini     deluxe_ctx->n_simple = 0;
3539bb4a8caSStefano Zampini     deluxe_ctx->idx_simple_B = 0;
354b1b3d7a2SStefano Zampini   }
35534a97f8cSStefano Zampini   PetscFunctionReturn(0);
35634a97f8cSStefano Zampini }
35734a97f8cSStefano Zampini 
35834a97f8cSStefano Zampini #undef __FUNCT__
3595a95e1ceSStefano Zampini #define __FUNCT__ "PCBDDCScalingSetUp_Deluxe_Private"
3605a95e1ceSStefano Zampini static PetscErrorCode PCBDDCScalingSetUp_Deluxe_Private(PC pc)
3615db18549SStefano Zampini {
3625db18549SStefano Zampini   PC_BDDC                *pcbddc=(PC_BDDC*)pc->data;
3635db18549SStefano Zampini   PCBDDCDeluxeScaling    deluxe_ctx=pcbddc->deluxe_ctx;
364b96c3477SStefano Zampini   PCBDDCSubSchurs        sub_schurs = pcbddc->sub_schurs;
3655db18549SStefano Zampini   PetscErrorCode         ierr;
3665db18549SStefano Zampini 
3675db18549SStefano Zampini   PetscFunctionBegin;
3685a95e1ceSStefano Zampini   if (!sub_schurs->n_subs) {
3699221af80SStefano Zampini     PetscFunctionReturn(0);
3709221af80SStefano Zampini   }
3719221af80SStefano Zampini 
3725db18549SStefano Zampini   /* Create work vectors for sequential part of deluxe */
373ac632422SStefano Zampini   ierr = MatCreateVecs(sub_schurs->S_Ej_all,&deluxe_ctx->seq_work1,&deluxe_ctx->seq_work2);CHKERRQ(ierr);
3745db18549SStefano Zampini 
3755db18549SStefano Zampini   /* Compute deluxe sequential scatter */
376d62866d3SStefano Zampini   if (sub_schurs->reuse_mumps && !sub_schurs->is_dir) {
377d62866d3SStefano Zampini     PCBDDCReuseMumps reuse_mumps = sub_schurs->reuse_mumps;
378d62866d3SStefano Zampini     ierr = PetscObjectReference((PetscObject)reuse_mumps->correction_scatter_B);CHKERRQ(ierr);
379d62866d3SStefano Zampini     deluxe_ctx->seq_scctx = reuse_mumps->correction_scatter_B;
380d62866d3SStefano Zampini   } else {
3815db18549SStefano Zampini     ierr = VecScatterCreate(pcbddc->work_scaling,sub_schurs->is_Ej_all,deluxe_ctx->seq_work1,NULL,&deluxe_ctx->seq_scctx);CHKERRQ(ierr);
382d62866d3SStefano Zampini   }
3835db18549SStefano Zampini 
3845a95e1ceSStefano Zampini   /* Create Mat object for deluxe scaling */
3855a95e1ceSStefano Zampini   ierr = PetscObjectReference((PetscObject)sub_schurs->S_Ej_all);CHKERRQ(ierr);
3865a95e1ceSStefano Zampini   deluxe_ctx->seq_mat = sub_schurs->S_Ej_all;
387ac632422SStefano Zampini   if (sub_schurs->sum_S_Ej_all) { /* if this matrix is present, then we need to create the KSP object to invert it */
388ac632422SStefano Zampini     PC               pc_temp;
389ac632422SStefano Zampini     MatSolverPackage solver=NULL;
390ac632422SStefano Zampini     char             ksp_prefix[256];
391ac632422SStefano Zampini     size_t           len;
392ac632422SStefano Zampini 
393ac632422SStefano Zampini     ierr = KSPCreate(PETSC_COMM_SELF,&deluxe_ctx->seq_ksp);CHKERRQ(ierr);
394ac632422SStefano Zampini     ierr = KSPSetOperators(deluxe_ctx->seq_ksp,sub_schurs->sum_S_Ej_all,sub_schurs->sum_S_Ej_all);CHKERRQ(ierr);
395ac632422SStefano Zampini     ierr = KSPSetType(deluxe_ctx->seq_ksp,KSPPREONLY);CHKERRQ(ierr);
396ac632422SStefano Zampini     ierr = KSPGetPC(deluxe_ctx->seq_ksp,&pc_temp);CHKERRQ(ierr);
397ac632422SStefano Zampini     ierr = PCSetType(pc_temp,PCLU);CHKERRQ(ierr);
398ac632422SStefano Zampini     ierr = KSPGetPC(pcbddc->ksp_D,&pc_temp);CHKERRQ(ierr);
399ac632422SStefano Zampini     ierr = PCFactorGetMatSolverPackage(pc_temp,(const MatSolverPackage*)&solver);CHKERRQ(ierr);
400ac632422SStefano Zampini     if (solver) {
401ac632422SStefano Zampini       PC     new_pc;
402ac632422SStefano Zampini       PCType type;
403ac632422SStefano Zampini 
404ac632422SStefano Zampini       ierr = PCGetType(pc_temp,&type);CHKERRQ(ierr);
405ac632422SStefano Zampini       ierr = KSPGetPC(deluxe_ctx->seq_ksp,&new_pc);CHKERRQ(ierr);
406ac632422SStefano Zampini       ierr = PCSetType(new_pc,type);CHKERRQ(ierr);
407ac632422SStefano Zampini       ierr = PCFactorSetMatSolverPackage(new_pc,solver);CHKERRQ(ierr);
408ac632422SStefano Zampini     }
409ac632422SStefano Zampini     ierr = PetscStrlen(((PetscObject)(pcbddc->ksp_D))->prefix,&len);CHKERRQ(ierr);
410ac632422SStefano Zampini     len -= 10; /* remove "dirichlet_" */
411ac632422SStefano Zampini     ierr = PetscStrncpy(ksp_prefix,((PetscObject)(pcbddc->ksp_D))->prefix,len+1);CHKERRQ(ierr);
412ac632422SStefano Zampini     ierr = PetscStrcat(ksp_prefix,"deluxe_");CHKERRQ(ierr);
413ac632422SStefano Zampini     ierr = KSPSetOptionsPrefix(deluxe_ctx->seq_ksp,ksp_prefix);CHKERRQ(ierr);
414ac632422SStefano Zampini     ierr = KSPSetFromOptions(deluxe_ctx->seq_ksp);CHKERRQ(ierr);
415ac632422SStefano Zampini     ierr = KSPSetUp(deluxe_ctx->seq_ksp);CHKERRQ(ierr);
416ac632422SStefano Zampini   }
4175db18549SStefano Zampini   PetscFunctionReturn(0);
4185db18549SStefano Zampini }
419