xref: /petsc/src/ksp/pc/impls/bddc/bddcscalingbasic.c (revision 5a95e1ceb00aaf913daa567a70a604cb4fd63728)
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);
8*5a95e1ceSStefano 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;
38*5a95e1ceSStefano Zampini   ierr = VecSet(pcbddc->work_scaling,0.0);CHKERRQ(ierr);
39*5a95e1ceSStefano Zampini   ierr = VecSet(y,0.0);CHKERRQ(ierr);
40*5a95e1ceSStefano Zampini   if (deluxe_ctx->n_simple) { /* scale deluxe vertices using diagonal scaling */
41*5a95e1ceSStefano 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   }
54*5a95e1ceSStefano Zampini   /* sequential part : all problems and Schur applications collapsed into a single matrix vector multiplication */
55*5a95e1ceSStefano Zampini   if (deluxe_ctx->seq_mat) {
56674ae819SStefano Zampini     ierr = VecScatterBegin(deluxe_ctx->seq_scctx,x,deluxe_ctx->seq_work1,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
57674ae819SStefano Zampini     ierr = VecScatterEnd(deluxe_ctx->seq_scctx,x,deluxe_ctx->seq_work1,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
58*5a95e1ceSStefano Zampini     ierr = MatMultTranspose(deluxe_ctx->seq_mat,deluxe_ctx->seq_work1,deluxe_ctx->seq_work2);CHKERRQ(ierr);
59*5a95e1ceSStefano Zampini     ierr = VecScatterBegin(deluxe_ctx->seq_scctx,deluxe_ctx->seq_work2,pcbddc->work_scaling,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
60*5a95e1ceSStefano Zampini     ierr = VecScatterEnd(deluxe_ctx->seq_scctx,deluxe_ctx->seq_work2,pcbddc->work_scaling,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
61674ae819SStefano Zampini   }
62674ae819SStefano Zampini   /* put local boundary part in global vector */
63674ae819SStefano Zampini   ierr = VecScatterBegin(pcis->global_to_B,pcbddc->work_scaling,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
64674ae819SStefano Zampini   ierr = VecScatterEnd(pcis->global_to_B,pcbddc->work_scaling,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
65674ae819SStefano Zampini   PetscFunctionReturn(0);
66674ae819SStefano Zampini }
67674ae819SStefano Zampini 
68674ae819SStefano Zampini #undef __FUNCT__
69674ae819SStefano Zampini #define __FUNCT__ "PCBDDCScalingExtension"
70674ae819SStefano Zampini PetscErrorCode PCBDDCScalingExtension(PC pc, Vec local_interface_vector, Vec global_vector)
71674ae819SStefano Zampini {
72674ae819SStefano Zampini   PC_BDDC *pcbddc=(PC_BDDC*)pc->data;
73674ae819SStefano Zampini   PetscErrorCode ierr;
74674ae819SStefano Zampini 
75674ae819SStefano Zampini   PetscFunctionBegin;
76674ae819SStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
77674ae819SStefano Zampini   PetscValidHeaderSpecific(local_interface_vector,VEC_CLASSID,2);
78674ae819SStefano Zampini   PetscValidHeaderSpecific(global_vector,VEC_CLASSID,3);
79674ae819SStefano Zampini   if (local_interface_vector == pcbddc->work_scaling) {
80674ae819SStefano Zampini     SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"Local vector cannot be pcbddc->work_scaling!\n");
81674ae819SStefano Zampini   }
82674ae819SStefano Zampini   ierr = PetscTryMethod(pc,"PCBDDCScalingExtension_C",(PC,Vec,Vec),(pc,local_interface_vector,global_vector));CHKERRQ(ierr);
83674ae819SStefano Zampini   PetscFunctionReturn(0);
84674ae819SStefano Zampini }
85674ae819SStefano Zampini 
86674ae819SStefano Zampini #undef __FUNCT__
87674ae819SStefano Zampini #define __FUNCT__ "PCBDDCScalingRestriction_Basic"
88674ae819SStefano Zampini static PetscErrorCode PCBDDCScalingRestriction_Basic(PC pc, Vec global_vector, Vec local_interface_vector)
89674ae819SStefano Zampini {
90674ae819SStefano Zampini   PetscErrorCode ierr;
91674ae819SStefano Zampini   PC_IS* pcis = (PC_IS*)pc->data;
92674ae819SStefano Zampini 
93674ae819SStefano Zampini   PetscFunctionBegin;
94674ae819SStefano Zampini   ierr = VecScatterBegin(pcis->global_to_B,global_vector,local_interface_vector,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
95674ae819SStefano Zampini   ierr = VecScatterEnd(pcis->global_to_B,global_vector,local_interface_vector,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
96674ae819SStefano Zampini   /* Apply partition of unity */
97674ae819SStefano Zampini   ierr = VecPointwiseMult(local_interface_vector,pcis->D,local_interface_vector);CHKERRQ(ierr);
98674ae819SStefano Zampini   PetscFunctionReturn(0);
99674ae819SStefano Zampini }
100674ae819SStefano Zampini 
101674ae819SStefano Zampini #undef __FUNCT__
102674ae819SStefano Zampini #define __FUNCT__ "PCBDDCScalingRestriction_Deluxe"
103674ae819SStefano Zampini static PetscErrorCode PCBDDCScalingRestriction_Deluxe(PC pc, Vec x, Vec y)
104674ae819SStefano Zampini {
105674ae819SStefano Zampini   PC_IS*              pcis=(PC_IS*)pc->data;
106674ae819SStefano Zampini   PC_BDDC*            pcbddc=(PC_BDDC*)pc->data;
107674ae819SStefano Zampini   PCBDDCDeluxeScaling deluxe_ctx = pcbddc->deluxe_ctx;
108674ae819SStefano Zampini   PetscErrorCode      ierr;
109674ae819SStefano Zampini 
110674ae819SStefano Zampini   PetscFunctionBegin;
111674ae819SStefano Zampini   /* get local boundary part of global vector */
112674ae819SStefano Zampini   ierr = VecScatterBegin(pcis->global_to_B,x,y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
113674ae819SStefano Zampini   ierr = VecScatterEnd(pcis->global_to_B,x,y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
114*5a95e1ceSStefano Zampini   if (deluxe_ctx->n_simple) { /* scale deluxe vertices using diagonal scaling */
115*5a95e1ceSStefano Zampini     PetscInt          i;
1162b095fd8SStefano Zampini     PetscScalar       *array_y;
1172b095fd8SStefano Zampini     const PetscScalar *array_D;
118674ae819SStefano Zampini     ierr = VecGetArray(y,&array_y);CHKERRQ(ierr);
1192b095fd8SStefano Zampini     ierr = VecGetArrayRead(pcis->D,&array_D);CHKERRQ(ierr);
120674ae819SStefano Zampini     for (i=0;i<deluxe_ctx->n_simple;i++) {
121674ae819SStefano Zampini       array_y[deluxe_ctx->idx_simple_B[i]] *= array_D[deluxe_ctx->idx_simple_B[i]];
122674ae819SStefano Zampini     }
1232b095fd8SStefano Zampini     ierr = VecRestoreArrayRead(pcis->D,&array_D);CHKERRQ(ierr);
124674ae819SStefano Zampini     ierr = VecRestoreArray(y,&array_y);CHKERRQ(ierr);
12534a97f8cSStefano Zampini   }
12634a97f8cSStefano Zampini   /* sequential part : all problems and Schur applications collapsed into a single matrix vector multiplication and ksp solution */
127*5a95e1ceSStefano Zampini   if (deluxe_ctx->seq_mat) {
128674ae819SStefano Zampini     ierr = VecScatterBegin(deluxe_ctx->seq_scctx,y,deluxe_ctx->seq_work1,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
129674ae819SStefano Zampini     ierr = VecScatterEnd(deluxe_ctx->seq_scctx,y,deluxe_ctx->seq_work1,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
130*5a95e1ceSStefano Zampini     ierr = MatMult(deluxe_ctx->seq_mat,deluxe_ctx->seq_work1,deluxe_ctx->seq_work2);CHKERRQ(ierr);
131*5a95e1ceSStefano Zampini     ierr = VecScatterBegin(deluxe_ctx->seq_scctx,deluxe_ctx->seq_work2,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
132*5a95e1ceSStefano Zampini     ierr = VecScatterEnd(deluxe_ctx->seq_scctx,deluxe_ctx->seq_work2,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
133674ae819SStefano Zampini   }
134674ae819SStefano Zampini   PetscFunctionReturn(0);
135674ae819SStefano Zampini }
136674ae819SStefano Zampini 
137674ae819SStefano Zampini #undef __FUNCT__
138674ae819SStefano Zampini #define __FUNCT__ "PCBDDCScalingRestriction"
139674ae819SStefano Zampini PetscErrorCode PCBDDCScalingRestriction(PC pc, Vec global_vector, Vec local_interface_vector)
140674ae819SStefano Zampini {
141674ae819SStefano Zampini   PC_BDDC        *pcbddc=(PC_BDDC*)pc->data;
142674ae819SStefano Zampini   PetscErrorCode ierr;
143674ae819SStefano Zampini 
144674ae819SStefano Zampini   PetscFunctionBegin;
145674ae819SStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
146674ae819SStefano Zampini   PetscValidHeaderSpecific(global_vector,VEC_CLASSID,2);
147674ae819SStefano Zampini   PetscValidHeaderSpecific(local_interface_vector,VEC_CLASSID,3);
148674ae819SStefano Zampini   if (local_interface_vector == pcbddc->work_scaling) {
149674ae819SStefano Zampini     SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"Local vector should cannot be pcbddc->work_scaling!\n");
150674ae819SStefano Zampini   }
151674ae819SStefano Zampini   ierr = PetscTryMethod(pc,"PCBDDCScalingRestriction_C",(PC,Vec,Vec),(pc,global_vector,local_interface_vector));CHKERRQ(ierr);
152674ae819SStefano Zampini   PetscFunctionReturn(0);
153674ae819SStefano Zampini }
154674ae819SStefano Zampini 
155674ae819SStefano Zampini #undef __FUNCT__
156674ae819SStefano Zampini #define __FUNCT__ "PCBDDCScalingSetUp"
157674ae819SStefano Zampini PetscErrorCode PCBDDCScalingSetUp(PC pc)
158674ae819SStefano Zampini {
159674ae819SStefano Zampini   PC_IS* pcis=(PC_IS*)pc->data;
160674ae819SStefano Zampini   PC_BDDC* pcbddc=(PC_BDDC*)pc->data;
161674ae819SStefano Zampini   PetscErrorCode ierr;
162674ae819SStefano Zampini 
163674ae819SStefano Zampini   PetscFunctionBegin;
164674ae819SStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
165674ae819SStefano Zampini   /* create work vector for the operator */
16634a97f8cSStefano Zampini   ierr = VecDestroy(&pcbddc->work_scaling);CHKERRQ(ierr);
167674ae819SStefano Zampini   ierr = VecDuplicate(pcis->vec1_B,&pcbddc->work_scaling);CHKERRQ(ierr);
16834a97f8cSStefano Zampini   /* always rebuild pcis->D */
16928d874f6SStefano Zampini   if (pcis->use_stiffness_scaling) {
170674ae819SStefano Zampini     ierr = MatGetDiagonal(pcbddc->local_mat,pcis->vec1_N);CHKERRQ(ierr);
171674ae819SStefano Zampini     ierr = VecScatterBegin(pcis->N_to_B,pcis->vec1_N,pcis->D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
172674ae819SStefano Zampini     ierr = VecScatterEnd(pcis->N_to_B,pcis->vec1_N,pcis->D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
173674ae819SStefano Zampini   }
174674ae819SStefano Zampini   ierr = VecCopy(pcis->D,pcis->vec1_B);CHKERRQ(ierr);
175674ae819SStefano Zampini   ierr = VecSet(pcis->vec1_global,0.0);CHKERRQ(ierr);
176674ae819SStefano Zampini   ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_B,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
177674ae819SStefano Zampini   ierr = VecScatterEnd(pcis->global_to_B,pcis->vec1_B,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
178674ae819SStefano Zampini   ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_global,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
179674ae819SStefano Zampini   ierr = VecScatterEnd(pcis->global_to_B,pcis->vec1_global,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
180674ae819SStefano Zampini   ierr = VecPointwiseDivide(pcis->D,pcis->D,pcis->vec1_B);CHKERRQ(ierr);
181674ae819SStefano Zampini   /* now setup */
182681e7c04SStefano Zampini   if (pcbddc->use_deluxe_scaling) {
18334a97f8cSStefano Zampini     if (!pcbddc->deluxe_ctx) {
18434a97f8cSStefano Zampini       ierr = PCBDDCScalingCreate_Deluxe(pc);CHKERRQ(ierr);
18534a97f8cSStefano Zampini     }
18634a97f8cSStefano Zampini     ierr = PCBDDCScalingSetUp_Deluxe(pc);CHKERRQ(ierr);
187674ae819SStefano Zampini     ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCScalingRestriction_C",PCBDDCScalingRestriction_Deluxe);CHKERRQ(ierr);
188674ae819SStefano Zampini     ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCScalingExtension_C",PCBDDCScalingExtension_Deluxe);CHKERRQ(ierr);
189674ae819SStefano Zampini   } else {
190674ae819SStefano Zampini     ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCScalingRestriction_C",PCBDDCScalingRestriction_Basic);CHKERRQ(ierr);
191674ae819SStefano Zampini     ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCScalingExtension_C",PCBDDCScalingExtension_Basic);CHKERRQ(ierr);
192674ae819SStefano Zampini   }
19334a97f8cSStefano Zampini 
194674ae819SStefano Zampini   /* test */
195674ae819SStefano Zampini   if (pcbddc->dbg_flag) {
19634a97f8cSStefano Zampini     Vec         vec2_global;
197674ae819SStefano Zampini     PetscViewer viewer=pcbddc->dbg_viewer;
19834a97f8cSStefano Zampini     PetscReal   error;
199674ae819SStefano Zampini 
200674ae819SStefano Zampini     /* extension -> from local to parallel */
20134a97f8cSStefano Zampini     ierr = VecSet(pcis->vec1_global,0.0);CHKERRQ(ierr);
20234a97f8cSStefano Zampini     ierr = VecSetRandom(pcis->vec1_B,NULL);CHKERRQ(ierr);
20334a97f8cSStefano Zampini     ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_B,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
20434a97f8cSStefano Zampini     ierr = VecScatterEnd(pcis->global_to_B,pcis->vec1_B,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
20534a97f8cSStefano Zampini     ierr = VecDuplicate(pcis->vec1_global,&vec2_global);CHKERRQ(ierr);
20634a97f8cSStefano Zampini     ierr = VecCopy(pcis->vec1_global,vec2_global);CHKERRQ(ierr);
20734a97f8cSStefano Zampini 
208674ae819SStefano Zampini     ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_global,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
209674ae819SStefano Zampini     ierr = VecScatterEnd(pcis->global_to_B,pcis->vec1_global,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
210674ae819SStefano Zampini     ierr = PCBDDCScalingExtension(pc,pcis->vec1_B,pcis->vec1_global);CHKERRQ(ierr);
21134a97f8cSStefano Zampini     ierr = VecAXPY(pcis->vec1_global,-1.0,vec2_global);CHKERRQ(ierr);
21234a97f8cSStefano Zampini     ierr = VecNorm(pcis->vec1_global,NORM_INFINITY,&error);CHKERRQ(ierr);
213674ae819SStefano Zampini     ierr = PetscViewerASCIIPrintf(viewer,"Error scaling extension %1.14e\n",error);CHKERRQ(ierr);
21434a97f8cSStefano Zampini     if (error>1.e-8 && pcbddc->dbg_flag>1) {
215674ae819SStefano Zampini       ierr = VecView(pcis->vec1_global,viewer);CHKERRQ(ierr);
216674ae819SStefano Zampini     }
21734a97f8cSStefano Zampini     ierr = VecDestroy(&vec2_global);CHKERRQ(ierr);
21834a97f8cSStefano Zampini 
219674ae819SStefano Zampini     /* restriction -> from parallel to local */
220674ae819SStefano Zampini     ierr = VecSet(pcis->vec1_global,0.0);CHKERRQ(ierr);
22134a97f8cSStefano Zampini     ierr = VecSetRandom(pcis->vec1_B,NULL);CHKERRQ(ierr);
222674ae819SStefano Zampini     ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_B,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
223674ae819SStefano Zampini     ierr = VecScatterEnd(pcis->global_to_B,pcis->vec1_B,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
22434a97f8cSStefano Zampini 
22534a97f8cSStefano Zampini     ierr = PCBDDCScalingRestriction(pc,pcis->vec1_global,pcis->vec1_B);CHKERRQ(ierr);
22634a97f8cSStefano Zampini     ierr = VecScale(pcis->vec1_B,-1.0);CHKERRQ(ierr);
22734a97f8cSStefano Zampini     ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_B,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
22834a97f8cSStefano Zampini     ierr = VecScatterEnd(pcis->global_to_B,pcis->vec1_B,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
22934a97f8cSStefano Zampini     ierr = VecNorm(pcis->vec1_global,NORM_INFINITY,&error);CHKERRQ(ierr);
23034a97f8cSStefano Zampini     ierr = PetscViewerASCIIPrintf(viewer,"Error scaling restriction %1.14e\n",error);CHKERRQ(ierr);
23134a97f8cSStefano Zampini     if (error>1.e-8 && pcbddc->dbg_flag>1) {
232674ae819SStefano Zampini       ierr = VecView(pcis->vec1_global,viewer);CHKERRQ(ierr);
233674ae819SStefano Zampini     }
234674ae819SStefano Zampini   }
235674ae819SStefano Zampini   PetscFunctionReturn(0);
236674ae819SStefano Zampini }
237674ae819SStefano Zampini 
238674ae819SStefano Zampini #undef __FUNCT__
239674ae819SStefano Zampini #define __FUNCT__ "PCBDDCScalingDestroy"
240674ae819SStefano Zampini PetscErrorCode PCBDDCScalingDestroy(PC pc)
241674ae819SStefano Zampini {
242674ae819SStefano Zampini   PC_BDDC* pcbddc=(PC_BDDC*)pc->data;
243674ae819SStefano Zampini   PetscErrorCode ierr;
244674ae819SStefano Zampini 
245674ae819SStefano Zampini   PetscFunctionBegin;
24634a97f8cSStefano Zampini   if (pcbddc->deluxe_ctx) {
24734a97f8cSStefano Zampini     ierr = PCBDDCScalingDestroy_Deluxe(pc);CHKERRQ(ierr);
248674ae819SStefano Zampini   }
249674ae819SStefano Zampini   ierr = VecDestroy(&pcbddc->work_scaling);CHKERRQ(ierr);
250674ae819SStefano Zampini   /* remove functions */
251674ae819SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCScalingRestriction_C",NULL);CHKERRQ(ierr);
252674ae819SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCScalingExtension_C",NULL);CHKERRQ(ierr);
253674ae819SStefano Zampini   PetscFunctionReturn(0);
254674ae819SStefano Zampini }
255674ae819SStefano Zampini 
25634a97f8cSStefano Zampini #undef __FUNCT__
25734a97f8cSStefano Zampini #define __FUNCT__ "PCBDDCScalingCreate_Deluxe"
25834a97f8cSStefano Zampini static PetscErrorCode PCBDDCScalingCreate_Deluxe(PC pc)
25934a97f8cSStefano Zampini {
26034a97f8cSStefano Zampini   PC_BDDC*            pcbddc=(PC_BDDC*)pc->data;
26134a97f8cSStefano Zampini   PCBDDCDeluxeScaling deluxe_ctx;
26234a97f8cSStefano Zampini   PetscErrorCode      ierr;
26334a97f8cSStefano Zampini 
26434a97f8cSStefano Zampini   PetscFunctionBegin;
26534a97f8cSStefano Zampini   ierr = PetscNew(&deluxe_ctx);CHKERRQ(ierr);
26634a97f8cSStefano Zampini   pcbddc->deluxe_ctx = deluxe_ctx;
26734a97f8cSStefano Zampini   PetscFunctionReturn(0);
26834a97f8cSStefano Zampini }
26934a97f8cSStefano Zampini 
27034a97f8cSStefano Zampini #undef __FUNCT__
27134a97f8cSStefano Zampini #define __FUNCT__ "PCBDDCScalingDestroy_Deluxe"
27234a97f8cSStefano Zampini static PetscErrorCode PCBDDCScalingDestroy_Deluxe(PC pc)
27334a97f8cSStefano Zampini {
27434a97f8cSStefano Zampini   PC_BDDC*            pcbddc=(PC_BDDC*)pc->data;
27534a97f8cSStefano Zampini   PetscErrorCode      ierr;
27634a97f8cSStefano Zampini 
27734a97f8cSStefano Zampini   PetscFunctionBegin;
27834a97f8cSStefano Zampini   ierr = PCBDDCScalingReset_Deluxe_Solvers(pcbddc->deluxe_ctx);CHKERRQ(ierr);
27934a97f8cSStefano Zampini   ierr = PetscFree(pcbddc->deluxe_ctx);CHKERRQ(ierr);
28034a97f8cSStefano Zampini   PetscFunctionReturn(0);
28134a97f8cSStefano Zampini }
28234a97f8cSStefano Zampini 
28334a97f8cSStefano Zampini #undef __FUNCT__
28434a97f8cSStefano Zampini #define __FUNCT__ "PCBDDCScalingReset_Deluxe_Solvers"
28534a97f8cSStefano Zampini static PetscErrorCode PCBDDCScalingReset_Deluxe_Solvers(PCBDDCDeluxeScaling deluxe_ctx)
28634a97f8cSStefano Zampini {
28734a97f8cSStefano Zampini   PetscErrorCode      ierr;
28834a97f8cSStefano Zampini 
28934a97f8cSStefano Zampini   PetscFunctionBegin;
29034a97f8cSStefano Zampini   ierr = PetscFree(deluxe_ctx->idx_simple_B);CHKERRQ(ierr);
29134a97f8cSStefano Zampini   deluxe_ctx->n_simple = 0;
29234a97f8cSStefano Zampini   ierr = VecScatterDestroy(&deluxe_ctx->seq_scctx);CHKERRQ(ierr);
29334a97f8cSStefano Zampini   ierr = VecDestroy(&deluxe_ctx->seq_work1);CHKERRQ(ierr);
29434a97f8cSStefano Zampini   ierr = VecDestroy(&deluxe_ctx->seq_work2);CHKERRQ(ierr);
295*5a95e1ceSStefano Zampini   ierr = MatDestroy(&deluxe_ctx->seq_mat);CHKERRQ(ierr);
29634a97f8cSStefano Zampini   PetscFunctionReturn(0);
29734a97f8cSStefano Zampini }
29834a97f8cSStefano Zampini 
29934a97f8cSStefano Zampini #undef __FUNCT__
30034a97f8cSStefano Zampini #define __FUNCT__ "PCBDDCScalingSetUp_Deluxe"
30134a97f8cSStefano Zampini static PetscErrorCode PCBDDCScalingSetUp_Deluxe(PC pc)
30234a97f8cSStefano Zampini {
30334a97f8cSStefano Zampini   PC_IS               *pcis=(PC_IS*)pc->data;
30434a97f8cSStefano Zampini   PC_BDDC             *pcbddc=(PC_BDDC*)pc->data;
30534a97f8cSStefano Zampini   PCBDDCDeluxeScaling deluxe_ctx=pcbddc->deluxe_ctx;
306b96c3477SStefano Zampini   PCBDDCSubSchurs     sub_schurs=pcbddc->sub_schurs;
3071b968477SStefano Zampini   IS                  dirIS;
30834a97f8cSStefano Zampini   PetscErrorCode      ierr;
30934a97f8cSStefano Zampini 
31034a97f8cSStefano Zampini   PetscFunctionBegin;
311b96c3477SStefano Zampini   /* (TODO: reuse) throw away the solvers */
31234a97f8cSStefano Zampini   ierr = PCBDDCScalingReset_Deluxe_Solvers(deluxe_ctx);CHKERRQ(ierr);
31334a97f8cSStefano Zampini 
314b1b3d7a2SStefano Zampini   /* Compute data structures to solve sequential problems */
315*5a95e1ceSStefano Zampini   ierr = PCBDDCScalingSetUp_Deluxe_Private(pc);CHKERRQ(ierr);
3165db18549SStefano Zampini 
317b1b3d7a2SStefano Zampini   /* diagonal scaling on interface dofs not contained in cc */
3181b968477SStefano Zampini   dirIS = NULL;
3191b968477SStefano Zampini   if (pcbddc->DirichletBoundariesLocal) {
3201b968477SStefano Zampini     ierr = PCBDDCGraphGetDirichletDofs(pcbddc->mat_graph,&dirIS);CHKERRQ(ierr);
3211b968477SStefano Zampini   }
3221b968477SStefano Zampini   if (sub_schurs->is_Ej_com || dirIS) {
3231b968477SStefano Zampini     PetscInt n_com,n_dir;
3241b968477SStefano Zampini     n_com = 0;
3251b968477SStefano Zampini     if (sub_schurs->is_Ej_com) {
3261b968477SStefano Zampini       ierr = ISGetLocalSize(sub_schurs->is_Ej_com,&n_com);CHKERRQ(ierr);
3271b968477SStefano Zampini     }
3281b968477SStefano Zampini     n_dir = 0;
3291b968477SStefano Zampini     if (dirIS) {
3301b968477SStefano Zampini       ierr = ISGetLocalSize(dirIS,&n_dir);CHKERRQ(ierr);
3311b968477SStefano Zampini     }
3321b968477SStefano Zampini     deluxe_ctx->n_simple = n_dir + n_com;
3331b968477SStefano Zampini     ierr = PetscMalloc1(deluxe_ctx->n_simple,&deluxe_ctx->idx_simple_B);CHKERRQ(ierr);
3349bb4a8caSStefano Zampini     if (sub_schurs->is_Ej_com) {
3359bb4a8caSStefano Zampini       PetscInt       nmap;
3369bb4a8caSStefano Zampini       const PetscInt *idxs;
3371b968477SStefano Zampini 
3389bb4a8caSStefano Zampini       ierr = ISGetIndices(sub_schurs->is_Ej_com,&idxs);CHKERRQ(ierr);
3391b968477SStefano Zampini       ierr = ISGlobalToLocalMappingApply(pcis->BtoNmap,IS_GTOLM_DROP,n_com,idxs,&nmap,deluxe_ctx->idx_simple_B);CHKERRQ(ierr);
3401b968477SStefano Zampini       if (nmap != n_com) {
3411b968477SStefano Zampini         SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Error when mapping simply scaled dofs (is_Ej_com)! %d != %d",nmap,n_com);
3429bb4a8caSStefano Zampini       }
3439bb4a8caSStefano Zampini       ierr = ISRestoreIndices(sub_schurs->is_Ej_com,&idxs);CHKERRQ(ierr);
3441b968477SStefano Zampini     }
3451b968477SStefano Zampini     if (dirIS) {
3461b968477SStefano Zampini       PetscInt       nmap;
3471b968477SStefano Zampini       const PetscInt *idxs;
3481b968477SStefano Zampini 
3491b968477SStefano Zampini       ierr = ISGetIndices(dirIS,&idxs);CHKERRQ(ierr);
3501b968477SStefano Zampini       ierr = ISGlobalToLocalMappingApply(pcis->BtoNmap,IS_GTOLM_DROP,n_dir,idxs,&nmap,deluxe_ctx->idx_simple_B+n_com);CHKERRQ(ierr);
3511b968477SStefano Zampini       ierr = ISRestoreIndices(dirIS,&idxs);CHKERRQ(ierr);
3521b968477SStefano Zampini     }
3531b968477SStefano Zampini     ierr = PetscSortInt(deluxe_ctx->n_simple,deluxe_ctx->idx_simple_B);CHKERRQ(ierr);
3549bb4a8caSStefano Zampini   } else {
355b1b3d7a2SStefano Zampini     deluxe_ctx->n_simple = 0;
3569bb4a8caSStefano Zampini     deluxe_ctx->idx_simple_B = 0;
357b1b3d7a2SStefano Zampini   }
3581b968477SStefano Zampini   ierr = ISDestroy(&dirIS);CHKERRQ(ierr);
3591b968477SStefano Zampini 
36034a97f8cSStefano Zampini   PetscFunctionReturn(0);
36134a97f8cSStefano Zampini }
36234a97f8cSStefano Zampini 
36334a97f8cSStefano Zampini #undef __FUNCT__
364*5a95e1ceSStefano Zampini #define __FUNCT__ "PCBDDCScalingSetUp_Deluxe_Private"
365*5a95e1ceSStefano Zampini static PetscErrorCode PCBDDCScalingSetUp_Deluxe_Private(PC pc)
3665db18549SStefano Zampini {
3675db18549SStefano Zampini   PC_BDDC                *pcbddc=(PC_BDDC*)pc->data;
3685db18549SStefano Zampini   PCBDDCDeluxeScaling    deluxe_ctx=pcbddc->deluxe_ctx;
369b96c3477SStefano Zampini   PCBDDCSubSchurs        sub_schurs = pcbddc->sub_schurs;
3705db18549SStefano Zampini   PetscErrorCode         ierr;
3715db18549SStefano Zampini 
3725db18549SStefano Zampini   PetscFunctionBegin;
373*5a95e1ceSStefano Zampini   if (!sub_schurs->n_subs) {
3749221af80SStefano Zampini     PetscFunctionReturn(0);
3759221af80SStefano Zampini   }
3769221af80SStefano Zampini 
3775db18549SStefano Zampini   /* Create work vectors for sequential part of deluxe */
37865d8bf0aSStefano Zampini   ierr = MatCreateVecs(sub_schurs->sum_S_Ej_all,&deluxe_ctx->seq_work1,&deluxe_ctx->seq_work2);CHKERRQ(ierr);
3795db18549SStefano Zampini 
3805db18549SStefano Zampini   /* Compute deluxe sequential scatter */
3815db18549SStefano Zampini   ierr = VecScatterCreate(pcbddc->work_scaling,sub_schurs->is_Ej_all,deluxe_ctx->seq_work1,NULL,&deluxe_ctx->seq_scctx);CHKERRQ(ierr);
3825db18549SStefano Zampini 
383*5a95e1ceSStefano Zampini   /* Create Mat object for deluxe scaling */
384*5a95e1ceSStefano Zampini   ierr = PetscObjectReference((PetscObject)sub_schurs->S_Ej_all);CHKERRQ(ierr);
385*5a95e1ceSStefano Zampini   deluxe_ctx->seq_mat = sub_schurs->S_Ej_all;
3865db18549SStefano Zampini   PetscFunctionReturn(0);
3875db18549SStefano Zampini }
388