xref: /petsc/src/ksp/pc/impls/bddc/bddcscalingbasic.c (revision 8dccbbd225dabd1ed59e17588f15bb79a24e0cb0)
1 #include "bddc.h"
2 #include "bddcprivate.h"
3 
4 /* prototypes for deluxe public functions */
5 extern PetscErrorCode PCBDDCScalingCreateDeluxe(PC);
6 extern PetscErrorCode PCBDDCScalingDestroyDeluxe(PC);
7 extern PetscErrorCode PCBDDCScalingSetUpDeluxe(PC);
8 
9 #undef __FUNCT__
10 #define __FUNCT__ "PCBDDCScalingExtension_Basic"
11 static PetscErrorCode PCBDDCScalingExtension_Basic(PC pc, Vec local_interface_vector, Vec global_vector)
12 {
13   PC_IS* pcis = (PC_IS*)pc->data;
14   PC_BDDC* pcbddc = (PC_BDDC*)pc->data;
15   PetscErrorCode ierr;
16 
17   PetscFunctionBegin;
18   /* Apply partition of unity */
19   ierr = VecPointwiseMult(pcbddc->work_scaling,pcis->D,local_interface_vector);CHKERRQ(ierr);
20   ierr = VecSet(global_vector,0.0);CHKERRQ(ierr);
21   ierr = VecScatterBegin(pcis->global_to_B,pcbddc->work_scaling,global_vector,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
22   ierr = VecScatterEnd(pcis->global_to_B,pcbddc->work_scaling,global_vector,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
23   PetscFunctionReturn(0);
24 }
25 
26 #undef __FUNCT__
27 #define __FUNCT__ "PCBDDCScalingExtension_Deluxe"
28 static PetscErrorCode PCBDDCScalingExtension_Deluxe(PC pc, Vec x, Vec y)
29 {
30   PC_IS*              pcis=(PC_IS*)pc->data;
31   PC_BDDC*            pcbddc=(PC_BDDC*)pc->data;
32   PCBDDCDeluxeScaling deluxe_ctx = pcbddc->deluxe_ctx;
33   PetscScalar         *array_x,*array_D,*array;
34   PetscScalar         zero=0.0;
35   PetscInt            i;
36   PetscMPIInt         color_rank;
37   PetscErrorCode      ierr;
38 
39   /* TODO CHECK STUFF RELATED WITH FAKE WORK */
40   PetscFunctionBegin;
41   ierr = VecSet(pcbddc->work_scaling,zero);CHKERRQ(ierr); /* needed by the fake work below */
42   /* scale deluxe vertices using diagonal scaling */
43   ierr = VecGetArray(x,&array_x);CHKERRQ(ierr);
44   ierr = VecGetArray(pcis->D,&array_D);CHKERRQ(ierr);
45   ierr = VecGetArray(pcbddc->work_scaling,&array);CHKERRQ(ierr);
46   for (i=0;i<deluxe_ctx->n_simple;i++) {
47     array[deluxe_ctx->idx_simple_B[i]] = array_x[deluxe_ctx->idx_simple_B[i]]*array_D[deluxe_ctx->idx_simple_B[i]];
48   }
49   ierr = VecRestoreArray(pcbddc->work_scaling,&array);CHKERRQ(ierr);
50   ierr = VecRestoreArray(pcis->D,&array_D);CHKERRQ(ierr);
51   ierr = VecRestoreArray(x,&array_x);CHKERRQ(ierr);
52   /* sequential part : all problems and Schur applications collapsed into seq_mat at setup phase */
53   if (deluxe_ctx->seq_mat) {
54     ierr = VecScatterBegin(deluxe_ctx->seq_scctx,x,deluxe_ctx->seq_work1,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
55     ierr = VecScatterEnd(deluxe_ctx->seq_scctx,x,deluxe_ctx->seq_work1,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
56     ierr = MatMult(deluxe_ctx->seq_mat,deluxe_ctx->seq_work1,deluxe_ctx->seq_work2);CHKERRQ(ierr);
57     ierr = KSPSolve(deluxe_ctx->seq_ksp,deluxe_ctx->seq_work2,deluxe_ctx->seq_work1);CHKERRQ(ierr);
58     /* fake work due to final ADD VALUES and vertices scaling needed? TODO: check it */
59     ierr = VecScatterBegin(deluxe_ctx->seq_scctx,deluxe_ctx->seq_work1,pcbddc->work_scaling,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
60     ierr = VecScatterEnd(deluxe_ctx->seq_scctx,deluxe_ctx->seq_work1,pcbddc->work_scaling,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
61   }
62   /* parallel part */
63   for (i=0;i<deluxe_ctx->par_colors;i++) {
64     if (deluxe_ctx->par_ksp[i]) {
65       ierr = MPI_Comm_rank(deluxe_ctx->par_subcomm[i]->comm,&color_rank);CHKERRQ(ierr);
66       ierr = VecSet(deluxe_ctx->work1_B,zero);CHKERRQ(ierr);
67       ierr = VecScatterBegin(deluxe_ctx->par_scctx_s[i],x,deluxe_ctx->work1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
68       ierr = VecScatterEnd(deluxe_ctx->par_scctx_s[i],x,deluxe_ctx->work1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
69       /* apply local schur on subset S_j^-1 */
70       ierr = PCBDDCApplySchur(pc,deluxe_ctx->work1_B,deluxe_ctx->work2_B,(Vec)0,deluxe_ctx->work1_D,deluxe_ctx->work2_D);CHKERRQ(ierr);
71       /* parallel transpose solve (\sum_j S_j)^-1 */
72       ierr = VecSet(deluxe_ctx->par_vec[i],zero);CHKERRQ(ierr);
73       ierr = VecScatterBegin(deluxe_ctx->par_scctx_p[i],deluxe_ctx->work2_B,deluxe_ctx->par_vec[i],ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
74       ierr = VecScatterEnd(deluxe_ctx->par_scctx_p[i],deluxe_ctx->work2_B,deluxe_ctx->par_vec[i],ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
75       ierr = KSPSolve(deluxe_ctx->par_ksp[i],deluxe_ctx->par_vec[i],deluxe_ctx->par_vec[i]);CHKERRQ(ierr);
76       if (!color_rank) {
77         ierr = VecScatterBegin(deluxe_ctx->par_scctx_p[i],deluxe_ctx->par_vec[i],pcbddc->work_scaling,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
78         ierr = VecScatterEnd(deluxe_ctx->par_scctx_p[i],deluxe_ctx->par_vec[i],pcbddc->work_scaling,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
79       } else { /* fake work due to final ADD VALUES and vertices scaling */
80         ierr = VecScatterBegin(deluxe_ctx->par_scctx_p[i],deluxe_ctx->par_vec[i],deluxe_ctx->work1_B,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
81         ierr = VecScatterEnd(deluxe_ctx->par_scctx_p[i],deluxe_ctx->par_vec[i],deluxe_ctx->work1_B,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
82       }
83     }
84   }
85   /* put local boundary part in global vector */
86   ierr = VecSet(y,zero);CHKERRQ(ierr);
87   ierr = VecScatterBegin(pcis->global_to_B,pcbddc->work_scaling,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
88   ierr = VecScatterEnd  (pcis->global_to_B,pcbddc->work_scaling,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
89   PetscFunctionReturn(0);
90 }
91 
92 #undef __FUNCT__
93 #define __FUNCT__ "PCBDDCScalingExtension"
94 PetscErrorCode PCBDDCScalingExtension(PC pc, Vec local_interface_vector, Vec global_vector)
95 {
96   PC_BDDC *pcbddc=(PC_BDDC*)pc->data;
97   PetscErrorCode ierr;
98 
99   PetscFunctionBegin;
100   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
101   PetscValidHeaderSpecific(local_interface_vector,VEC_CLASSID,2);
102   PetscValidHeaderSpecific(global_vector,VEC_CLASSID,3);
103   if (local_interface_vector == pcbddc->work_scaling) {
104     SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"Local vector cannot be pcbddc->work_scaling!\n");
105   }
106   ierr = PetscTryMethod(pc,"PCBDDCScalingExtension_C",(PC,Vec,Vec),(pc,local_interface_vector,global_vector));CHKERRQ(ierr);
107   PetscFunctionReturn(0);
108 }
109 
110 #undef __FUNCT__
111 #define __FUNCT__ "PCBDDCScalingRestriction_Basic"
112 static PetscErrorCode PCBDDCScalingRestriction_Basic(PC pc, Vec global_vector, Vec local_interface_vector)
113 {
114   PetscErrorCode ierr;
115   PC_IS* pcis = (PC_IS*)pc->data;
116 
117   PetscFunctionBegin;
118   ierr = VecScatterBegin(pcis->global_to_B,global_vector,local_interface_vector,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
119   ierr = VecScatterEnd(pcis->global_to_B,global_vector,local_interface_vector,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
120   /* Apply partition of unity */
121   ierr = VecPointwiseMult(local_interface_vector,pcis->D,local_interface_vector);CHKERRQ(ierr);
122   PetscFunctionReturn(0);
123 }
124 
125 #undef __FUNCT__
126 #define __FUNCT__ "PCBDDCScalingRestriction_Deluxe"
127 static PetscErrorCode PCBDDCScalingRestriction_Deluxe(PC pc, Vec x, Vec y)
128 {
129   PC_IS*              pcis=(PC_IS*)pc->data;
130   PC_BDDC*            pcbddc=(PC_BDDC*)pc->data;
131   PCBDDCDeluxeScaling deluxe_ctx = pcbddc->deluxe_ctx;
132   PetscScalar         *array_y,*array_D,zero=0.0;
133   PetscInt            i;
134   PetscErrorCode      ierr;
135 
136   PetscFunctionBegin;
137   /* get local boundary part of global vector */
138   ierr = VecScatterBegin(pcis->global_to_B,x,y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
139   ierr = VecScatterEnd  (pcis->global_to_B,x,y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
140   /* scale vertices using diagonal scaling -> every scaling perform the same */
141   ierr = VecGetArray(y,&array_y);CHKERRQ(ierr);
142   ierr = VecGetArray(pcis->D,&array_D);CHKERRQ(ierr);
143   for (i=0;i<deluxe_ctx->n_simple;i++) {
144     array_y[deluxe_ctx->idx_simple_B[i]] *= array_D[deluxe_ctx->idx_simple_B[i]];
145   }
146   ierr = VecRestoreArray(pcis->D,&array_D);CHKERRQ(ierr);
147   ierr = VecRestoreArray(y,&array_y);CHKERRQ(ierr);
148   /* sequential part : all problems and Schur applications collapsed into seq_mat at setup phase */
149   if (deluxe_ctx->seq_mat) {
150     ierr = VecScatterBegin(deluxe_ctx->seq_scctx,y,deluxe_ctx->seq_work1,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
151     ierr = VecScatterEnd(deluxe_ctx->seq_scctx,y,deluxe_ctx->seq_work1,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
152     ierr = KSPSolveTranspose(deluxe_ctx->seq_ksp,deluxe_ctx->seq_work1,deluxe_ctx->seq_work2);CHKERRQ(ierr);
153     ierr = MatMultTranspose(deluxe_ctx->seq_mat,deluxe_ctx->seq_work2,deluxe_ctx->seq_work1);CHKERRQ(ierr);
154     ierr = VecScatterBegin(deluxe_ctx->seq_scctx,deluxe_ctx->seq_work1,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
155     ierr = VecScatterEnd(deluxe_ctx->seq_scctx,deluxe_ctx->seq_work1,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
156   }
157   /* parallel part */
158   for (i=0;i<deluxe_ctx->par_colors;i++) {
159     if (deluxe_ctx->par_ksp[i]) {
160       /* parallel solve (\sum_j S_j)^-1 */
161       ierr = VecSet(deluxe_ctx->par_vec[i],zero);CHKERRQ(ierr);
162       ierr = VecScatterBegin(deluxe_ctx->par_scctx_p[i],y,deluxe_ctx->par_vec[i],INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
163       ierr = VecScatterEnd(deluxe_ctx->par_scctx_p[i],y,deluxe_ctx->par_vec[i],INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
164       ierr = KSPSolveTranspose(deluxe_ctx->par_ksp[i],deluxe_ctx->par_vec[i],deluxe_ctx->par_vec[i]);CHKERRQ(ierr);
165       /* apply local schur S_j^-1 */
166       ierr = VecSet(deluxe_ctx->work1_B,zero);CHKERRQ(ierr);
167       ierr = VecScatterBegin(deluxe_ctx->par_scctx_p[i],deluxe_ctx->par_vec[i],deluxe_ctx->work1_B,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
168       ierr = VecScatterEnd(deluxe_ctx->par_scctx_p[i],deluxe_ctx->par_vec[i],deluxe_ctx->work1_B,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
169       ierr = PCBDDCApplySchurTranspose(pc,deluxe_ctx->work1_B,deluxe_ctx->work2_B,(Vec)0,deluxe_ctx->work1_D,deluxe_ctx->work2_D);CHKERRQ(ierr);
170       ierr = VecScatterBegin(deluxe_ctx->par_scctx_s[i],deluxe_ctx->work2_B,y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
171       ierr = VecScatterEnd(deluxe_ctx->par_scctx_s[i],deluxe_ctx->work2_B,y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
172     }
173   }
174   PetscFunctionReturn(0);
175 }
176 
177 #undef __FUNCT__
178 #define __FUNCT__ "PCBDDCScalingRestriction"
179 PetscErrorCode PCBDDCScalingRestriction(PC pc, Vec global_vector, Vec local_interface_vector)
180 {
181   PC_BDDC        *pcbddc=(PC_BDDC*)pc->data;
182   PetscErrorCode ierr;
183 
184   PetscFunctionBegin;
185   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
186   PetscValidHeaderSpecific(global_vector,VEC_CLASSID,2);
187   PetscValidHeaderSpecific(local_interface_vector,VEC_CLASSID,3);
188   if (local_interface_vector == pcbddc->work_scaling) {
189     SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"Local vector should cannot be pcbddc->work_scaling!\n");
190   }
191   ierr = PetscTryMethod(pc,"PCBDDCScalingRestriction_C",(PC,Vec,Vec),(pc,global_vector,local_interface_vector));CHKERRQ(ierr);
192   PetscFunctionReturn(0);
193 }
194 
195 #undef __FUNCT__
196 #define __FUNCT__ "PCBDDCScalingSetUp"
197 PetscErrorCode PCBDDCScalingSetUp(PC pc)
198 {
199   PC_IS* pcis=(PC_IS*)pc->data;
200   PC_BDDC* pcbddc=(PC_BDDC*)pc->data;
201   PetscErrorCode ierr;
202 
203   PetscFunctionBegin;
204   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
205   /* create work vector for the operator */
206   ierr = VecDuplicate(pcis->vec1_B,&pcbddc->work_scaling);CHKERRQ(ierr);
207   /* rebuild pcis->D if stiffness scaling has been requested */
208   if (pcis->use_stiffness_scaling) {
209     ierr = MatGetDiagonal(pcbddc->local_mat,pcis->vec1_N);CHKERRQ(ierr);
210     ierr = VecScatterBegin(pcis->N_to_B,pcis->vec1_N,pcis->D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
211     ierr = VecScatterEnd(pcis->N_to_B,pcis->vec1_N,pcis->D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
212   }
213   ierr = VecCopy(pcis->D,pcis->vec1_B);CHKERRQ(ierr);
214   ierr = VecSet(pcis->vec1_global,0.0);CHKERRQ(ierr);
215   ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_B,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
216   ierr = VecScatterEnd(pcis->global_to_B,pcis->vec1_B,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
217   ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_global,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
218   ierr = VecScatterEnd(pcis->global_to_B,pcis->vec1_global,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
219   ierr = VecPointwiseDivide(pcis->D,pcis->D,pcis->vec1_B);CHKERRQ(ierr);
220   /* now setup */
221   /* if (pcbddc->use_deluxe_scaling) { */
222   if (PETSC_FALSE) {
223     ierr = PCBDDCScalingCreateDeluxe(pc);CHKERRQ(ierr);
224     ierr = PCBDDCScalingSetUpDeluxe(pc);CHKERRQ(ierr);
225     ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCScalingRestriction_C",PCBDDCScalingRestriction_Deluxe);CHKERRQ(ierr);
226     ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCScalingExtension_C",PCBDDCScalingExtension_Deluxe);CHKERRQ(ierr);
227   } else {
228     ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCScalingRestriction_C",PCBDDCScalingRestriction_Basic);CHKERRQ(ierr);
229     ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCScalingExtension_C",PCBDDCScalingExtension_Basic);CHKERRQ(ierr);
230   }
231   /* test */
232   if (pcbddc->dbg_flag) {
233     PetscViewer viewer=pcbddc->dbg_viewer;
234     PetscReal   error,gerror;
235     MPI_Comm    test_comm;
236 
237     /* extension -> from local to parallel */
238     ierr = PetscObjectGetComm((PetscObject)pc,&test_comm);CHKERRQ(ierr);
239     ierr = VecSetRandom(pcis->vec1_global,NULL);CHKERRQ(ierr);
240     ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_global,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
241     ierr = VecScatterEnd(pcis->global_to_B,pcis->vec1_global,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
242     ierr = PCBDDCScalingExtension(pc,pcis->vec1_B,pcis->vec1_global);CHKERRQ(ierr);
243     ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_global,pcbddc->work_scaling,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
244     ierr = VecScatterEnd(pcis->global_to_B,pcis->vec1_global,pcbddc->work_scaling,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
245     ierr = VecAXPY(pcis->vec1_B,-1.0,pcbddc->work_scaling);CHKERRQ(ierr);
246     ierr = VecNorm(pcis->vec1_B,NORM_INFINITY,&error);CHKERRQ(ierr);
247     ierr = MPI_Allreduce(&error,&gerror,1,MPIU_REAL,MPI_MAX,test_comm);CHKERRQ(ierr);
248     ierr = PetscViewerASCIIPrintf(viewer,"Error scaling extension %1.14e\n",error);CHKERRQ(ierr);
249     if (PetscAbsReal(gerror)>1.e-8) {
250       ierr = VecSet(pcis->vec1_global,0.0);CHKERRQ(ierr);
251       ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_B,pcis->vec1_global,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
252       ierr = VecScatterEnd(pcis->global_to_B,pcis->vec1_B,pcis->vec1_global,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
253       ierr = VecView(pcis->vec1_global,viewer);CHKERRQ(ierr);
254     }
255     /* restriction -> from parallel to local */
256     ierr = VecSetRandom(pcis->vec1_global,NULL);CHKERRQ(ierr);
257     ierr = PCBDDCScalingRestriction(pc,pcis->vec1_global,pcis->vec1_B);CHKERRQ(ierr);
258     ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_global,pcbddc->work_scaling,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
259     ierr = VecScatterEnd(pcis->global_to_B,pcis->vec1_global,pcbddc->work_scaling,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
260     ierr = VecSet(pcis->vec1_global,0.0);CHKERRQ(ierr);
261     ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_B,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
262     ierr = VecScatterEnd(pcis->global_to_B,pcis->vec1_B,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
263     ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_global,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
264     ierr = VecScatterEnd(pcis->global_to_B,pcis->vec1_global,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
265     ierr = VecAXPY(pcis->vec1_B,-1.0,pcbddc->work_scaling);CHKERRQ(ierr);
266     ierr = VecNorm(pcis->vec1_B,NORM_INFINITY,&error);CHKERRQ(ierr);
267     ierr = MPI_Allreduce(&error,&gerror,1,MPIU_REAL,MPI_MAX,test_comm);CHKERRQ(ierr);
268     ierr = PetscViewerASCIIPrintf(viewer,"Error scaling restriction %1.14e\n",gerror);CHKERRQ(ierr);
269     if (PetscAbsReal(gerror)>1.e-8) {
270       ierr = VecSet(pcis->vec1_global,0.0);CHKERRQ(ierr);
271       ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_B,pcis->vec1_global,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
272       ierr = VecScatterEnd(pcis->global_to_B,pcis->vec1_B,pcis->vec1_global,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
273       ierr = VecView(pcis->vec1_global,viewer);CHKERRQ(ierr);
274     }
275   }
276   PetscFunctionReturn(0);
277 }
278 
279 #undef __FUNCT__
280 #define __FUNCT__ "PCBDDCScalingDestroy"
281 PetscErrorCode PCBDDCScalingDestroy(PC pc)
282 {
283   PC_BDDC* pcbddc=(PC_BDDC*)pc->data;
284   PetscErrorCode ierr;
285 
286   PetscFunctionBegin;
287   /* if (pcbddc->use_deluxe_scaling) { */
288   if (PETSC_FALSE) {
289     ierr = PCBDDCScalingDestroyDeluxe(pc);CHKERRQ(ierr);
290   }
291   ierr = VecDestroy(&pcbddc->work_scaling);CHKERRQ(ierr);
292   /* remove functions */
293   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCScalingRestriction_C",NULL);CHKERRQ(ierr);
294   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCScalingExtension_C",NULL);CHKERRQ(ierr);
295   PetscFunctionReturn(0);
296 }
297 
298