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