xref: /petsc/src/ksp/pc/impls/bddc/bddcscalingbasic.c (revision 2b095fd868400130684ee72a1c39d4c0d6f4ebf1)
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);
834a97f8cSStefano Zampini static PetscErrorCode PCBDDCScalingSetUp_Deluxe_Par(PC,PetscInt,PetscInt,PetscInt[],PetscInt[]);
9883469d8SStefano Zampini static PetscErrorCode PCBDDCScalingSetUp_Deluxe_Seq(PC);
1034a97f8cSStefano Zampini static PetscErrorCode PCBDDCScalingReset_Deluxe_Solvers(PCBDDCDeluxeScaling);
11674ae819SStefano Zampini 
12674ae819SStefano Zampini #undef __FUNCT__
13674ae819SStefano Zampini #define __FUNCT__ "PCBDDCScalingExtension_Basic"
14674ae819SStefano Zampini static PetscErrorCode PCBDDCScalingExtension_Basic(PC pc, Vec local_interface_vector, Vec global_vector)
15674ae819SStefano Zampini {
16674ae819SStefano Zampini   PC_IS* pcis = (PC_IS*)pc->data;
17674ae819SStefano Zampini   PC_BDDC* pcbddc = (PC_BDDC*)pc->data;
18674ae819SStefano Zampini   PetscErrorCode ierr;
19674ae819SStefano Zampini 
20674ae819SStefano Zampini   PetscFunctionBegin;
21674ae819SStefano Zampini   /* Apply partition of unity */
22674ae819SStefano Zampini   ierr = VecPointwiseMult(pcbddc->work_scaling,pcis->D,local_interface_vector);CHKERRQ(ierr);
23674ae819SStefano Zampini   ierr = VecSet(global_vector,0.0);CHKERRQ(ierr);
24674ae819SStefano Zampini   ierr = VecScatterBegin(pcis->global_to_B,pcbddc->work_scaling,global_vector,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
25674ae819SStefano Zampini   ierr = VecScatterEnd(pcis->global_to_B,pcbddc->work_scaling,global_vector,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
26674ae819SStefano Zampini   PetscFunctionReturn(0);
27674ae819SStefano Zampini }
28674ae819SStefano Zampini 
29674ae819SStefano Zampini #undef __FUNCT__
30674ae819SStefano Zampini #define __FUNCT__ "PCBDDCScalingExtension_Deluxe"
31674ae819SStefano Zampini static PetscErrorCode PCBDDCScalingExtension_Deluxe(PC pc, Vec x, Vec y)
32674ae819SStefano Zampini {
33674ae819SStefano Zampini   PC_IS*              pcis=(PC_IS*)pc->data;
34674ae819SStefano Zampini   PC_BDDC*            pcbddc=(PC_BDDC*)pc->data;
35674ae819SStefano Zampini   PCBDDCDeluxeScaling deluxe_ctx = pcbddc->deluxe_ctx;
36b96c3477SStefano Zampini   PCBDDCSubSchurs     sub_schurs = pcbddc->sub_schurs;
37674ae819SStefano Zampini   PetscInt            i;
38674ae819SStefano Zampini   PetscErrorCode      ierr;
39674ae819SStefano Zampini 
40674ae819SStefano Zampini   /* TODO CHECK STUFF RELATED WITH FAKE WORK */
41674ae819SStefano Zampini   PetscFunctionBegin;
4234a97f8cSStefano Zampini   ierr = VecSet(pcbddc->work_scaling,0.0);CHKERRQ(ierr); /* needed by the fake work below */
4334a97f8cSStefano Zampini   if (deluxe_ctx->n_simple) {
44674ae819SStefano Zampini     /* scale deluxe vertices using diagonal scaling */
45*2b095fd8SStefano Zampini     const PetscScalar *array_x,*array_D;
46*2b095fd8SStefano Zampini     PetscScalar       *array;
47*2b095fd8SStefano Zampini     ierr = VecGetArrayRead(x,&array_x);CHKERRQ(ierr);
48*2b095fd8SStefano Zampini     ierr = VecGetArrayRead(pcis->D,&array_D);CHKERRQ(ierr);
49674ae819SStefano Zampini     ierr = VecGetArray(pcbddc->work_scaling,&array);CHKERRQ(ierr);
50674ae819SStefano Zampini     for (i=0;i<deluxe_ctx->n_simple;i++) {
51674ae819SStefano Zampini       array[deluxe_ctx->idx_simple_B[i]] = array_x[deluxe_ctx->idx_simple_B[i]]*array_D[deluxe_ctx->idx_simple_B[i]];
52674ae819SStefano Zampini     }
53674ae819SStefano Zampini     ierr = VecRestoreArray(pcbddc->work_scaling,&array);CHKERRQ(ierr);
54*2b095fd8SStefano Zampini     ierr = VecRestoreArrayRead(pcis->D,&array_D);CHKERRQ(ierr);
55*2b095fd8SStefano Zampini     ierr = VecRestoreArrayRead(x,&array_x);CHKERRQ(ierr);
5634a97f8cSStefano Zampini   }
5734a97f8cSStefano Zampini   /* sequential part : all problems and Schur applications collapsed into a single matrix vector multiplication and ksp solution */
5841c3ba1bSStefano Zampini   if (deluxe_ctx->seq_ksp) {
59674ae819SStefano Zampini     ierr = VecScatterBegin(deluxe_ctx->seq_scctx,x,deluxe_ctx->seq_work1,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
60674ae819SStefano Zampini     ierr = VecScatterEnd(deluxe_ctx->seq_scctx,x,deluxe_ctx->seq_work1,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
611d850880SStefano Zampini     ierr = MatMultTranspose(sub_schurs->S_Ej_all,deluxe_ctx->seq_work1,deluxe_ctx->seq_work2);CHKERRQ(ierr);
621d850880SStefano Zampini     ierr = KSPSolveTranspose(deluxe_ctx->seq_ksp,deluxe_ctx->seq_work2,deluxe_ctx->seq_work1);CHKERRQ(ierr);
63674ae819SStefano Zampini     /* fake work due to final ADD VALUES and vertices scaling needed? TODO: check it */
64674ae819SStefano Zampini     ierr = VecScatterBegin(deluxe_ctx->seq_scctx,deluxe_ctx->seq_work1,pcbddc->work_scaling,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
65674ae819SStefano Zampini     ierr = VecScatterEnd(deluxe_ctx->seq_scctx,deluxe_ctx->seq_work1,pcbddc->work_scaling,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
66674ae819SStefano Zampini   }
67674ae819SStefano Zampini   /* parallel part */
68674ae819SStefano Zampini   for (i=0;i<deluxe_ctx->par_colors;i++) {
69674ae819SStefano Zampini     if (deluxe_ctx->par_ksp[i]) {
7034a97f8cSStefano Zampini       PetscMPIInt color_rank;
7134a97f8cSStefano Zampini       PetscInt    subidx = deluxe_ctx->par_col2sub[i];
7234a97f8cSStefano Zampini       /* restrict on subset */
73b96c3477SStefano Zampini       ierr = VecScatterBegin(deluxe_ctx->par_scctx_s[i],x,deluxe_ctx->par_work1[i],INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
74b96c3477SStefano Zampini       ierr = VecScatterEnd(deluxe_ctx->par_scctx_s[i],x,deluxe_ctx->par_work1[i],INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
7534a97f8cSStefano Zampini       /* S_Ej */
76b96c3477SStefano Zampini       ierr = MatMult(sub_schurs->S_Ej[subidx],deluxe_ctx->par_work1[i],deluxe_ctx->par_work2[i]);CHKERRQ(ierr);
7734a97f8cSStefano Zampini       /* (\sum_j S_Ej)^-1 */
7834a97f8cSStefano Zampini       ierr = VecSet(deluxe_ctx->par_vec[i],0.0);CHKERRQ(ierr);
79b96c3477SStefano Zampini       ierr = VecScatterBegin(deluxe_ctx->par_scctx_p[i],deluxe_ctx->par_work2[i],deluxe_ctx->par_vec[i],ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
80b96c3477SStefano Zampini       ierr = VecScatterEnd(deluxe_ctx->par_scctx_p[i],deluxe_ctx->par_work2[i],deluxe_ctx->par_vec[i],ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
81674ae819SStefano Zampini       ierr = KSPSolve(deluxe_ctx->par_ksp[i],deluxe_ctx->par_vec[i],deluxe_ctx->par_vec[i]);CHKERRQ(ierr);
8234a97f8cSStefano Zampini       ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)deluxe_ctx->par_ksp[i]),&color_rank);CHKERRQ(ierr);
8334a97f8cSStefano Zampini       /* get back solution on subset */
84b96c3477SStefano Zampini       ierr = VecScatterBegin(deluxe_ctx->par_scctx_p[i],deluxe_ctx->par_vec[i],deluxe_ctx->par_work1[i],INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
85b96c3477SStefano Zampini       ierr = VecScatterEnd(deluxe_ctx->par_scctx_p[i],deluxe_ctx->par_vec[i],deluxe_ctx->par_work1[i],INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
8634a97f8cSStefano Zampini       if (!color_rank) { /* only the master process in coloured comm copies the computed values */
87b96c3477SStefano Zampini         ierr = VecScatterBegin(deluxe_ctx->par_scctx_s[i],deluxe_ctx->par_work1[i],pcbddc->work_scaling,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
88b96c3477SStefano Zampini         ierr = VecScatterEnd(deluxe_ctx->par_scctx_s[i],deluxe_ctx->par_work1[i],pcbddc->work_scaling,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
89674ae819SStefano Zampini       }
90674ae819SStefano Zampini     }
91674ae819SStefano Zampini   }
92674ae819SStefano Zampini   /* put local boundary part in global vector */
9334a97f8cSStefano Zampini   ierr = VecSet(y,0.0);CHKERRQ(ierr);
94674ae819SStefano Zampini   ierr = VecScatterBegin(pcis->global_to_B,pcbddc->work_scaling,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
95674ae819SStefano Zampini   ierr = VecScatterEnd(pcis->global_to_B,pcbddc->work_scaling,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
96674ae819SStefano Zampini   PetscFunctionReturn(0);
97674ae819SStefano Zampini }
98674ae819SStefano Zampini 
99674ae819SStefano Zampini #undef __FUNCT__
100674ae819SStefano Zampini #define __FUNCT__ "PCBDDCScalingExtension"
101674ae819SStefano Zampini PetscErrorCode PCBDDCScalingExtension(PC pc, Vec local_interface_vector, Vec global_vector)
102674ae819SStefano Zampini {
103674ae819SStefano Zampini   PC_BDDC *pcbddc=(PC_BDDC*)pc->data;
104674ae819SStefano Zampini   PetscErrorCode ierr;
105674ae819SStefano Zampini 
106674ae819SStefano Zampini   PetscFunctionBegin;
107674ae819SStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
108674ae819SStefano Zampini   PetscValidHeaderSpecific(local_interface_vector,VEC_CLASSID,2);
109674ae819SStefano Zampini   PetscValidHeaderSpecific(global_vector,VEC_CLASSID,3);
110674ae819SStefano Zampini   if (local_interface_vector == pcbddc->work_scaling) {
111674ae819SStefano Zampini     SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"Local vector cannot be pcbddc->work_scaling!\n");
112674ae819SStefano Zampini   }
113674ae819SStefano Zampini   ierr = PetscTryMethod(pc,"PCBDDCScalingExtension_C",(PC,Vec,Vec),(pc,local_interface_vector,global_vector));CHKERRQ(ierr);
114674ae819SStefano Zampini   PetscFunctionReturn(0);
115674ae819SStefano Zampini }
116674ae819SStefano Zampini 
117674ae819SStefano Zampini #undef __FUNCT__
118674ae819SStefano Zampini #define __FUNCT__ "PCBDDCScalingRestriction_Basic"
119674ae819SStefano Zampini static PetscErrorCode PCBDDCScalingRestriction_Basic(PC pc, Vec global_vector, Vec local_interface_vector)
120674ae819SStefano Zampini {
121674ae819SStefano Zampini   PetscErrorCode ierr;
122674ae819SStefano Zampini   PC_IS* pcis = (PC_IS*)pc->data;
123674ae819SStefano Zampini 
124674ae819SStefano Zampini   PetscFunctionBegin;
125674ae819SStefano Zampini   ierr = VecScatterBegin(pcis->global_to_B,global_vector,local_interface_vector,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
126674ae819SStefano Zampini   ierr = VecScatterEnd(pcis->global_to_B,global_vector,local_interface_vector,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
127674ae819SStefano Zampini   /* Apply partition of unity */
128674ae819SStefano Zampini   ierr = VecPointwiseMult(local_interface_vector,pcis->D,local_interface_vector);CHKERRQ(ierr);
129674ae819SStefano Zampini   PetscFunctionReturn(0);
130674ae819SStefano Zampini }
131674ae819SStefano Zampini 
132674ae819SStefano Zampini #undef __FUNCT__
133674ae819SStefano Zampini #define __FUNCT__ "PCBDDCScalingRestriction_Deluxe"
134674ae819SStefano Zampini static PetscErrorCode PCBDDCScalingRestriction_Deluxe(PC pc, Vec x, Vec y)
135674ae819SStefano Zampini {
136674ae819SStefano Zampini   PC_IS*              pcis=(PC_IS*)pc->data;
137674ae819SStefano Zampini   PC_BDDC*            pcbddc=(PC_BDDC*)pc->data;
138674ae819SStefano Zampini   PCBDDCDeluxeScaling deluxe_ctx = pcbddc->deluxe_ctx;
139b96c3477SStefano Zampini   PCBDDCSubSchurs     sub_schurs = pcbddc->sub_schurs;
140674ae819SStefano Zampini   PetscInt            i;
141674ae819SStefano Zampini   PetscErrorCode      ierr;
142674ae819SStefano Zampini 
143674ae819SStefano Zampini   PetscFunctionBegin;
144674ae819SStefano Zampini   /* get local boundary part of global vector */
145674ae819SStefano Zampini   ierr = VecScatterBegin(pcis->global_to_B,x,y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
146674ae819SStefano Zampini   ierr = VecScatterEnd(pcis->global_to_B,x,y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
14734a97f8cSStefano Zampini   if (deluxe_ctx->n_simple) {
14834a97f8cSStefano Zampini     /* scale deluxe vertices using diagonal scaling */
149*2b095fd8SStefano Zampini     PetscScalar       *array_y;
150*2b095fd8SStefano Zampini     const PetscScalar *array_D;
151674ae819SStefano Zampini     ierr = VecGetArray(y,&array_y);CHKERRQ(ierr);
152*2b095fd8SStefano Zampini     ierr = VecGetArrayRead(pcis->D,&array_D);CHKERRQ(ierr);
153674ae819SStefano Zampini     for (i=0;i<deluxe_ctx->n_simple;i++) {
154674ae819SStefano Zampini       array_y[deluxe_ctx->idx_simple_B[i]] *= array_D[deluxe_ctx->idx_simple_B[i]];
155674ae819SStefano Zampini     }
156*2b095fd8SStefano Zampini     ierr = VecRestoreArrayRead(pcis->D,&array_D);CHKERRQ(ierr);
157674ae819SStefano Zampini     ierr = VecRestoreArray(y,&array_y);CHKERRQ(ierr);
15834a97f8cSStefano Zampini   }
15934a97f8cSStefano Zampini   /* sequential part : all problems and Schur applications collapsed into a single matrix vector multiplication and ksp solution */
16041c3ba1bSStefano Zampini   if (deluxe_ctx->seq_ksp) {
161674ae819SStefano Zampini     ierr = VecScatterBegin(deluxe_ctx->seq_scctx,y,deluxe_ctx->seq_work1,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
162674ae819SStefano Zampini     ierr = VecScatterEnd(deluxe_ctx->seq_scctx,y,deluxe_ctx->seq_work1,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1631d850880SStefano Zampini     ierr = KSPSolve(deluxe_ctx->seq_ksp,deluxe_ctx->seq_work1,deluxe_ctx->seq_work2);CHKERRQ(ierr);
1641d850880SStefano Zampini     ierr = MatMult(sub_schurs->S_Ej_all,deluxe_ctx->seq_work2,deluxe_ctx->seq_work1);CHKERRQ(ierr);
165674ae819SStefano Zampini     ierr = VecScatterBegin(deluxe_ctx->seq_scctx,deluxe_ctx->seq_work1,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
166674ae819SStefano Zampini     ierr = VecScatterEnd(deluxe_ctx->seq_scctx,deluxe_ctx->seq_work1,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
167674ae819SStefano Zampini   }
168674ae819SStefano Zampini   /* parallel part */
169674ae819SStefano Zampini   for (i=0;i<deluxe_ctx->par_colors;i++) {
170674ae819SStefano Zampini     if (deluxe_ctx->par_ksp[i]) {
17134a97f8cSStefano Zampini       PetscInt subidx = deluxe_ctx->par_col2sub[i];
17234a97f8cSStefano Zampini       /* restrict on subset */
173b96c3477SStefano Zampini       ierr = VecScatterBegin(deluxe_ctx->par_scctx_s[i],y,deluxe_ctx->par_work1[i],INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
174b96c3477SStefano Zampini       ierr = VecScatterEnd(deluxe_ctx->par_scctx_s[i],y,deluxe_ctx->par_work1[i],INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
17534a97f8cSStefano Zampini       /* (\sum_j S_Ej)^-T */
17634a97f8cSStefano Zampini       ierr = VecSet(deluxe_ctx->par_vec[i],0.0);CHKERRQ(ierr);
177b96c3477SStefano Zampini       ierr = VecScatterBegin(deluxe_ctx->par_scctx_p[i],deluxe_ctx->par_work1[i],deluxe_ctx->par_vec[i],INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
178b96c3477SStefano Zampini       ierr = VecScatterEnd(deluxe_ctx->par_scctx_p[i],deluxe_ctx->par_work1[i],deluxe_ctx->par_vec[i],INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
179674ae819SStefano Zampini       ierr = KSPSolveTranspose(deluxe_ctx->par_ksp[i],deluxe_ctx->par_vec[i],deluxe_ctx->par_vec[i]);CHKERRQ(ierr);
180b96c3477SStefano Zampini       ierr = VecScatterBegin(deluxe_ctx->par_scctx_p[i],deluxe_ctx->par_vec[i],deluxe_ctx->par_work1[i],INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
181b96c3477SStefano Zampini       ierr = VecScatterEnd(deluxe_ctx->par_scctx_p[i],deluxe_ctx->par_vec[i],deluxe_ctx->par_work1[i],INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
18234a97f8cSStefano Zampini       /* S_Ej^T */
183b96c3477SStefano Zampini       ierr = MatMultTranspose(sub_schurs->S_Ej[subidx],deluxe_ctx->par_work1[i],deluxe_ctx->par_work2[i]);CHKERRQ(ierr);
18434a97f8cSStefano Zampini       /* extend to boundary */
185b96c3477SStefano Zampini       ierr = VecScatterBegin(deluxe_ctx->par_scctx_s[i],deluxe_ctx->par_work2[i],y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
186b96c3477SStefano Zampini       ierr = VecScatterEnd(deluxe_ctx->par_scctx_s[i],deluxe_ctx->par_work2[i],y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
187674ae819SStefano Zampini     }
188674ae819SStefano Zampini   }
189674ae819SStefano Zampini   PetscFunctionReturn(0);
190674ae819SStefano Zampini }
191674ae819SStefano Zampini 
192674ae819SStefano Zampini #undef __FUNCT__
193674ae819SStefano Zampini #define __FUNCT__ "PCBDDCScalingRestriction"
194674ae819SStefano Zampini PetscErrorCode PCBDDCScalingRestriction(PC pc, Vec global_vector, Vec local_interface_vector)
195674ae819SStefano Zampini {
196674ae819SStefano Zampini   PC_BDDC        *pcbddc=(PC_BDDC*)pc->data;
197674ae819SStefano Zampini   PetscErrorCode ierr;
198674ae819SStefano Zampini 
199674ae819SStefano Zampini   PetscFunctionBegin;
200674ae819SStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
201674ae819SStefano Zampini   PetscValidHeaderSpecific(global_vector,VEC_CLASSID,2);
202674ae819SStefano Zampini   PetscValidHeaderSpecific(local_interface_vector,VEC_CLASSID,3);
203674ae819SStefano Zampini   if (local_interface_vector == pcbddc->work_scaling) {
204674ae819SStefano Zampini     SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"Local vector should cannot be pcbddc->work_scaling!\n");
205674ae819SStefano Zampini   }
206674ae819SStefano Zampini   ierr = PetscTryMethod(pc,"PCBDDCScalingRestriction_C",(PC,Vec,Vec),(pc,global_vector,local_interface_vector));CHKERRQ(ierr);
207674ae819SStefano Zampini   PetscFunctionReturn(0);
208674ae819SStefano Zampini }
209674ae819SStefano Zampini 
210674ae819SStefano Zampini #undef __FUNCT__
211674ae819SStefano Zampini #define __FUNCT__ "PCBDDCScalingSetUp"
212674ae819SStefano Zampini PetscErrorCode PCBDDCScalingSetUp(PC pc)
213674ae819SStefano Zampini {
214674ae819SStefano Zampini   PC_IS* pcis=(PC_IS*)pc->data;
215674ae819SStefano Zampini   PC_BDDC* pcbddc=(PC_BDDC*)pc->data;
216674ae819SStefano Zampini   PetscErrorCode ierr;
217674ae819SStefano Zampini 
218674ae819SStefano Zampini   PetscFunctionBegin;
219674ae819SStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
220674ae819SStefano Zampini   /* create work vector for the operator */
22134a97f8cSStefano Zampini   ierr = VecDestroy(&pcbddc->work_scaling);CHKERRQ(ierr);
222674ae819SStefano Zampini   ierr = VecDuplicate(pcis->vec1_B,&pcbddc->work_scaling);CHKERRQ(ierr);
22334a97f8cSStefano Zampini   /* always rebuild pcis->D */
22428d874f6SStefano Zampini   if (pcis->use_stiffness_scaling) {
225674ae819SStefano Zampini     ierr = MatGetDiagonal(pcbddc->local_mat,pcis->vec1_N);CHKERRQ(ierr);
226674ae819SStefano Zampini     ierr = VecScatterBegin(pcis->N_to_B,pcis->vec1_N,pcis->D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
227674ae819SStefano Zampini     ierr = VecScatterEnd(pcis->N_to_B,pcis->vec1_N,pcis->D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
228674ae819SStefano Zampini   }
229674ae819SStefano Zampini   ierr = VecCopy(pcis->D,pcis->vec1_B);CHKERRQ(ierr);
230674ae819SStefano Zampini   ierr = VecSet(pcis->vec1_global,0.0);CHKERRQ(ierr);
231674ae819SStefano Zampini   ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_B,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
232674ae819SStefano Zampini   ierr = VecScatterEnd(pcis->global_to_B,pcis->vec1_B,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
233674ae819SStefano Zampini   ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_global,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
234674ae819SStefano Zampini   ierr = VecScatterEnd(pcis->global_to_B,pcis->vec1_global,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
235674ae819SStefano Zampini   ierr = VecPointwiseDivide(pcis->D,pcis->D,pcis->vec1_B);CHKERRQ(ierr);
236674ae819SStefano Zampini   /* now setup */
237681e7c04SStefano Zampini   if (pcbddc->use_deluxe_scaling) {
23834a97f8cSStefano Zampini     if (!pcbddc->deluxe_ctx) {
23934a97f8cSStefano Zampini       ierr = PCBDDCScalingCreate_Deluxe(pc);CHKERRQ(ierr);
24034a97f8cSStefano Zampini     }
24134a97f8cSStefano Zampini     ierr = PCBDDCScalingSetUp_Deluxe(pc);CHKERRQ(ierr);
242674ae819SStefano Zampini     ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCScalingRestriction_C",PCBDDCScalingRestriction_Deluxe);CHKERRQ(ierr);
243674ae819SStefano Zampini     ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCScalingExtension_C",PCBDDCScalingExtension_Deluxe);CHKERRQ(ierr);
244674ae819SStefano Zampini   } else {
245674ae819SStefano Zampini     ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCScalingRestriction_C",PCBDDCScalingRestriction_Basic);CHKERRQ(ierr);
246674ae819SStefano Zampini     ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCScalingExtension_C",PCBDDCScalingExtension_Basic);CHKERRQ(ierr);
247674ae819SStefano Zampini   }
24834a97f8cSStefano Zampini 
249674ae819SStefano Zampini   /* test */
250674ae819SStefano Zampini   if (pcbddc->dbg_flag) {
25134a97f8cSStefano Zampini     Vec         vec2_global;
252674ae819SStefano Zampini     PetscViewer viewer=pcbddc->dbg_viewer;
25334a97f8cSStefano Zampini     PetscReal   error;
254674ae819SStefano Zampini 
255674ae819SStefano Zampini     /* extension -> from local to parallel */
25634a97f8cSStefano Zampini     ierr = VecSet(pcis->vec1_global,0.0);CHKERRQ(ierr);
25734a97f8cSStefano Zampini     ierr = VecSetRandom(pcis->vec1_B,NULL);CHKERRQ(ierr);
25834a97f8cSStefano Zampini     ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_B,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
25934a97f8cSStefano Zampini     ierr = VecScatterEnd(pcis->global_to_B,pcis->vec1_B,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
26034a97f8cSStefano Zampini     ierr = VecDuplicate(pcis->vec1_global,&vec2_global);CHKERRQ(ierr);
26134a97f8cSStefano Zampini     ierr = VecCopy(pcis->vec1_global,vec2_global);CHKERRQ(ierr);
26234a97f8cSStefano Zampini 
263674ae819SStefano Zampini     ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_global,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
264674ae819SStefano Zampini     ierr = VecScatterEnd(pcis->global_to_B,pcis->vec1_global,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
265674ae819SStefano Zampini     ierr = PCBDDCScalingExtension(pc,pcis->vec1_B,pcis->vec1_global);CHKERRQ(ierr);
26634a97f8cSStefano Zampini     ierr = VecAXPY(pcis->vec1_global,-1.0,vec2_global);CHKERRQ(ierr);
26734a97f8cSStefano Zampini     ierr = VecNorm(pcis->vec1_global,NORM_INFINITY,&error);CHKERRQ(ierr);
268674ae819SStefano Zampini     ierr = PetscViewerASCIIPrintf(viewer,"Error scaling extension %1.14e\n",error);CHKERRQ(ierr);
26934a97f8cSStefano Zampini     if (error>1.e-8 && pcbddc->dbg_flag>1) {
270674ae819SStefano Zampini       ierr = VecView(pcis->vec1_global,viewer);CHKERRQ(ierr);
271674ae819SStefano Zampini     }
27234a97f8cSStefano Zampini     ierr = VecDestroy(&vec2_global);CHKERRQ(ierr);
27334a97f8cSStefano Zampini 
274674ae819SStefano Zampini     /* restriction -> from parallel to local */
275674ae819SStefano Zampini     ierr = VecSet(pcis->vec1_global,0.0);CHKERRQ(ierr);
27634a97f8cSStefano Zampini     ierr = VecSetRandom(pcis->vec1_B,NULL);CHKERRQ(ierr);
277674ae819SStefano Zampini     ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_B,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
278674ae819SStefano Zampini     ierr = VecScatterEnd(pcis->global_to_B,pcis->vec1_B,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
27934a97f8cSStefano Zampini 
28034a97f8cSStefano Zampini     ierr = PCBDDCScalingRestriction(pc,pcis->vec1_global,pcis->vec1_B);CHKERRQ(ierr);
28134a97f8cSStefano Zampini     ierr = VecScale(pcis->vec1_B,-1.0);CHKERRQ(ierr);
28234a97f8cSStefano Zampini     ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_B,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
28334a97f8cSStefano Zampini     ierr = VecScatterEnd(pcis->global_to_B,pcis->vec1_B,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
28434a97f8cSStefano Zampini     ierr = VecNorm(pcis->vec1_global,NORM_INFINITY,&error);CHKERRQ(ierr);
28534a97f8cSStefano Zampini     ierr = PetscViewerASCIIPrintf(viewer,"Error scaling restriction %1.14e\n",error);CHKERRQ(ierr);
28634a97f8cSStefano Zampini     if (error>1.e-8 && pcbddc->dbg_flag>1) {
287674ae819SStefano Zampini       ierr = VecView(pcis->vec1_global,viewer);CHKERRQ(ierr);
288674ae819SStefano Zampini     }
289674ae819SStefano Zampini   }
290674ae819SStefano Zampini   PetscFunctionReturn(0);
291674ae819SStefano Zampini }
292674ae819SStefano Zampini 
293674ae819SStefano Zampini #undef __FUNCT__
294674ae819SStefano Zampini #define __FUNCT__ "PCBDDCScalingDestroy"
295674ae819SStefano Zampini PetscErrorCode PCBDDCScalingDestroy(PC pc)
296674ae819SStefano Zampini {
297674ae819SStefano Zampini   PC_BDDC* pcbddc=(PC_BDDC*)pc->data;
298674ae819SStefano Zampini   PetscErrorCode ierr;
299674ae819SStefano Zampini 
300674ae819SStefano Zampini   PetscFunctionBegin;
30134a97f8cSStefano Zampini   if (pcbddc->deluxe_ctx) {
30234a97f8cSStefano Zampini     ierr = PCBDDCScalingDestroy_Deluxe(pc);CHKERRQ(ierr);
303674ae819SStefano Zampini   }
304674ae819SStefano Zampini   ierr = VecDestroy(&pcbddc->work_scaling);CHKERRQ(ierr);
305674ae819SStefano Zampini   /* remove functions */
306674ae819SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCScalingRestriction_C",NULL);CHKERRQ(ierr);
307674ae819SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCScalingExtension_C",NULL);CHKERRQ(ierr);
308674ae819SStefano Zampini   PetscFunctionReturn(0);
309674ae819SStefano Zampini }
310674ae819SStefano Zampini 
31134a97f8cSStefano Zampini #undef __FUNCT__
31234a97f8cSStefano Zampini #define __FUNCT__ "PCBDDCScalingCreate_Deluxe"
31334a97f8cSStefano Zampini static PetscErrorCode PCBDDCScalingCreate_Deluxe(PC pc)
31434a97f8cSStefano Zampini {
31534a97f8cSStefano Zampini   PC_BDDC*            pcbddc=(PC_BDDC*)pc->data;
31634a97f8cSStefano Zampini   PCBDDCDeluxeScaling deluxe_ctx;
31734a97f8cSStefano Zampini   PetscErrorCode      ierr;
31834a97f8cSStefano Zampini 
31934a97f8cSStefano Zampini   PetscFunctionBegin;
32034a97f8cSStefano Zampini   ierr = PetscNew(&deluxe_ctx);CHKERRQ(ierr);
32134a97f8cSStefano Zampini   pcbddc->deluxe_ctx = deluxe_ctx;
32234a97f8cSStefano Zampini   PetscFunctionReturn(0);
32334a97f8cSStefano Zampini }
32434a97f8cSStefano Zampini 
32534a97f8cSStefano Zampini #undef __FUNCT__
32634a97f8cSStefano Zampini #define __FUNCT__ "PCBDDCScalingDestroy_Deluxe"
32734a97f8cSStefano Zampini static PetscErrorCode PCBDDCScalingDestroy_Deluxe(PC pc)
32834a97f8cSStefano Zampini {
32934a97f8cSStefano Zampini   PC_BDDC*            pcbddc=(PC_BDDC*)pc->data;
33034a97f8cSStefano Zampini   PetscErrorCode      ierr;
33134a97f8cSStefano Zampini 
33234a97f8cSStefano Zampini   PetscFunctionBegin;
33334a97f8cSStefano Zampini   ierr = PCBDDCScalingReset_Deluxe_Solvers(pcbddc->deluxe_ctx);CHKERRQ(ierr);
33434a97f8cSStefano Zampini   ierr = PetscFree(pcbddc->deluxe_ctx);CHKERRQ(ierr);
33534a97f8cSStefano Zampini   PetscFunctionReturn(0);
33634a97f8cSStefano Zampini }
33734a97f8cSStefano Zampini 
33834a97f8cSStefano Zampini #undef __FUNCT__
33934a97f8cSStefano Zampini #define __FUNCT__ "PCBDDCScalingReset_Deluxe_Solvers"
34034a97f8cSStefano Zampini static PetscErrorCode PCBDDCScalingReset_Deluxe_Solvers(PCBDDCDeluxeScaling deluxe_ctx)
34134a97f8cSStefano Zampini {
34234a97f8cSStefano Zampini   PetscErrorCode      ierr;
34334a97f8cSStefano Zampini 
34434a97f8cSStefano Zampini   PetscFunctionBegin;
34534a97f8cSStefano Zampini   ierr = PetscFree(deluxe_ctx->idx_simple_B);CHKERRQ(ierr);
34634a97f8cSStefano Zampini   deluxe_ctx->n_simple = 0;
34741c3ba1bSStefano Zampini   if (deluxe_ctx->seq_ksp) {
34834a97f8cSStefano Zampini     ierr = VecScatterDestroy(&deluxe_ctx->seq_scctx);CHKERRQ(ierr);
34934a97f8cSStefano Zampini     ierr = VecDestroy(&deluxe_ctx->seq_work1);CHKERRQ(ierr);
35034a97f8cSStefano Zampini     ierr = VecDestroy(&deluxe_ctx->seq_work2);CHKERRQ(ierr);
35134a97f8cSStefano Zampini     ierr = KSPDestroy(&deluxe_ctx->seq_ksp);CHKERRQ(ierr);
35234a97f8cSStefano Zampini   }
35334a97f8cSStefano Zampini   if (deluxe_ctx->par_colors) {
35434a97f8cSStefano Zampini     PetscInt i;
35534a97f8cSStefano Zampini     for (i=0;i<deluxe_ctx->par_colors;i++) {
35634a97f8cSStefano Zampini       ierr = VecScatterDestroy(&deluxe_ctx->par_scctx_s[i]);CHKERRQ(ierr);
35734a97f8cSStefano Zampini       ierr = VecScatterDestroy(&deluxe_ctx->par_scctx_p[i]);CHKERRQ(ierr);
35834a97f8cSStefano Zampini       ierr = VecDestroy(&deluxe_ctx->par_vec[i]);CHKERRQ(ierr);
359b96c3477SStefano Zampini       ierr = VecDestroy(&deluxe_ctx->par_work1[i]);CHKERRQ(ierr);
360b96c3477SStefano Zampini       ierr = VecDestroy(&deluxe_ctx->par_work2[i]);CHKERRQ(ierr);
36134a97f8cSStefano Zampini       ierr = KSPDestroy(&deluxe_ctx->par_ksp[i]);CHKERRQ(ierr);
36234a97f8cSStefano Zampini     }
363b96c3477SStefano Zampini     ierr = PetscFree7(deluxe_ctx->par_ksp,
36434a97f8cSStefano Zampini                       deluxe_ctx->par_scctx_s,
36534a97f8cSStefano Zampini                       deluxe_ctx->par_scctx_p,
36634a97f8cSStefano Zampini                       deluxe_ctx->par_vec,
367b96c3477SStefano Zampini                       deluxe_ctx->par_work1,
368b96c3477SStefano Zampini                       deluxe_ctx->par_work2,
36934a97f8cSStefano Zampini                       deluxe_ctx->par_col2sub);CHKERRQ(ierr);
37034a97f8cSStefano Zampini   }
37134a97f8cSStefano Zampini   deluxe_ctx->par_colors = 0;
37234a97f8cSStefano Zampini   PetscFunctionReturn(0);
37334a97f8cSStefano Zampini }
37434a97f8cSStefano Zampini 
37534a97f8cSStefano Zampini #undef __FUNCT__
37634a97f8cSStefano Zampini #define __FUNCT__ "PCBDDCScalingSetUp_Deluxe"
37734a97f8cSStefano Zampini static PetscErrorCode PCBDDCScalingSetUp_Deluxe(PC pc)
37834a97f8cSStefano Zampini {
37934a97f8cSStefano Zampini   PC_IS               *pcis=(PC_IS*)pc->data;
38034a97f8cSStefano Zampini   PC_BDDC             *pcbddc=(PC_BDDC*)pc->data;
38134a97f8cSStefano Zampini   PCBDDCDeluxeScaling deluxe_ctx=pcbddc->deluxe_ctx;
382b96c3477SStefano Zampini   PCBDDCSubSchurs     sub_schurs=pcbddc->sub_schurs;
38334a97f8cSStefano Zampini   PetscErrorCode      ierr;
38434a97f8cSStefano Zampini 
38534a97f8cSStefano Zampini   PetscFunctionBegin;
386b96c3477SStefano Zampini   /* (TODO: reuse) throw away the solvers */
38734a97f8cSStefano Zampini   ierr = PCBDDCScalingReset_Deluxe_Solvers(deluxe_ctx);CHKERRQ(ierr);
38834a97f8cSStefano Zampini 
389b1b3d7a2SStefano Zampini   /* Compute data structures to solve parallel problems */
390b1b3d7a2SStefano Zampini   ierr = PCBDDCScalingSetUp_Deluxe_Par(pc,sub_schurs->n_subs_par,sub_schurs->n_subs_par_g,
391b1b3d7a2SStefano Zampini                                        sub_schurs->auxglobal_parallel,
392b1b3d7a2SStefano Zampini                                        sub_schurs->index_parallel);CHKERRQ(ierr);
393b96c3477SStefano Zampini 
394b1b3d7a2SStefano Zampini   /* Compute data structures to solve sequential problems */
395883469d8SStefano Zampini   ierr = PCBDDCScalingSetUp_Deluxe_Seq(pc);CHKERRQ(ierr);
3965db18549SStefano Zampini 
397b1b3d7a2SStefano Zampini   /* diagonal scaling on interface dofs not contained in cc */
3989bb4a8caSStefano Zampini   if (sub_schurs->is_Ej_com) {
3999bb4a8caSStefano Zampini     PetscInt       nmap;
4009bb4a8caSStefano Zampini     const PetscInt *idxs;
4019bb4a8caSStefano Zampini     ierr = ISGetLocalSize(sub_schurs->is_Ej_com,&deluxe_ctx->n_simple);CHKERRQ(ierr);
4029bb4a8caSStefano Zampini     ierr = ISGetIndices(sub_schurs->is_Ej_com,&idxs);CHKERRQ(ierr);
403b1b3d7a2SStefano Zampini     ierr = PetscMalloc1(deluxe_ctx->n_simple,&deluxe_ctx->idx_simple_B);CHKERRQ(ierr);
4049bb4a8caSStefano Zampini     ierr = ISGlobalToLocalMappingApply(pcis->BtoNmap,IS_GTOLM_DROP,deluxe_ctx->n_simple,idxs,&nmap,deluxe_ctx->idx_simple_B);CHKERRQ(ierr);
4059bb4a8caSStefano Zampini     if (nmap != deluxe_ctx->n_simple) {
4069bb4a8caSStefano Zampini       SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Error when mapping simply scaled dofs! %d != %d",nmap,deluxe_ctx->n_simple);
4079bb4a8caSStefano Zampini     }
4089bb4a8caSStefano Zampini     ierr = ISRestoreIndices(sub_schurs->is_Ej_com,&idxs);CHKERRQ(ierr);
4099bb4a8caSStefano Zampini   } else {
410b1b3d7a2SStefano Zampini     deluxe_ctx->n_simple = 0;
4119bb4a8caSStefano Zampini     deluxe_ctx->idx_simple_B = 0;
412b1b3d7a2SStefano Zampini   }
41334a97f8cSStefano Zampini   PetscFunctionReturn(0);
41434a97f8cSStefano Zampini }
41534a97f8cSStefano Zampini 
41634a97f8cSStefano Zampini #undef __FUNCT__
41734a97f8cSStefano Zampini #define __FUNCT__ "PCBDDCScalingSetUp_Deluxe_Par"
41834a97f8cSStefano Zampini static PetscErrorCode PCBDDCScalingSetUp_Deluxe_Par(PC pc, PetscInt n_local_parallel_problems,PetscInt n_parallel_problems,PetscInt global_parallel[],PetscInt index_parallel[])
41934a97f8cSStefano Zampini {
42034a97f8cSStefano Zampini   PC_BDDC                *pcbddc=(PC_BDDC*)pc->data;
42134a97f8cSStefano Zampini   PCBDDCDeluxeScaling    deluxe_ctx=pcbddc->deluxe_ctx;
42234a97f8cSStefano Zampini   /* coloring */
42334a97f8cSStefano Zampini   Mat                    parallel_problems;
42434a97f8cSStefano Zampini   MatColoring            coloring_obj;
42534a97f8cSStefano Zampini   ISColoring             coloring_parallel_problems;
42634a97f8cSStefano Zampini   IS                     *par_is_colors,*is_colors;
42734a97f8cSStefano Zampini   /* working stuff */
42834a97f8cSStefano Zampini   PetscInt               i,j;
42934a97f8cSStefano Zampini   PetscErrorCode         ierr;
43034a97f8cSStefano Zampini 
43134a97f8cSStefano Zampini   PetscFunctionBegin;
43234a97f8cSStefano Zampini   if (!n_parallel_problems) {
43334a97f8cSStefano Zampini     PetscFunctionReturn(0);
43434a97f8cSStefano Zampini   }
43534a97f8cSStefano Zampini   /* Color parallel subproblems */
43634a97f8cSStefano Zampini   ierr = MatCreate(PetscObjectComm((PetscObject)pc),&parallel_problems);CHKERRQ(ierr);
43734a97f8cSStefano Zampini   ierr = MatSetSizes(parallel_problems,PETSC_DECIDE,PETSC_DECIDE,n_parallel_problems,n_parallel_problems);CHKERRQ(ierr);
43834a97f8cSStefano Zampini   ierr = MatSetType(parallel_problems,MATAIJ);CHKERRQ(ierr);
43934a97f8cSStefano Zampini   ierr = MatSetUp(parallel_problems);CHKERRQ(ierr);
44034a97f8cSStefano Zampini   ierr = MatSetOption(parallel_problems,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_FALSE);CHKERRQ(ierr);
44134a97f8cSStefano Zampini   ierr = MatSetOption(parallel_problems,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_FALSE);CHKERRQ(ierr);
44234a97f8cSStefano Zampini   for (i=0;i<n_local_parallel_problems;i++) {
44334a97f8cSStefano Zampini     PetscInt row = global_parallel[i];
44434a97f8cSStefano Zampini     for (j=0;j<n_local_parallel_problems;j++) {
44534a97f8cSStefano Zampini       PetscInt col = global_parallel[j];
44634a97f8cSStefano Zampini       if (row != col) {
44734a97f8cSStefano Zampini         ierr = MatSetValue(parallel_problems,row,col,1.0,INSERT_VALUES);CHKERRQ(ierr);
44834a97f8cSStefano Zampini       }
44934a97f8cSStefano Zampini     }
45034a97f8cSStefano Zampini   }
45134a97f8cSStefano Zampini   ierr = MatAssemblyBegin(parallel_problems,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
45234a97f8cSStefano Zampini   ierr = MatAssemblyEnd(parallel_problems,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
45334a97f8cSStefano Zampini   if (pcbddc->dbg_flag > 1) {
45434a97f8cSStefano Zampini     ierr = PetscViewerFlush(pcbddc->dbg_viewer);CHKERRQ(ierr);
45534a97f8cSStefano Zampini     ierr = PetscViewerASCIIPrintf(pcbddc->dbg_viewer,"Adj matrix for deluxe parallel problems\n");CHKERRQ(ierr);
45634a97f8cSStefano Zampini     ierr = MatView(parallel_problems,pcbddc->dbg_viewer);CHKERRQ(ierr);
45734a97f8cSStefano Zampini   }
45834a97f8cSStefano Zampini   ierr = MatColoringCreate(parallel_problems,&coloring_obj);CHKERRQ(ierr);
45934a97f8cSStefano Zampini   ierr = MatColoringSetDistance(coloring_obj,1);CHKERRQ(ierr);
46034a97f8cSStefano Zampini   ierr = MatColoringSetType(coloring_obj,MATCOLORINGJP);CHKERRQ(ierr);
46134a97f8cSStefano Zampini   ierr = MatColoringApply(coloring_obj,&coloring_parallel_problems);CHKERRQ(ierr);
46234a97f8cSStefano Zampini   ierr = ISColoringGetIS(coloring_parallel_problems,&deluxe_ctx->par_colors,&par_is_colors);CHKERRQ(ierr);
46334a97f8cSStefano Zampini   if (pcbddc->dbg_flag) {
46434a97f8cSStefano Zampini     ierr = PetscViewerFlush(pcbddc->dbg_viewer);CHKERRQ(ierr);
46534a97f8cSStefano Zampini     ierr = PetscViewerASCIIPrintf(pcbddc->dbg_viewer,"Number of colors %d for parallel part of deluxe\n",deluxe_ctx->par_colors);CHKERRQ(ierr);
46634a97f8cSStefano Zampini   }
46734a97f8cSStefano Zampini 
46834a97f8cSStefano Zampini   /* all procs should know the color distribution */
46934a97f8cSStefano Zampini   ierr = PetscMalloc1(deluxe_ctx->par_colors,&is_colors);CHKERRQ(ierr);
47034a97f8cSStefano Zampini   for (i=0;i<deluxe_ctx->par_colors;i++) {
47134a97f8cSStefano Zampini     if (pcbddc->dbg_flag) {
47234a97f8cSStefano Zampini       ierr = PetscViewerASCIIPrintf(pcbddc->dbg_viewer,"Global problem indexes for color %d\n",i);CHKERRQ(ierr);
47334a97f8cSStefano Zampini       ierr = ISView(par_is_colors[i],pcbddc->dbg_viewer);CHKERRQ(ierr);
47434a97f8cSStefano Zampini       ierr = PetscViewerFlush(pcbddc->dbg_viewer);CHKERRQ(ierr);
47534a97f8cSStefano Zampini     }
47634a97f8cSStefano Zampini     ierr = ISAllGather(par_is_colors[i],&is_colors[i]);CHKERRQ(ierr);
47734a97f8cSStefano Zampini   }
47834a97f8cSStefano Zampini 
47934a97f8cSStefano Zampini   /* free unneeded objects */
48034a97f8cSStefano Zampini   ierr = ISColoringRestoreIS(coloring_parallel_problems,&par_is_colors);CHKERRQ(ierr);
48134a97f8cSStefano Zampini   ierr = ISColoringDestroy(&coloring_parallel_problems);CHKERRQ(ierr);
48234a97f8cSStefano Zampini   ierr = MatColoringDestroy(&coloring_obj);CHKERRQ(ierr);
48334a97f8cSStefano Zampini   ierr = MatDestroy(&parallel_problems);CHKERRQ(ierr);
48434a97f8cSStefano Zampini 
48534a97f8cSStefano Zampini   /* allocate deluxe arrays for parallel problems */
486b96c3477SStefano Zampini   ierr = PetscCalloc7(deluxe_ctx->par_colors,&deluxe_ctx->par_ksp,
48734a97f8cSStefano Zampini                       deluxe_ctx->par_colors,&deluxe_ctx->par_scctx_s,
48834a97f8cSStefano Zampini                       deluxe_ctx->par_colors,&deluxe_ctx->par_scctx_p,
48934a97f8cSStefano Zampini                       deluxe_ctx->par_colors,&deluxe_ctx->par_vec,
490b96c3477SStefano Zampini                       deluxe_ctx->par_colors,&deluxe_ctx->par_work1,
491b96c3477SStefano Zampini                       deluxe_ctx->par_colors,&deluxe_ctx->par_work2,
49234a97f8cSStefano Zampini                       deluxe_ctx->par_colors,&deluxe_ctx->par_col2sub);CHKERRQ(ierr);
49334a97f8cSStefano Zampini 
49434a97f8cSStefano Zampini   /* cycle on colors */
49534a97f8cSStefano Zampini   for (i=0;i<deluxe_ctx->par_colors;i++) {
49634a97f8cSStefano Zampini     PetscSubcomm    par_subcomm;
49734a97f8cSStefano Zampini     const PetscInt* idxs_subproblems;
49834a97f8cSStefano Zampini     PetscInt        color_size;
49934a97f8cSStefano Zampini     PetscMPIInt     rank,active_color;
50034a97f8cSStefano Zampini 
50134a97f8cSStefano Zampini     /* get local index of i-th parallel colored problem */
50234a97f8cSStefano Zampini     ierr = ISGetLocalSize(is_colors[i],&color_size);CHKERRQ(ierr);
50334a97f8cSStefano Zampini     ierr = ISGetIndices(is_colors[i],&idxs_subproblems);CHKERRQ(ierr);
50434a97f8cSStefano Zampini     /* split comm for computing parallel problems for this color */
50534a97f8cSStefano Zampini     /* Processes not partecipating at this stage will have color = color_size */
50634a97f8cSStefano Zampini     /* because PetscCommDuplicate does not handle MPI_COMM_NULL */
50734a97f8cSStefano Zampini     active_color = color_size;
50834a97f8cSStefano Zampini     deluxe_ctx->par_col2sub[i] = -1;
50934a97f8cSStefano Zampini     for (j=0;j<n_local_parallel_problems;j++) {
51034a97f8cSStefano Zampini       PetscInt local_idx;
51134a97f8cSStefano Zampini       ierr = PetscFindInt(global_parallel[j],color_size,idxs_subproblems,&local_idx);CHKERRQ(ierr);
51234a97f8cSStefano Zampini       if (local_idx > -1) {
51334a97f8cSStefano Zampini         ierr = PetscMPIIntCast(local_idx,&active_color);CHKERRQ(ierr);
51434a97f8cSStefano Zampini         deluxe_ctx->par_col2sub[i] = index_parallel[j];
51534a97f8cSStefano Zampini         break;
51634a97f8cSStefano Zampini       }
51734a97f8cSStefano Zampini     }
51834a97f8cSStefano Zampini     ierr = ISRestoreIndices(is_colors[i],&idxs_subproblems);CHKERRQ(ierr);
51934a97f8cSStefano Zampini     ierr = PetscSubcommCreate(PetscObjectComm((PetscObject)pc),&par_subcomm);CHKERRQ(ierr);
52034a97f8cSStefano Zampini     ierr = PetscSubcommSetNumber(par_subcomm,color_size+1);CHKERRQ(ierr);
52134a97f8cSStefano Zampini     ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)pc),&rank);CHKERRQ(ierr);
52234a97f8cSStefano Zampini     ierr = PetscSubcommSetTypeGeneral(par_subcomm,active_color,rank);CHKERRQ(ierr);
52334a97f8cSStefano Zampini     /* print debug info */
52434a97f8cSStefano Zampini     if (pcbddc->dbg_flag) {
52534a97f8cSStefano Zampini       PetscMPIInt crank,csize;
52634a97f8cSStefano Zampini       ierr = MPI_Comm_rank(par_subcomm->comm,&crank);CHKERRQ(ierr);
52734a97f8cSStefano Zampini       ierr = MPI_Comm_size(par_subcomm->comm,&csize);CHKERRQ(ierr);
52834a97f8cSStefano Zampini       ierr = PetscViewerASCIIPrintf(pcbddc->dbg_viewer,"Color %d: size %d, details follows.\n",i,color_size);CHKERRQ(ierr);
52934a97f8cSStefano Zampini       ierr = PetscViewerFlush(pcbddc->dbg_viewer);CHKERRQ(ierr);
53034a97f8cSStefano Zampini       ierr = PetscViewerASCIISynchronizedAllow(pcbddc->dbg_viewer,PETSC_TRUE);CHKERRQ(ierr);
53134a97f8cSStefano Zampini       ierr = PetscViewerASCIISynchronizedPrintf(pcbddc->dbg_viewer,"  Subdomain %d: color in subcomm %d (rank %d out of %d) (lidx %d)\n",PetscGlobalRank,par_subcomm->color,crank,csize,deluxe_ctx->par_col2sub[i]);CHKERRQ(ierr);
53234a97f8cSStefano Zampini       ierr = PetscViewerFlush(pcbddc->dbg_viewer);CHKERRQ(ierr);
53334a97f8cSStefano Zampini     }
53434a97f8cSStefano Zampini 
53534a97f8cSStefano Zampini     if (deluxe_ctx->par_col2sub[i] >= 0) {
5365e8657edSStefano Zampini       PC                     pctemp;
5375e8657edSStefano Zampini       PC_IS                  *pcis=(PC_IS*)pc->data;
53834a97f8cSStefano Zampini       Mat                    color_mat,color_mat_is,temp_mat;
53934a97f8cSStefano Zampini       ISLocalToGlobalMapping WtoNmap,l2gmap_subset;
54034a97f8cSStefano Zampini       IS                     is_local_numbering,isB_local,isW_local,isW;
541b96c3477SStefano Zampini       PCBDDCSubSchurs        sub_schurs = pcbddc->sub_schurs;
54234a97f8cSStefano Zampini       PetscInt               subidx,n_local_dofs,n_global_dofs;
54334a97f8cSStefano Zampini       PetscInt               *global_numbering,*local_numbering;
54434a97f8cSStefano Zampini       char                   ksp_prefix[256];
54534a97f8cSStefano Zampini       size_t                 len;
54634a97f8cSStefano Zampini 
54734a97f8cSStefano Zampini       /* Local index for schur complement on subset */
54834a97f8cSStefano Zampini       subidx = deluxe_ctx->par_col2sub[i];
54934a97f8cSStefano Zampini 
55034a97f8cSStefano Zampini       /* Parallel numbering for dofs in colored subset */
551b96c3477SStefano Zampini       ierr = ISSum(sub_schurs->is_I_layer,sub_schurs->is_subs[subidx],&is_local_numbering);CHKERRQ(ierr);
55234a97f8cSStefano Zampini       ierr = ISGetLocalSize(is_local_numbering,&n_local_dofs);CHKERRQ(ierr);
55334a97f8cSStefano Zampini       ierr = ISGetIndices(is_local_numbering,(const PetscInt **)&local_numbering);CHKERRQ(ierr);
55434a97f8cSStefano Zampini       ierr = PCBDDCSubsetNumbering(par_subcomm->comm,pcbddc->mat_graph->l2gmap,n_local_dofs,local_numbering,PETSC_NULL,&n_global_dofs,&global_numbering);CHKERRQ(ierr);
55534a97f8cSStefano Zampini       ierr = ISRestoreIndices(is_local_numbering,(const PetscInt **)&local_numbering);CHKERRQ(ierr);
55634a97f8cSStefano Zampini 
55734a97f8cSStefano Zampini       /* L2Gmap from relevant dofs to local dofs */
55834a97f8cSStefano Zampini       ierr = ISLocalToGlobalMappingCreateIS(is_local_numbering,&WtoNmap);CHKERRQ(ierr);
55934a97f8cSStefano Zampini 
56034a97f8cSStefano Zampini       /* L2Gmap from local to global dofs */
56134a97f8cSStefano Zampini       ierr = ISLocalToGlobalMappingCreate(par_subcomm->comm,1,n_local_dofs,global_numbering,PETSC_COPY_VALUES,&l2gmap_subset);CHKERRQ(ierr);
56234a97f8cSStefano Zampini 
56334a97f8cSStefano Zampini       /* compute parallel matrix (extended dirichlet problem on subset) */
56434a97f8cSStefano Zampini       ierr = MatCreateIS(par_subcomm->comm,1,PETSC_DECIDE,PETSC_DECIDE,n_global_dofs,n_global_dofs,l2gmap_subset,&color_mat_is);CHKERRQ(ierr);
56534a97f8cSStefano Zampini       ierr = MatGetSubMatrix(pcbddc->local_mat,is_local_numbering,is_local_numbering,MAT_INITIAL_MATRIX,&temp_mat);CHKERRQ(ierr);
56634a97f8cSStefano Zampini       ierr = MatISSetLocalMat(color_mat_is,temp_mat);CHKERRQ(ierr);
56734a97f8cSStefano Zampini       ierr = MatDestroy(&temp_mat);CHKERRQ(ierr);
56834a97f8cSStefano Zampini       ierr = MatISGetMPIXAIJ(color_mat_is,MAT_INITIAL_MATRIX,&color_mat);CHKERRQ(ierr);
56934a97f8cSStefano Zampini       ierr = MatDestroy(&color_mat_is);CHKERRQ(ierr);
57034a97f8cSStefano Zampini 
57134a97f8cSStefano Zampini       /* work vector for (parallel) extended dirichlet problem */
5728a26ef87SStefano Zampini       ierr = MatCreateVecs(color_mat,&deluxe_ctx->par_vec[i],NULL);CHKERRQ(ierr);
57334a97f8cSStefano Zampini 
57434a97f8cSStefano Zampini       /* compute scatters */
57534a97f8cSStefano Zampini       /* deluxe_ctx->par_scctx_p[i] extension from local subset to extended dirichlet problem
57634a97f8cSStefano Zampini          deluxe_ctx->par_scctx_s[i] restriction from local boundary to subset -> simple copy of selected values */
577b96c3477SStefano Zampini       ierr = ISGlobalToLocalMappingApplyIS(pcis->BtoNmap,IS_GTOLM_DROP,sub_schurs->is_subs[subidx],&isB_local);CHKERRQ(ierr);
578b96c3477SStefano Zampini       ierr = MatCreateVecs(sub_schurs->S_Ej[subidx],&deluxe_ctx->par_work1[i],&deluxe_ctx->par_work2[i]);CHKERRQ(ierr);
579b96c3477SStefano Zampini       ierr = VecScatterCreate(pcbddc->work_scaling,isB_local,deluxe_ctx->par_work1[i],NULL,&deluxe_ctx->par_scctx_s[i]);CHKERRQ(ierr);
580b96c3477SStefano Zampini       ierr = ISGlobalToLocalMappingApplyIS(WtoNmap,IS_GTOLM_DROP,sub_schurs->is_subs[subidx],&isW_local);CHKERRQ(ierr);
58134a97f8cSStefano Zampini       ierr = ISLocalToGlobalMappingApplyIS(l2gmap_subset,isW_local,&isW);CHKERRQ(ierr);
582b96c3477SStefano Zampini       ierr = VecScatterCreate(deluxe_ctx->par_work1[i],NULL,deluxe_ctx->par_vec[i],isW,&deluxe_ctx->par_scctx_p[i]);CHKERRQ(ierr);
58334a97f8cSStefano Zampini 
58434a97f8cSStefano Zampini       /* free objects no longer neeeded */
58534a97f8cSStefano Zampini       ierr = ISDestroy(&isW);CHKERRQ(ierr);
58634a97f8cSStefano Zampini       ierr = ISDestroy(&isW_local);CHKERRQ(ierr);
58734a97f8cSStefano Zampini       ierr = ISDestroy(&isB_local);CHKERRQ(ierr);
58834a97f8cSStefano Zampini       ierr = ISLocalToGlobalMappingDestroy(&WtoNmap);CHKERRQ(ierr);
58934a97f8cSStefano Zampini       ierr = ISLocalToGlobalMappingDestroy(&l2gmap_subset);CHKERRQ(ierr);
59034a97f8cSStefano Zampini       ierr = ISDestroy(&is_local_numbering);CHKERRQ(ierr);
59134a97f8cSStefano Zampini       ierr = PetscFree(global_numbering);CHKERRQ(ierr);
59234a97f8cSStefano Zampini 
59334a97f8cSStefano Zampini       /* KSP for extended dirichlet problem */
59434a97f8cSStefano Zampini       ierr = KSPCreate(par_subcomm->comm,&deluxe_ctx->par_ksp[i]);CHKERRQ(ierr);
59534a97f8cSStefano Zampini       ierr = KSPSetOperators(deluxe_ctx->par_ksp[i],color_mat,color_mat);CHKERRQ(ierr);
59634a97f8cSStefano Zampini       ierr = KSPSetTolerances(deluxe_ctx->par_ksp[i],1.e-12,1.e-12,1.e10,10000);CHKERRQ(ierr);
59734a97f8cSStefano Zampini       ierr = KSPSetType(deluxe_ctx->par_ksp[i],KSPPREONLY);CHKERRQ(ierr);
5985e8657edSStefano Zampini       ierr = KSPGetPC(deluxe_ctx->par_ksp[i],&pctemp);CHKERRQ(ierr);
5995e8657edSStefano Zampini       ierr = PCSetType(pctemp,PCREDUNDANT);CHKERRQ(ierr);
60034a97f8cSStefano Zampini       ierr = PetscStrlen(((PetscObject)(pcbddc->ksp_D))->prefix,&len);CHKERRQ(ierr);
6018856534dSStefano Zampini       len -= 10; /* remove "dirichlet_" */
6028856534dSStefano Zampini       ierr = PetscStrncpy(ksp_prefix,((PetscObject)(pcbddc->ksp_D))->prefix,len+1);CHKERRQ(ierr); /* PetscStrncpy puts a terminating char at the end */
60334a97f8cSStefano Zampini       ierr = PetscStrcat(ksp_prefix,"deluxe_par_");CHKERRQ(ierr);
60434a97f8cSStefano Zampini       ierr = KSPSetOptionsPrefix(deluxe_ctx->par_ksp[i],ksp_prefix);CHKERRQ(ierr);
60534a97f8cSStefano Zampini       ierr = KSPSetFromOptions(deluxe_ctx->par_ksp[i]);CHKERRQ(ierr);
60634a97f8cSStefano Zampini       ierr = KSPSetUp(deluxe_ctx->par_ksp[i]);CHKERRQ(ierr);
60734a97f8cSStefano Zampini       ierr = MatDestroy(&color_mat);CHKERRQ(ierr);
60834a97f8cSStefano Zampini     }
60934a97f8cSStefano Zampini     ierr = PetscSubcommDestroy(&par_subcomm);CHKERRQ(ierr);
61034a97f8cSStefano Zampini   }
61134a97f8cSStefano Zampini   for (i=0;i<deluxe_ctx->par_colors;i++) {
61234a97f8cSStefano Zampini     ierr = ISDestroy(&is_colors[i]);CHKERRQ(ierr);
61334a97f8cSStefano Zampini   }
61434a97f8cSStefano Zampini   ierr = PetscFree(is_colors);CHKERRQ(ierr);
61534a97f8cSStefano Zampini 
61634a97f8cSStefano Zampini   if (pcbddc->dbg_flag) {
61734a97f8cSStefano Zampini     Vec test_vec;
61834a97f8cSStefano Zampini     PetscReal error;
619b96c3477SStefano Zampini     PCBDDCSubSchurs sub_schurs = pcbddc->sub_schurs;
62034a97f8cSStefano Zampini     /* test partition of unity of coloured schur complements  */
62134a97f8cSStefano Zampini     for (i=0;i<deluxe_ctx->par_colors;i++) {
62234a97f8cSStefano Zampini       PetscInt  subidx = deluxe_ctx->par_col2sub[i];
62334a97f8cSStefano Zampini       PetscBool error_found = PETSC_FALSE;
62434a97f8cSStefano Zampini       ierr = PetscViewerASCIISynchronizedAllow(pcbddc->dbg_viewer,PETSC_TRUE);CHKERRQ(ierr);
62534a97f8cSStefano Zampini 
62634a97f8cSStefano Zampini       if (deluxe_ctx->par_ksp[i]) {
62734a97f8cSStefano Zampini         /* create random test vec being zero on internal nodes of the extende dirichlet problem */
62834a97f8cSStefano Zampini         ierr = VecDuplicate(deluxe_ctx->par_vec[i],&test_vec);CHKERRQ(ierr);
629b96c3477SStefano Zampini         ierr = VecSetRandom(deluxe_ctx->par_work1[i],PETSC_NULL);CHKERRQ(ierr);
63034a97f8cSStefano Zampini         ierr = VecSet(test_vec,0.0);CHKERRQ(ierr);
631b96c3477SStefano Zampini         ierr = VecScatterBegin(deluxe_ctx->par_scctx_p[i],deluxe_ctx->par_work1[i],test_vec,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
632b96c3477SStefano Zampini         ierr = VecScatterEnd(deluxe_ctx->par_scctx_p[i],deluxe_ctx->par_work1[i],test_vec,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
63334a97f8cSStefano Zampini         /* w_j */
634b96c3477SStefano Zampini         ierr = VecScatterBegin(deluxe_ctx->par_scctx_p[i],test_vec,deluxe_ctx->par_work1[i],INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
635b96c3477SStefano Zampini         ierr = VecScatterEnd(deluxe_ctx->par_scctx_p[i],test_vec,deluxe_ctx->par_work1[i],INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
63634a97f8cSStefano Zampini         /* S_j*w_j */
637b96c3477SStefano Zampini         ierr = MatMult(sub_schurs->S_Ej[subidx],deluxe_ctx->par_work1[i],deluxe_ctx->par_work2[i]);CHKERRQ(ierr);
63834a97f8cSStefano Zampini         /* \sum_j S_j*w_j */
63934a97f8cSStefano Zampini         ierr = VecSet(deluxe_ctx->par_vec[i],0.0);CHKERRQ(ierr);
640b96c3477SStefano Zampini         ierr = VecScatterBegin(deluxe_ctx->par_scctx_p[i],deluxe_ctx->par_work2[i],deluxe_ctx->par_vec[i],ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
641b96c3477SStefano Zampini         ierr = VecScatterEnd(deluxe_ctx->par_scctx_p[i],deluxe_ctx->par_work2[i],deluxe_ctx->par_vec[i],ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
64234a97f8cSStefano Zampini         /* (\sum_j S_j)^(-1)(\sum_j S_j*w_j) */
64334a97f8cSStefano Zampini         ierr = KSPSolve(deluxe_ctx->par_ksp[i],deluxe_ctx->par_vec[i],deluxe_ctx->par_vec[i]);CHKERRQ(ierr);
644b96c3477SStefano Zampini         ierr = VecScatterBegin(deluxe_ctx->par_scctx_p[i],deluxe_ctx->par_vec[i],deluxe_ctx->par_work1[i],INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
645b96c3477SStefano Zampini         ierr = VecScatterEnd(deluxe_ctx->par_scctx_p[i],deluxe_ctx->par_vec[i],deluxe_ctx->par_work1[i],INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
64634a97f8cSStefano Zampini         ierr = VecSet(deluxe_ctx->par_vec[i],0.0);CHKERRQ(ierr);
647b96c3477SStefano Zampini         ierr = VecScatterBegin(deluxe_ctx->par_scctx_p[i],deluxe_ctx->par_work1[i],deluxe_ctx->par_vec[i],INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
648b96c3477SStefano Zampini         ierr = VecScatterEnd(deluxe_ctx->par_scctx_p[i],deluxe_ctx->par_work1[i],deluxe_ctx->par_vec[i],INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
64934a97f8cSStefano Zampini         /* test partition of unity */
65034a97f8cSStefano Zampini         ierr = VecAXPY(test_vec,-1.0,deluxe_ctx->par_vec[i]);CHKERRQ(ierr);
65134a97f8cSStefano Zampini         ierr = VecNorm(test_vec,NORM_INFINITY,&error);CHKERRQ(ierr);
652c63f45b2SStefano Zampini         if (PetscAbsReal(error) > 1.e-2) {
65334a97f8cSStefano Zampini           /* ierr = VecView(test_vec,0);CHKERRQ(ierr); */
65434a97f8cSStefano Zampini           error_found = PETSC_TRUE;
65534a97f8cSStefano Zampini         }
65634a97f8cSStefano Zampini         ierr = VecDestroy(&test_vec);CHKERRQ(ierr);
65734a97f8cSStefano Zampini       }
65834a97f8cSStefano Zampini       if (error_found) {
65934a97f8cSStefano Zampini         ierr = PetscViewerASCIISynchronizedPrintf(pcbddc->dbg_viewer,"Error testing local schur for color %d and subdomain %d\n",i,PetscGlobalRank);CHKERRQ(ierr);
66034a97f8cSStefano Zampini       }
66134a97f8cSStefano Zampini       ierr = PetscViewerFlush(pcbddc->dbg_viewer);CHKERRQ(ierr);
66234a97f8cSStefano Zampini     }
66334a97f8cSStefano Zampini   }
66434a97f8cSStefano Zampini   PetscFunctionReturn(0);
66534a97f8cSStefano Zampini }
66634a97f8cSStefano Zampini 
6675db18549SStefano Zampini #undef __FUNCT__
668883469d8SStefano Zampini #define __FUNCT__ "PCBDDCScalingSetUp_Deluxe_Seq"
669883469d8SStefano Zampini static PetscErrorCode PCBDDCScalingSetUp_Deluxe_Seq(PC pc)
6705db18549SStefano Zampini {
6715db18549SStefano Zampini   PC_BDDC                *pcbddc=(PC_BDDC*)pc->data;
6725db18549SStefano Zampini   PCBDDCDeluxeScaling    deluxe_ctx=pcbddc->deluxe_ctx;
673b96c3477SStefano Zampini   PCBDDCSubSchurs        sub_schurs = pcbddc->sub_schurs;
6745db18549SStefano Zampini   PC                     pc_temp;
6755db18549SStefano Zampini   MatSolverPackage       solver=NULL;
6765db18549SStefano Zampini   char                   ksp_prefix[256];
6775db18549SStefano Zampini   size_t                 len;
6785db18549SStefano Zampini   PetscInt               local_size;
6795db18549SStefano Zampini   PetscErrorCode         ierr;
6805db18549SStefano Zampini 
6815db18549SStefano Zampini   PetscFunctionBegin;
6829221af80SStefano Zampini   if (!sub_schurs->n_subs_seq_g) {
6839221af80SStefano Zampini     PetscFunctionReturn(0);
6849221af80SStefano Zampini   }
6859221af80SStefano Zampini 
6865db18549SStefano Zampini   /* Create work vectors for sequential part of deluxe */
68765d8bf0aSStefano Zampini   ierr = MatCreateVecs(sub_schurs->sum_S_Ej_all,&deluxe_ctx->seq_work1,&deluxe_ctx->seq_work2);CHKERRQ(ierr);
6885db18549SStefano Zampini 
6895db18549SStefano Zampini   /* Compute deluxe sequential scatter */
6905db18549SStefano Zampini   ierr = VecScatterCreate(pcbddc->work_scaling,sub_schurs->is_Ej_all,deluxe_ctx->seq_work1,NULL,&deluxe_ctx->seq_scctx);CHKERRQ(ierr);
6915db18549SStefano Zampini 
6925db18549SStefano Zampini   /* Create KSP object for sequential part of deluxe scaling */
6935db18549SStefano Zampini   ierr = KSPCreate(PETSC_COMM_SELF,&deluxe_ctx->seq_ksp);CHKERRQ(ierr);
6945db18549SStefano Zampini   ierr = KSPSetOperators(deluxe_ctx->seq_ksp,sub_schurs->sum_S_Ej_all,sub_schurs->sum_S_Ej_all);CHKERRQ(ierr);
6955db18549SStefano Zampini   ierr = KSPSetType(deluxe_ctx->seq_ksp,KSPPREONLY);CHKERRQ(ierr);
6965db18549SStefano Zampini   ierr = KSPGetPC(deluxe_ctx->seq_ksp,&pc_temp);CHKERRQ(ierr);
6975db18549SStefano Zampini   ierr = PCSetType(pc_temp,PCLU);CHKERRQ(ierr);
6985db18549SStefano Zampini   ierr = KSPGetPC(pcbddc->ksp_D,&pc_temp);CHKERRQ(ierr);
6995db18549SStefano Zampini   ierr = PCFactorGetMatSolverPackage(pc_temp,(const MatSolverPackage*)&solver);CHKERRQ(ierr);
7005db18549SStefano Zampini   ierr = MatGetSize(sub_schurs->sum_S_Ej_all,&local_size,NULL);CHKERRQ(ierr);
7015db18549SStefano Zampini   if (solver && local_size) { /* if local_size is null, some external packages will report errors */
7025db18549SStefano Zampini     PC     new_pc;
7035db18549SStefano Zampini     PCType type;
7045db18549SStefano Zampini     ierr = PCGetType(pc_temp,&type);CHKERRQ(ierr);
7055db18549SStefano Zampini     ierr = KSPGetPC(deluxe_ctx->seq_ksp,&new_pc);CHKERRQ(ierr);
7065db18549SStefano Zampini     ierr = PCSetType(new_pc,type);CHKERRQ(ierr);
7075db18549SStefano Zampini     ierr = PCFactorSetMatSolverPackage(new_pc,solver);CHKERRQ(ierr);
7085db18549SStefano Zampini   }
7095db18549SStefano Zampini   ierr = PetscStrlen(((PetscObject)(pcbddc->ksp_D))->prefix,&len);CHKERRQ(ierr);
7105db18549SStefano Zampini   len -= 10; /* remove "dirichlet_" */
7115db18549SStefano Zampini   ierr = PetscStrncpy(ksp_prefix,((PetscObject)(pcbddc->ksp_D))->prefix,len+1);CHKERRQ(ierr);
7125db18549SStefano Zampini   ierr = PetscStrcat(ksp_prefix,"deluxe_seq_");CHKERRQ(ierr);
7135db18549SStefano Zampini   ierr = KSPSetOptionsPrefix(deluxe_ctx->seq_ksp,ksp_prefix);CHKERRQ(ierr);
7145db18549SStefano Zampini   if (local_size) {
7155db18549SStefano Zampini     ierr = KSPSetFromOptions(deluxe_ctx->seq_ksp);CHKERRQ(ierr);
7165db18549SStefano Zampini   }
7175db18549SStefano Zampini   ierr = KSPSetUp(deluxe_ctx->seq_ksp);CHKERRQ(ierr);
7185db18549SStefano Zampini   PetscFunctionReturn(0);
7195db18549SStefano Zampini }
720