xref: /petsc/src/ksp/pc/impls/bddc/bddcscalingbasic.c (revision 65d8bf0a82b01a6bbcb5d0fd403f7ff82932c095)
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;
36b1b3d7a2SStefano Zampini   PCBDDCSubSchurs     sub_schurs = pcbddc->sub_schurs[1];
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 */
4534a97f8cSStefano Zampini     PetscScalar *array_x,*array_D,*array;
46674ae819SStefano Zampini     ierr = VecGetArray(x,&array_x);CHKERRQ(ierr);
47674ae819SStefano Zampini     ierr = VecGetArray(pcis->D,&array_D);CHKERRQ(ierr);
48674ae819SStefano Zampini     ierr = VecGetArray(pcbddc->work_scaling,&array);CHKERRQ(ierr);
49674ae819SStefano Zampini     for (i=0;i<deluxe_ctx->n_simple;i++) {
50674ae819SStefano Zampini       array[deluxe_ctx->idx_simple_B[i]] = array_x[deluxe_ctx->idx_simple_B[i]]*array_D[deluxe_ctx->idx_simple_B[i]];
51674ae819SStefano Zampini     }
52674ae819SStefano Zampini     ierr = VecRestoreArray(pcbddc->work_scaling,&array);CHKERRQ(ierr);
53674ae819SStefano Zampini     ierr = VecRestoreArray(pcis->D,&array_D);CHKERRQ(ierr);
54674ae819SStefano Zampini     ierr = VecRestoreArray(x,&array_x);CHKERRQ(ierr);
5534a97f8cSStefano Zampini   }
5634a97f8cSStefano Zampini   /* sequential part : all problems and Schur applications collapsed into a single matrix vector multiplication and ksp solution */
5741c3ba1bSStefano Zampini   if (deluxe_ctx->seq_ksp) {
58674ae819SStefano Zampini     ierr = VecScatterBegin(deluxe_ctx->seq_scctx,x,deluxe_ctx->seq_work1,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
59674ae819SStefano Zampini     ierr = VecScatterEnd(deluxe_ctx->seq_scctx,x,deluxe_ctx->seq_work1,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
601d850880SStefano Zampini     ierr = MatMultTranspose(sub_schurs->S_Ej_all,deluxe_ctx->seq_work1,deluxe_ctx->seq_work2);CHKERRQ(ierr);
611d850880SStefano Zampini     ierr = KSPSolveTranspose(deluxe_ctx->seq_ksp,deluxe_ctx->seq_work2,deluxe_ctx->seq_work1);CHKERRQ(ierr);
62674ae819SStefano Zampini     /* fake work due to final ADD VALUES and vertices scaling needed? TODO: check it */
63674ae819SStefano Zampini     ierr = VecScatterBegin(deluxe_ctx->seq_scctx,deluxe_ctx->seq_work1,pcbddc->work_scaling,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
64674ae819SStefano Zampini     ierr = VecScatterEnd(deluxe_ctx->seq_scctx,deluxe_ctx->seq_work1,pcbddc->work_scaling,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
65674ae819SStefano Zampini   }
66674ae819SStefano Zampini   /* parallel part */
67674ae819SStefano Zampini   for (i=0;i<deluxe_ctx->par_colors;i++) {
68674ae819SStefano Zampini     if (deluxe_ctx->par_ksp[i]) {
6934a97f8cSStefano Zampini       PetscMPIInt color_rank;
7034a97f8cSStefano Zampini       PetscInt    subidx = deluxe_ctx->par_col2sub[i];
7134a97f8cSStefano Zampini       /* restrict on subset */
7234a97f8cSStefano Zampini       ierr = VecScatterBegin(deluxe_ctx->par_scctx_s[i],x,sub_schurs->work1[subidx],INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
7334a97f8cSStefano Zampini       ierr = VecScatterEnd(deluxe_ctx->par_scctx_s[i],x,sub_schurs->work1[subidx],INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
7434a97f8cSStefano Zampini       /* S_Ej */
7534a97f8cSStefano Zampini       ierr = MatMult(sub_schurs->S_Ej[subidx],sub_schurs->work1[subidx],sub_schurs->work2[subidx]);CHKERRQ(ierr);
7634a97f8cSStefano Zampini       /* (\sum_j S_Ej)^-1 */
7734a97f8cSStefano Zampini       ierr = VecSet(deluxe_ctx->par_vec[i],0.0);CHKERRQ(ierr);
7834a97f8cSStefano Zampini       ierr = VecScatterBegin(deluxe_ctx->par_scctx_p[i],sub_schurs->work2[subidx],deluxe_ctx->par_vec[i],ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
7934a97f8cSStefano Zampini       ierr = VecScatterEnd(deluxe_ctx->par_scctx_p[i],sub_schurs->work2[subidx],deluxe_ctx->par_vec[i],ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
80674ae819SStefano Zampini       ierr = KSPSolve(deluxe_ctx->par_ksp[i],deluxe_ctx->par_vec[i],deluxe_ctx->par_vec[i]);CHKERRQ(ierr);
8134a97f8cSStefano Zampini       ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)deluxe_ctx->par_ksp[i]),&color_rank);CHKERRQ(ierr);
8234a97f8cSStefano Zampini       /* get back solution on subset */
8334a97f8cSStefano Zampini       ierr = VecScatterBegin(deluxe_ctx->par_scctx_p[i],deluxe_ctx->par_vec[i],sub_schurs->work1[subidx],INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
8434a97f8cSStefano Zampini       ierr = VecScatterEnd(deluxe_ctx->par_scctx_p[i],deluxe_ctx->par_vec[i],sub_schurs->work1[subidx],INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
8534a97f8cSStefano Zampini       if (!color_rank) { /* only the master process in coloured comm copies the computed values */
8634a97f8cSStefano Zampini         ierr = VecScatterBegin(deluxe_ctx->par_scctx_s[i],sub_schurs->work1[subidx],pcbddc->work_scaling,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
8734a97f8cSStefano Zampini         ierr = VecScatterEnd(deluxe_ctx->par_scctx_s[i],sub_schurs->work1[subidx],pcbddc->work_scaling,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
88674ae819SStefano Zampini       }
89674ae819SStefano Zampini     }
90674ae819SStefano Zampini   }
91674ae819SStefano Zampini   /* put local boundary part in global vector */
9234a97f8cSStefano Zampini   ierr = VecSet(y,0.0);CHKERRQ(ierr);
93674ae819SStefano Zampini   ierr = VecScatterBegin(pcis->global_to_B,pcbddc->work_scaling,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
94674ae819SStefano Zampini   ierr = VecScatterEnd(pcis->global_to_B,pcbddc->work_scaling,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
95674ae819SStefano Zampini   PetscFunctionReturn(0);
96674ae819SStefano Zampini }
97674ae819SStefano Zampini 
98674ae819SStefano Zampini #undef __FUNCT__
99674ae819SStefano Zampini #define __FUNCT__ "PCBDDCScalingExtension"
100674ae819SStefano Zampini PetscErrorCode PCBDDCScalingExtension(PC pc, Vec local_interface_vector, Vec global_vector)
101674ae819SStefano Zampini {
102674ae819SStefano Zampini   PC_BDDC *pcbddc=(PC_BDDC*)pc->data;
103674ae819SStefano Zampini   PetscErrorCode ierr;
104674ae819SStefano Zampini 
105674ae819SStefano Zampini   PetscFunctionBegin;
106674ae819SStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
107674ae819SStefano Zampini   PetscValidHeaderSpecific(local_interface_vector,VEC_CLASSID,2);
108674ae819SStefano Zampini   PetscValidHeaderSpecific(global_vector,VEC_CLASSID,3);
109674ae819SStefano Zampini   if (local_interface_vector == pcbddc->work_scaling) {
110674ae819SStefano Zampini     SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"Local vector cannot be pcbddc->work_scaling!\n");
111674ae819SStefano Zampini   }
112674ae819SStefano Zampini   ierr = PetscTryMethod(pc,"PCBDDCScalingExtension_C",(PC,Vec,Vec),(pc,local_interface_vector,global_vector));CHKERRQ(ierr);
113674ae819SStefano Zampini   PetscFunctionReturn(0);
114674ae819SStefano Zampini }
115674ae819SStefano Zampini 
116674ae819SStefano Zampini #undef __FUNCT__
117674ae819SStefano Zampini #define __FUNCT__ "PCBDDCScalingRestriction_Basic"
118674ae819SStefano Zampini static PetscErrorCode PCBDDCScalingRestriction_Basic(PC pc, Vec global_vector, Vec local_interface_vector)
119674ae819SStefano Zampini {
120674ae819SStefano Zampini   PetscErrorCode ierr;
121674ae819SStefano Zampini   PC_IS* pcis = (PC_IS*)pc->data;
122674ae819SStefano Zampini 
123674ae819SStefano Zampini   PetscFunctionBegin;
124674ae819SStefano Zampini   ierr = VecScatterBegin(pcis->global_to_B,global_vector,local_interface_vector,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
125674ae819SStefano Zampini   ierr = VecScatterEnd(pcis->global_to_B,global_vector,local_interface_vector,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
126674ae819SStefano Zampini   /* Apply partition of unity */
127674ae819SStefano Zampini   ierr = VecPointwiseMult(local_interface_vector,pcis->D,local_interface_vector);CHKERRQ(ierr);
128674ae819SStefano Zampini   PetscFunctionReturn(0);
129674ae819SStefano Zampini }
130674ae819SStefano Zampini 
131674ae819SStefano Zampini #undef __FUNCT__
132674ae819SStefano Zampini #define __FUNCT__ "PCBDDCScalingRestriction_Deluxe"
133674ae819SStefano Zampini static PetscErrorCode PCBDDCScalingRestriction_Deluxe(PC pc, Vec x, Vec y)
134674ae819SStefano Zampini {
135674ae819SStefano Zampini   PC_IS*              pcis=(PC_IS*)pc->data;
136674ae819SStefano Zampini   PC_BDDC*            pcbddc=(PC_BDDC*)pc->data;
137674ae819SStefano Zampini   PCBDDCDeluxeScaling deluxe_ctx = pcbddc->deluxe_ctx;
138b1b3d7a2SStefano Zampini   PCBDDCSubSchurs     sub_schurs = pcbddc->sub_schurs[1];
139674ae819SStefano Zampini   PetscInt            i;
140674ae819SStefano Zampini   PetscErrorCode      ierr;
141674ae819SStefano Zampini 
142674ae819SStefano Zampini   PetscFunctionBegin;
143674ae819SStefano Zampini   /* get local boundary part of global vector */
144674ae819SStefano Zampini   ierr = VecScatterBegin(pcis->global_to_B,x,y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
145674ae819SStefano Zampini   ierr = VecScatterEnd(pcis->global_to_B,x,y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
14634a97f8cSStefano Zampini   if (deluxe_ctx->n_simple) {
14734a97f8cSStefano Zampini     /* scale deluxe vertices using diagonal scaling */
14834a97f8cSStefano Zampini     PetscScalar *array_y,*array_D;
149674ae819SStefano Zampini     ierr = VecGetArray(y,&array_y);CHKERRQ(ierr);
150674ae819SStefano Zampini     ierr = VecGetArray(pcis->D,&array_D);CHKERRQ(ierr);
151674ae819SStefano Zampini     for (i=0;i<deluxe_ctx->n_simple;i++) {
152674ae819SStefano Zampini       array_y[deluxe_ctx->idx_simple_B[i]] *= array_D[deluxe_ctx->idx_simple_B[i]];
153674ae819SStefano Zampini     }
154674ae819SStefano Zampini     ierr = VecRestoreArray(pcis->D,&array_D);CHKERRQ(ierr);
155674ae819SStefano Zampini     ierr = VecRestoreArray(y,&array_y);CHKERRQ(ierr);
15634a97f8cSStefano Zampini   }
15734a97f8cSStefano Zampini   /* sequential part : all problems and Schur applications collapsed into a single matrix vector multiplication and ksp solution */
15841c3ba1bSStefano Zampini   if (deluxe_ctx->seq_ksp) {
159674ae819SStefano Zampini     ierr = VecScatterBegin(deluxe_ctx->seq_scctx,y,deluxe_ctx->seq_work1,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
160674ae819SStefano Zampini     ierr = VecScatterEnd(deluxe_ctx->seq_scctx,y,deluxe_ctx->seq_work1,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1611d850880SStefano Zampini     ierr = KSPSolve(deluxe_ctx->seq_ksp,deluxe_ctx->seq_work1,deluxe_ctx->seq_work2);CHKERRQ(ierr);
1621d850880SStefano Zampini     ierr = MatMult(sub_schurs->S_Ej_all,deluxe_ctx->seq_work2,deluxe_ctx->seq_work1);CHKERRQ(ierr);
163674ae819SStefano Zampini     ierr = VecScatterBegin(deluxe_ctx->seq_scctx,deluxe_ctx->seq_work1,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
164674ae819SStefano Zampini     ierr = VecScatterEnd(deluxe_ctx->seq_scctx,deluxe_ctx->seq_work1,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
165674ae819SStefano Zampini   }
166674ae819SStefano Zampini   /* parallel part */
167674ae819SStefano Zampini   for (i=0;i<deluxe_ctx->par_colors;i++) {
168674ae819SStefano Zampini     if (deluxe_ctx->par_ksp[i]) {
16934a97f8cSStefano Zampini       PetscInt subidx = deluxe_ctx->par_col2sub[i];
17034a97f8cSStefano Zampini       /* restrict on subset */
17134a97f8cSStefano Zampini       ierr = VecScatterBegin(deluxe_ctx->par_scctx_s[i],y,sub_schurs->work1[subidx],INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
17234a97f8cSStefano Zampini       ierr = VecScatterEnd(deluxe_ctx->par_scctx_s[i],y,sub_schurs->work1[subidx],INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
17334a97f8cSStefano Zampini       /* (\sum_j S_Ej)^-T */
17434a97f8cSStefano Zampini       ierr = VecSet(deluxe_ctx->par_vec[i],0.0);CHKERRQ(ierr);
17534a97f8cSStefano Zampini       ierr = VecScatterBegin(deluxe_ctx->par_scctx_p[i],sub_schurs->work1[subidx],deluxe_ctx->par_vec[i],INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
17634a97f8cSStefano Zampini       ierr = VecScatterEnd(deluxe_ctx->par_scctx_p[i],sub_schurs->work1[subidx],deluxe_ctx->par_vec[i],INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
177674ae819SStefano Zampini       ierr = KSPSolveTranspose(deluxe_ctx->par_ksp[i],deluxe_ctx->par_vec[i],deluxe_ctx->par_vec[i]);CHKERRQ(ierr);
17834a97f8cSStefano Zampini       ierr = VecScatterBegin(deluxe_ctx->par_scctx_p[i],deluxe_ctx->par_vec[i],sub_schurs->work1[subidx],INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
17934a97f8cSStefano Zampini       ierr = VecScatterEnd(deluxe_ctx->par_scctx_p[i],deluxe_ctx->par_vec[i],sub_schurs->work1[subidx],INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
18034a97f8cSStefano Zampini       /* S_Ej^T */
18134a97f8cSStefano Zampini       ierr = MatMultTranspose(sub_schurs->S_Ej[subidx],sub_schurs->work1[subidx],sub_schurs->work2[subidx]);CHKERRQ(ierr);
18234a97f8cSStefano Zampini       /* extend to boundary */
18334a97f8cSStefano Zampini       ierr = VecScatterBegin(deluxe_ctx->par_scctx_s[i],sub_schurs->work2[subidx],y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
18434a97f8cSStefano Zampini       ierr = VecScatterEnd(deluxe_ctx->par_scctx_s[i],sub_schurs->work2[subidx],y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
185674ae819SStefano Zampini     }
186674ae819SStefano Zampini   }
187674ae819SStefano Zampini   PetscFunctionReturn(0);
188674ae819SStefano Zampini }
189674ae819SStefano Zampini 
190674ae819SStefano Zampini #undef __FUNCT__
191674ae819SStefano Zampini #define __FUNCT__ "PCBDDCScalingRestriction"
192674ae819SStefano Zampini PetscErrorCode PCBDDCScalingRestriction(PC pc, Vec global_vector, Vec local_interface_vector)
193674ae819SStefano Zampini {
194674ae819SStefano Zampini   PC_BDDC        *pcbddc=(PC_BDDC*)pc->data;
195674ae819SStefano Zampini   PetscErrorCode ierr;
196674ae819SStefano Zampini 
197674ae819SStefano Zampini   PetscFunctionBegin;
198674ae819SStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
199674ae819SStefano Zampini   PetscValidHeaderSpecific(global_vector,VEC_CLASSID,2);
200674ae819SStefano Zampini   PetscValidHeaderSpecific(local_interface_vector,VEC_CLASSID,3);
201674ae819SStefano Zampini   if (local_interface_vector == pcbddc->work_scaling) {
202674ae819SStefano Zampini     SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"Local vector should cannot be pcbddc->work_scaling!\n");
203674ae819SStefano Zampini   }
204674ae819SStefano Zampini   ierr = PetscTryMethod(pc,"PCBDDCScalingRestriction_C",(PC,Vec,Vec),(pc,global_vector,local_interface_vector));CHKERRQ(ierr);
205674ae819SStefano Zampini   PetscFunctionReturn(0);
206674ae819SStefano Zampini }
207674ae819SStefano Zampini 
208674ae819SStefano Zampini #undef __FUNCT__
209674ae819SStefano Zampini #define __FUNCT__ "PCBDDCScalingSetUp"
210674ae819SStefano Zampini PetscErrorCode PCBDDCScalingSetUp(PC pc)
211674ae819SStefano Zampini {
212674ae819SStefano Zampini   PC_IS* pcis=(PC_IS*)pc->data;
213674ae819SStefano Zampini   PC_BDDC* pcbddc=(PC_BDDC*)pc->data;
214674ae819SStefano Zampini   PetscErrorCode ierr;
215674ae819SStefano Zampini 
216674ae819SStefano Zampini   PetscFunctionBegin;
217674ae819SStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
218674ae819SStefano Zampini   /* create work vector for the operator */
21934a97f8cSStefano Zampini   ierr = VecDestroy(&pcbddc->work_scaling);CHKERRQ(ierr);
220674ae819SStefano Zampini   ierr = VecDuplicate(pcis->vec1_B,&pcbddc->work_scaling);CHKERRQ(ierr);
22134a97f8cSStefano Zampini   /* always rebuild pcis->D */
22228d874f6SStefano Zampini   if (pcis->use_stiffness_scaling) {
223674ae819SStefano Zampini     ierr = MatGetDiagonal(pcbddc->local_mat,pcis->vec1_N);CHKERRQ(ierr);
224674ae819SStefano Zampini     ierr = VecScatterBegin(pcis->N_to_B,pcis->vec1_N,pcis->D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
225674ae819SStefano Zampini     ierr = VecScatterEnd(pcis->N_to_B,pcis->vec1_N,pcis->D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
226674ae819SStefano Zampini   }
227674ae819SStefano Zampini   ierr = VecCopy(pcis->D,pcis->vec1_B);CHKERRQ(ierr);
228674ae819SStefano Zampini   ierr = VecSet(pcis->vec1_global,0.0);CHKERRQ(ierr);
229674ae819SStefano Zampini   ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_B,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
230674ae819SStefano Zampini   ierr = VecScatterEnd(pcis->global_to_B,pcis->vec1_B,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
231674ae819SStefano Zampini   ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_global,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
232674ae819SStefano Zampini   ierr = VecScatterEnd(pcis->global_to_B,pcis->vec1_global,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
233674ae819SStefano Zampini   ierr = VecPointwiseDivide(pcis->D,pcis->D,pcis->vec1_B);CHKERRQ(ierr);
234674ae819SStefano Zampini   /* now setup */
235681e7c04SStefano Zampini   if (pcbddc->use_deluxe_scaling) {
23634a97f8cSStefano Zampini     if (!pcbddc->deluxe_ctx) {
23734a97f8cSStefano Zampini       ierr = PCBDDCScalingCreate_Deluxe(pc);CHKERRQ(ierr);
23834a97f8cSStefano Zampini     }
23934a97f8cSStefano Zampini     ierr = PCBDDCScalingSetUp_Deluxe(pc);CHKERRQ(ierr);
240674ae819SStefano Zampini     ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCScalingRestriction_C",PCBDDCScalingRestriction_Deluxe);CHKERRQ(ierr);
241674ae819SStefano Zampini     ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCScalingExtension_C",PCBDDCScalingExtension_Deluxe);CHKERRQ(ierr);
242674ae819SStefano Zampini   } else {
243674ae819SStefano Zampini     ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCScalingRestriction_C",PCBDDCScalingRestriction_Basic);CHKERRQ(ierr);
244674ae819SStefano Zampini     ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCScalingExtension_C",PCBDDCScalingExtension_Basic);CHKERRQ(ierr);
245674ae819SStefano Zampini   }
24634a97f8cSStefano Zampini 
247674ae819SStefano Zampini   /* test */
248674ae819SStefano Zampini   if (pcbddc->dbg_flag) {
24934a97f8cSStefano Zampini     Vec         vec2_global;
250674ae819SStefano Zampini     PetscViewer viewer=pcbddc->dbg_viewer;
25134a97f8cSStefano Zampini     PetscReal   error;
252674ae819SStefano Zampini 
253674ae819SStefano Zampini     /* extension -> from local to parallel */
25434a97f8cSStefano Zampini     ierr = VecSet(pcis->vec1_global,0.0);CHKERRQ(ierr);
25534a97f8cSStefano Zampini     ierr = VecSetRandom(pcis->vec1_B,NULL);CHKERRQ(ierr);
25634a97f8cSStefano Zampini     ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_B,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
25734a97f8cSStefano Zampini     ierr = VecScatterEnd(pcis->global_to_B,pcis->vec1_B,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
25834a97f8cSStefano Zampini     ierr = VecDuplicate(pcis->vec1_global,&vec2_global);CHKERRQ(ierr);
25934a97f8cSStefano Zampini     ierr = VecCopy(pcis->vec1_global,vec2_global);CHKERRQ(ierr);
26034a97f8cSStefano Zampini 
261674ae819SStefano Zampini     ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_global,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
262674ae819SStefano Zampini     ierr = VecScatterEnd(pcis->global_to_B,pcis->vec1_global,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
263674ae819SStefano Zampini     ierr = PCBDDCScalingExtension(pc,pcis->vec1_B,pcis->vec1_global);CHKERRQ(ierr);
26434a97f8cSStefano Zampini     ierr = VecAXPY(pcis->vec1_global,-1.0,vec2_global);CHKERRQ(ierr);
26534a97f8cSStefano Zampini     ierr = VecNorm(pcis->vec1_global,NORM_INFINITY,&error);CHKERRQ(ierr);
266674ae819SStefano Zampini     ierr = PetscViewerASCIIPrintf(viewer,"Error scaling extension %1.14e\n",error);CHKERRQ(ierr);
26734a97f8cSStefano Zampini     if (error>1.e-8 && pcbddc->dbg_flag>1) {
268674ae819SStefano Zampini       ierr = VecView(pcis->vec1_global,viewer);CHKERRQ(ierr);
269674ae819SStefano Zampini     }
27034a97f8cSStefano Zampini     ierr = VecDestroy(&vec2_global);CHKERRQ(ierr);
27134a97f8cSStefano Zampini 
272674ae819SStefano Zampini     /* restriction -> from parallel to local */
273674ae819SStefano Zampini     ierr = VecSet(pcis->vec1_global,0.0);CHKERRQ(ierr);
27434a97f8cSStefano Zampini     ierr = VecSetRandom(pcis->vec1_B,NULL);CHKERRQ(ierr);
275674ae819SStefano Zampini     ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_B,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
276674ae819SStefano Zampini     ierr = VecScatterEnd(pcis->global_to_B,pcis->vec1_B,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
27734a97f8cSStefano Zampini 
27834a97f8cSStefano Zampini     ierr = PCBDDCScalingRestriction(pc,pcis->vec1_global,pcis->vec1_B);CHKERRQ(ierr);
27934a97f8cSStefano Zampini     ierr = VecScale(pcis->vec1_B,-1.0);CHKERRQ(ierr);
28034a97f8cSStefano Zampini     ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_B,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
28134a97f8cSStefano Zampini     ierr = VecScatterEnd(pcis->global_to_B,pcis->vec1_B,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
28234a97f8cSStefano Zampini     ierr = VecNorm(pcis->vec1_global,NORM_INFINITY,&error);CHKERRQ(ierr);
28334a97f8cSStefano Zampini     ierr = PetscViewerASCIIPrintf(viewer,"Error scaling restriction %1.14e\n",error);CHKERRQ(ierr);
28434a97f8cSStefano Zampini     if (error>1.e-8 && pcbddc->dbg_flag>1) {
285674ae819SStefano Zampini       ierr = VecView(pcis->vec1_global,viewer);CHKERRQ(ierr);
286674ae819SStefano Zampini     }
287674ae819SStefano Zampini   }
288674ae819SStefano Zampini   PetscFunctionReturn(0);
289674ae819SStefano Zampini }
290674ae819SStefano Zampini 
291674ae819SStefano Zampini #undef __FUNCT__
292674ae819SStefano Zampini #define __FUNCT__ "PCBDDCScalingDestroy"
293674ae819SStefano Zampini PetscErrorCode PCBDDCScalingDestroy(PC pc)
294674ae819SStefano Zampini {
295674ae819SStefano Zampini   PC_BDDC* pcbddc=(PC_BDDC*)pc->data;
296674ae819SStefano Zampini   PetscErrorCode ierr;
297674ae819SStefano Zampini 
298674ae819SStefano Zampini   PetscFunctionBegin;
29934a97f8cSStefano Zampini   if (pcbddc->deluxe_ctx) {
30034a97f8cSStefano Zampini     ierr = PCBDDCScalingDestroy_Deluxe(pc);CHKERRQ(ierr);
301674ae819SStefano Zampini   }
302674ae819SStefano Zampini   ierr = VecDestroy(&pcbddc->work_scaling);CHKERRQ(ierr);
303674ae819SStefano Zampini   /* remove functions */
304674ae819SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCScalingRestriction_C",NULL);CHKERRQ(ierr);
305674ae819SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCScalingExtension_C",NULL);CHKERRQ(ierr);
306674ae819SStefano Zampini   PetscFunctionReturn(0);
307674ae819SStefano Zampini }
308674ae819SStefano Zampini 
30934a97f8cSStefano Zampini #undef __FUNCT__
31034a97f8cSStefano Zampini #define __FUNCT__ "PCBDDCScalingCreate_Deluxe"
31134a97f8cSStefano Zampini static PetscErrorCode PCBDDCScalingCreate_Deluxe(PC pc)
31234a97f8cSStefano Zampini {
31334a97f8cSStefano Zampini   PC_BDDC*            pcbddc=(PC_BDDC*)pc->data;
31434a97f8cSStefano Zampini   PCBDDCDeluxeScaling deluxe_ctx;
31534a97f8cSStefano Zampini   PetscErrorCode      ierr;
31634a97f8cSStefano Zampini 
31734a97f8cSStefano Zampini   PetscFunctionBegin;
31834a97f8cSStefano Zampini   ierr = PetscNew(&deluxe_ctx);CHKERRQ(ierr);
31934a97f8cSStefano Zampini   pcbddc->deluxe_ctx = deluxe_ctx;
32034a97f8cSStefano Zampini   PetscFunctionReturn(0);
32134a97f8cSStefano Zampini }
32234a97f8cSStefano Zampini 
32334a97f8cSStefano Zampini #undef __FUNCT__
32434a97f8cSStefano Zampini #define __FUNCT__ "PCBDDCScalingDestroy_Deluxe"
32534a97f8cSStefano Zampini static PetscErrorCode PCBDDCScalingDestroy_Deluxe(PC pc)
32634a97f8cSStefano Zampini {
32734a97f8cSStefano Zampini   PC_BDDC*            pcbddc=(PC_BDDC*)pc->data;
32834a97f8cSStefano Zampini   PetscErrorCode      ierr;
32934a97f8cSStefano Zampini 
33034a97f8cSStefano Zampini   PetscFunctionBegin;
33134a97f8cSStefano Zampini   ierr = PCBDDCScalingReset_Deluxe_Solvers(pcbddc->deluxe_ctx);CHKERRQ(ierr);
33234a97f8cSStefano Zampini   ierr = PetscFree(pcbddc->deluxe_ctx);CHKERRQ(ierr);
33334a97f8cSStefano Zampini   PetscFunctionReturn(0);
33434a97f8cSStefano Zampini }
33534a97f8cSStefano Zampini 
33634a97f8cSStefano Zampini #undef __FUNCT__
33734a97f8cSStefano Zampini #define __FUNCT__ "PCBDDCScalingReset_Deluxe_Solvers"
33834a97f8cSStefano Zampini static PetscErrorCode PCBDDCScalingReset_Deluxe_Solvers(PCBDDCDeluxeScaling deluxe_ctx)
33934a97f8cSStefano Zampini {
34034a97f8cSStefano Zampini   PetscErrorCode      ierr;
34134a97f8cSStefano Zampini 
34234a97f8cSStefano Zampini   PetscFunctionBegin;
34334a97f8cSStefano Zampini   ierr = PetscFree(deluxe_ctx->idx_simple_B);CHKERRQ(ierr);
34434a97f8cSStefano Zampini   deluxe_ctx->n_simple = 0;
34541c3ba1bSStefano Zampini   if (deluxe_ctx->seq_ksp) {
34634a97f8cSStefano Zampini     ierr = VecScatterDestroy(&deluxe_ctx->seq_scctx);CHKERRQ(ierr);
34734a97f8cSStefano Zampini     ierr = VecDestroy(&deluxe_ctx->seq_work1);CHKERRQ(ierr);
34834a97f8cSStefano Zampini     ierr = VecDestroy(&deluxe_ctx->seq_work2);CHKERRQ(ierr);
34934a97f8cSStefano Zampini     ierr = KSPDestroy(&deluxe_ctx->seq_ksp);CHKERRQ(ierr);
35034a97f8cSStefano Zampini   }
35134a97f8cSStefano Zampini   if (deluxe_ctx->par_colors) {
35234a97f8cSStefano Zampini     PetscInt i;
35334a97f8cSStefano Zampini     for (i=0;i<deluxe_ctx->par_colors;i++) {
35434a97f8cSStefano Zampini       ierr = VecScatterDestroy(&deluxe_ctx->par_scctx_s[i]);CHKERRQ(ierr);
35534a97f8cSStefano Zampini       ierr = VecScatterDestroy(&deluxe_ctx->par_scctx_p[i]);CHKERRQ(ierr);
35634a97f8cSStefano Zampini       ierr = VecDestroy(&deluxe_ctx->par_vec[i]);CHKERRQ(ierr);
35734a97f8cSStefano Zampini       ierr = KSPDestroy(&deluxe_ctx->par_ksp[i]);CHKERRQ(ierr);
35834a97f8cSStefano Zampini     }
35934a97f8cSStefano Zampini     ierr = PetscFree5(deluxe_ctx->par_ksp,
36034a97f8cSStefano Zampini                       deluxe_ctx->par_scctx_s,
36134a97f8cSStefano Zampini                       deluxe_ctx->par_scctx_p,
36234a97f8cSStefano Zampini                       deluxe_ctx->par_vec,
36334a97f8cSStefano Zampini                       deluxe_ctx->par_col2sub);CHKERRQ(ierr);
36434a97f8cSStefano Zampini   }
36534a97f8cSStefano Zampini   deluxe_ctx->par_colors = 0;
36634a97f8cSStefano Zampini   PetscFunctionReturn(0);
36734a97f8cSStefano Zampini }
36834a97f8cSStefano Zampini 
36934a97f8cSStefano Zampini #undef __FUNCT__
37034a97f8cSStefano Zampini #define __FUNCT__ "PCBDDCScalingSetUp_Deluxe"
37134a97f8cSStefano Zampini static PetscErrorCode PCBDDCScalingSetUp_Deluxe(PC pc)
37234a97f8cSStefano Zampini {
37334a97f8cSStefano Zampini   PC_IS               *pcis=(PC_IS*)pc->data;
37434a97f8cSStefano Zampini   PC_BDDC             *pcbddc=(PC_BDDC*)pc->data;
37534a97f8cSStefano Zampini   PCBDDCDeluxeScaling deluxe_ctx=pcbddc->deluxe_ctx;
376b1b3d7a2SStefano Zampini   PCBDDCSubSchurs     sub_schurs=pcbddc->sub_schurs[1];
37734a97f8cSStefano Zampini   PCBDDCGraph         graph;
378b1b3d7a2SStefano Zampini   Mat                 S_j;
379883469d8SStefano Zampini   PetscBT             bitmask;
380883469d8SStefano Zampini   PetscInt            i,*used_xadj,*used_adjncy;
381883469d8SStefano Zampini   const PetscInt      *idxs;
382b1b3d7a2SStefano Zampini   PetscBool           free_used_adj;
38334a97f8cSStefano Zampini   PetscErrorCode      ierr;
38434a97f8cSStefano Zampini 
38534a97f8cSStefano Zampini   PetscFunctionBegin;
38634a97f8cSStefano Zampini   /* throw away the solvers */
38734a97f8cSStefano Zampini   ierr = PCBDDCScalingReset_Deluxe_Solvers(deluxe_ctx);CHKERRQ(ierr);
38834a97f8cSStefano Zampini 
38934a97f8cSStefano Zampini   /* attach interface graph for determining subsets */
39034a97f8cSStefano Zampini   if (pcbddc->deluxe_rebuild) { /* in case rebuild has been requested, it uses a graph generated only by the neighbouring information */
39134a97f8cSStefano Zampini     PetscInt *idx_V_N;
39234a97f8cSStefano Zampini     IS       verticesIS;
39334a97f8cSStefano Zampini     ierr = PCBDDCGetPrimalVerticesLocalIdx(pc,&i,&idx_V_N);CHKERRQ(ierr);
39434a97f8cSStefano Zampini     ierr = ISCreateGeneral(PETSC_COMM_SELF,i,idx_V_N,PETSC_OWN_POINTER,&verticesIS);CHKERRQ(ierr);
39534a97f8cSStefano Zampini     ierr = PCBDDCGraphCreate(&graph);CHKERRQ(ierr);
39634a97f8cSStefano Zampini     ierr = PCBDDCGraphInit(graph,pcbddc->mat_graph->l2gmap);CHKERRQ(ierr);
3970a95f60dSStefano Zampini     ierr = PCBDDCGraphSetUp(graph,0,NULL,pcbddc->DirichletBoundariesLocal,0,NULL,verticesIS);CHKERRQ(ierr);
39834a97f8cSStefano Zampini     ierr = PCBDDCGraphComputeConnectedComponents(graph);CHKERRQ(ierr);
39934a97f8cSStefano Zampini     ierr = ISDestroy(&verticesIS);CHKERRQ(ierr);
40034a97f8cSStefano Zampini /*
40134a97f8cSStefano Zampini     if (pcbddc->dbg_flag) {
40234a97f8cSStefano Zampini       ierr = PCBDDCGraphASCIIView(graph,pcbddc->dbg_flag,pcbddc->dbg_viewer);CHKERRQ(ierr);
40334a97f8cSStefano Zampini     }
40434a97f8cSStefano Zampini */
40534a97f8cSStefano Zampini   } else {
40634a97f8cSStefano Zampini     graph = pcbddc->mat_graph;
40734a97f8cSStefano Zampini   }
40834a97f8cSStefano Zampini 
40934a97f8cSStefano Zampini   /* decide the adjacency to be used for determining internal problems for local schur on subsets */
41034a97f8cSStefano Zampini   free_used_adj = PETSC_FALSE;
41134a97f8cSStefano Zampini   if (pcbddc->deluxe_layers == -1) {
41234a97f8cSStefano Zampini     used_xadj = NULL;
41334a97f8cSStefano Zampini     used_adjncy = NULL;
41434a97f8cSStefano Zampini   } else {
41534a97f8cSStefano Zampini     if ((pcbddc->deluxe_use_useradj && pcbddc->mat_graph->xadj) || !pcbddc->deluxe_compute_rowadj) {
41634a97f8cSStefano Zampini       used_xadj = pcbddc->mat_graph->xadj;
41734a97f8cSStefano Zampini       used_adjncy = pcbddc->mat_graph->adjncy;
41834a97f8cSStefano Zampini     } else {
41934a97f8cSStefano Zampini       Mat            mat_adj;
42034a97f8cSStefano Zampini       PetscBool      flg_row=PETSC_TRUE;
42134a97f8cSStefano Zampini       const PetscInt *xadj,*adjncy;
42234a97f8cSStefano Zampini       PetscInt       nvtxs;
42334a97f8cSStefano Zampini 
42434a97f8cSStefano Zampini       ierr = MatConvert(pcbddc->local_mat,MATMPIADJ,MAT_INITIAL_MATRIX,&mat_adj);CHKERRQ(ierr);
42534a97f8cSStefano Zampini       ierr = MatGetRowIJ(mat_adj,0,PETSC_TRUE,PETSC_FALSE,&nvtxs,&xadj,&adjncy,&flg_row);CHKERRQ(ierr);
42634a97f8cSStefano Zampini       if (!flg_row) {
42734a97f8cSStefano Zampini         SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Error in MatGetRowIJ called in %s\n",__FUNCT__);
42834a97f8cSStefano Zampini       }
42934a97f8cSStefano Zampini       ierr = PetscMalloc2(nvtxs+1,&used_xadj,xadj[nvtxs],&used_adjncy);CHKERRQ(ierr);
43034a97f8cSStefano Zampini       ierr = PetscMemcpy(used_xadj,xadj,(nvtxs+1)*sizeof(*xadj));CHKERRQ(ierr);
43134a97f8cSStefano Zampini       ierr = PetscMemcpy(used_adjncy,adjncy,(xadj[nvtxs])*sizeof(*adjncy));CHKERRQ(ierr);
43234a97f8cSStefano Zampini       ierr = MatRestoreRowIJ(mat_adj,0,PETSC_TRUE,PETSC_FALSE,&nvtxs,&xadj,&adjncy,&flg_row);CHKERRQ(ierr);
43334a97f8cSStefano Zampini       if (!flg_row) {
43434a97f8cSStefano Zampini         SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Error in MatRestoreRowIJ called in %s\n",__FUNCT__);
43534a97f8cSStefano Zampini       }
43634a97f8cSStefano Zampini       ierr = MatDestroy(&mat_adj);CHKERRQ(ierr);
43734a97f8cSStefano Zampini       free_used_adj = PETSC_TRUE;
43834a97f8cSStefano Zampini     }
43934a97f8cSStefano Zampini   }
44034a97f8cSStefano Zampini 
44134a97f8cSStefano Zampini   /* Create Schur complement matrix */
44234a97f8cSStefano Zampini   ierr = MatCreateSchurComplement(pcis->A_II,pcis->A_II,pcis->A_IB,pcis->A_BI,pcis->A_BB,&S_j);CHKERRQ(ierr);
44334a97f8cSStefano Zampini   ierr = MatSchurComplementSetKSP(S_j,pcbddc->ksp_D);CHKERRQ(ierr);
4447821e9e7SStefano Zampini 
4451e9c79c2SStefano Zampini   /* sub_schurs init */ /* TODO reuse adaptive one if valid (i.e. pcbddc->local_mat == matis->A and same graph info (HOW?) ) */
4465db18549SStefano Zampini   ierr = PCBDDCSubSchursInit(sub_schurs,pcbddc->local_mat,S_j,pcis->is_I_local,pcis->is_B_local,graph,pcis->BtoNmap,pcbddc->deluxe_threshold);CHKERRQ(ierr);
447b1b3d7a2SStefano Zampini   ierr = MatDestroy(&S_j);CHKERRQ(ierr);
4481580ed26SStefano Zampini   ierr = PCBDDCSubSchursSetUp(sub_schurs,used_xadj,used_adjncy,pcbddc->deluxe_layers);CHKERRQ(ierr);
449b1b3d7a2SStefano Zampini 
450b1b3d7a2SStefano Zampini   /* Compute data structures to solve parallel problems */
451b1b3d7a2SStefano Zampini   ierr = PCBDDCScalingSetUp_Deluxe_Par(pc,sub_schurs->n_subs_par,sub_schurs->n_subs_par_g,
452b1b3d7a2SStefano Zampini                                        sub_schurs->auxglobal_parallel,
453b1b3d7a2SStefano Zampini                                        sub_schurs->index_parallel);CHKERRQ(ierr);
454b1b3d7a2SStefano Zampini   /* Compute data structures to solve sequential problems */
455883469d8SStefano Zampini   ierr = PCBDDCScalingSetUp_Deluxe_Seq(pc);CHKERRQ(ierr);
4565db18549SStefano Zampini 
457b1b3d7a2SStefano Zampini   /* free adjacency */
458b1b3d7a2SStefano Zampini   if (free_used_adj) {
459b1b3d7a2SStefano Zampini     ierr = PetscFree2(used_xadj,used_adjncy);CHKERRQ(ierr);
460b1b3d7a2SStefano Zampini   }
461b1b3d7a2SStefano Zampini 
462b1b3d7a2SStefano Zampini   /* diagonal scaling on interface dofs not contained in cc */
463b1b3d7a2SStefano Zampini   ierr = PetscBTCreate(pcis->n,&bitmask);CHKERRQ(ierr);
464b1b3d7a2SStefano Zampini   ierr = ISGetIndices(pcis->is_I_local,&idxs);CHKERRQ(ierr);
465b1b3d7a2SStefano Zampini   for (i=0;i<pcis->n-pcis->n_B;i++) {
466b1b3d7a2SStefano Zampini     ierr = PetscBTSet(bitmask,idxs[i]);CHKERRQ(ierr);
467b1b3d7a2SStefano Zampini   }
468b1b3d7a2SStefano Zampini   ierr = ISRestoreIndices(pcis->is_I_local,&idxs);CHKERRQ(ierr);
469b1b3d7a2SStefano Zampini 
470b1b3d7a2SStefano Zampini   for (i=0;i<sub_schurs->n_subs;i++) {
471b1b3d7a2SStefano Zampini     PetscInt subset_size,j;
472b1b3d7a2SStefano Zampini 
473b1b3d7a2SStefano Zampini     ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr);
474b1b3d7a2SStefano Zampini     ierr = ISGetIndices(sub_schurs->is_subs[i],&idxs);CHKERRQ(ierr);
475b1b3d7a2SStefano Zampini     for (j=0;j<subset_size;j++) {
476b1b3d7a2SStefano Zampini       ierr = PetscBTSet(bitmask,idxs[j]);CHKERRQ(ierr);
477b1b3d7a2SStefano Zampini     }
478b1b3d7a2SStefano Zampini     ierr = ISRestoreIndices(sub_schurs->is_subs[i],&idxs);CHKERRQ(ierr);
479b1b3d7a2SStefano Zampini   }
480b1b3d7a2SStefano Zampini 
481b1b3d7a2SStefano Zampini   deluxe_ctx->n_simple = 0;
482b1b3d7a2SStefano Zampini   for (i=0;i<pcis->n;i++) {
483b1b3d7a2SStefano Zampini     if (!PetscBTLookup(bitmask,i)) {
484b1b3d7a2SStefano Zampini       deluxe_ctx->n_simple++;
485b1b3d7a2SStefano Zampini     }
486b1b3d7a2SStefano Zampini   }
487b1b3d7a2SStefano Zampini   ierr = PetscMalloc1(deluxe_ctx->n_simple,&deluxe_ctx->idx_simple_B);CHKERRQ(ierr);
488b1b3d7a2SStefano Zampini   deluxe_ctx->n_simple = 0;
489b1b3d7a2SStefano Zampini   for (i=0;i<pcis->n;i++) {
490b1b3d7a2SStefano Zampini     if (!PetscBTLookup(bitmask,i)) {
491b1b3d7a2SStefano Zampini       deluxe_ctx->idx_simple_B[deluxe_ctx->n_simple++] = i;
492b1b3d7a2SStefano Zampini     }
493b1b3d7a2SStefano Zampini   }
494b1b3d7a2SStefano Zampini   ierr = ISGlobalToLocalMappingApply(pcis->BtoNmap,IS_GTOLM_DROP,deluxe_ctx->n_simple,deluxe_ctx->idx_simple_B,&i,deluxe_ctx->idx_simple_B);CHKERRQ(ierr);
495b1b3d7a2SStefano Zampini   if (i != deluxe_ctx->n_simple) {
496b1b3d7a2SStefano Zampini     SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Error when mapping simple scaling dofs! %d != %d",i,deluxe_ctx->n_simple);
497b1b3d7a2SStefano Zampini   }
498b1b3d7a2SStefano Zampini   ierr = PetscBTDestroy(&bitmask);CHKERRQ(ierr);
499b1b3d7a2SStefano Zampini 
50034a97f8cSStefano Zampini   /* free graph struct */
50134a97f8cSStefano Zampini   if (pcbddc->deluxe_rebuild) {
50234a97f8cSStefano Zampini     ierr = PCBDDCGraphDestroy(&graph);CHKERRQ(ierr);
50334a97f8cSStefano Zampini   }
50434a97f8cSStefano Zampini   PetscFunctionReturn(0);
50534a97f8cSStefano Zampini }
50634a97f8cSStefano Zampini 
50734a97f8cSStefano Zampini #undef __FUNCT__
50834a97f8cSStefano Zampini #define __FUNCT__ "PCBDDCScalingSetUp_Deluxe_Par"
50934a97f8cSStefano Zampini static PetscErrorCode PCBDDCScalingSetUp_Deluxe_Par(PC pc, PetscInt n_local_parallel_problems,PetscInt n_parallel_problems,PetscInt global_parallel[],PetscInt index_parallel[])
51034a97f8cSStefano Zampini {
51134a97f8cSStefano Zampini   PC_BDDC                *pcbddc=(PC_BDDC*)pc->data;
51234a97f8cSStefano Zampini   PCBDDCDeluxeScaling    deluxe_ctx=pcbddc->deluxe_ctx;
51334a97f8cSStefano Zampini   /* coloring */
51434a97f8cSStefano Zampini   Mat                    parallel_problems;
51534a97f8cSStefano Zampini   MatColoring            coloring_obj;
51634a97f8cSStefano Zampini   ISColoring             coloring_parallel_problems;
51734a97f8cSStefano Zampini   IS                     *par_is_colors,*is_colors;
51834a97f8cSStefano Zampini   /* working stuff */
51934a97f8cSStefano Zampini   PetscInt               i,j;
52034a97f8cSStefano Zampini   PetscErrorCode         ierr;
52134a97f8cSStefano Zampini 
52234a97f8cSStefano Zampini   PetscFunctionBegin;
52334a97f8cSStefano Zampini   if (!n_parallel_problems) {
52434a97f8cSStefano Zampini     PetscFunctionReturn(0);
52534a97f8cSStefano Zampini   }
52634a97f8cSStefano Zampini   /* Color parallel subproblems */
52734a97f8cSStefano Zampini   ierr = MatCreate(PetscObjectComm((PetscObject)pc),&parallel_problems);CHKERRQ(ierr);
52834a97f8cSStefano Zampini   ierr = MatSetSizes(parallel_problems,PETSC_DECIDE,PETSC_DECIDE,n_parallel_problems,n_parallel_problems);CHKERRQ(ierr);
52934a97f8cSStefano Zampini   ierr = MatSetType(parallel_problems,MATAIJ);CHKERRQ(ierr);
53034a97f8cSStefano Zampini   ierr = MatSetUp(parallel_problems);CHKERRQ(ierr);
53134a97f8cSStefano Zampini   ierr = MatSetOption(parallel_problems,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_FALSE);CHKERRQ(ierr);
53234a97f8cSStefano Zampini   ierr = MatSetOption(parallel_problems,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_FALSE);CHKERRQ(ierr);
53334a97f8cSStefano Zampini   for (i=0;i<n_local_parallel_problems;i++) {
53434a97f8cSStefano Zampini     PetscInt row = global_parallel[i];
53534a97f8cSStefano Zampini     for (j=0;j<n_local_parallel_problems;j++) {
53634a97f8cSStefano Zampini       PetscInt col = global_parallel[j];
53734a97f8cSStefano Zampini       if (row != col) {
53834a97f8cSStefano Zampini         ierr = MatSetValue(parallel_problems,row,col,1.0,INSERT_VALUES);CHKERRQ(ierr);
53934a97f8cSStefano Zampini       }
54034a97f8cSStefano Zampini     }
54134a97f8cSStefano Zampini   }
54234a97f8cSStefano Zampini   ierr = MatAssemblyBegin(parallel_problems,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
54334a97f8cSStefano Zampini   ierr = MatAssemblyEnd(parallel_problems,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
54434a97f8cSStefano Zampini   if (pcbddc->dbg_flag > 1) {
54534a97f8cSStefano Zampini     ierr = PetscViewerFlush(pcbddc->dbg_viewer);CHKERRQ(ierr);
54634a97f8cSStefano Zampini     ierr = PetscViewerASCIIPrintf(pcbddc->dbg_viewer,"Adj matrix for deluxe parallel problems\n");CHKERRQ(ierr);
54734a97f8cSStefano Zampini     ierr = MatView(parallel_problems,pcbddc->dbg_viewer);CHKERRQ(ierr);
54834a97f8cSStefano Zampini   }
54934a97f8cSStefano Zampini   ierr = MatColoringCreate(parallel_problems,&coloring_obj);CHKERRQ(ierr);
55034a97f8cSStefano Zampini   ierr = MatColoringSetDistance(coloring_obj,1);CHKERRQ(ierr);
55134a97f8cSStefano Zampini   ierr = MatColoringSetType(coloring_obj,MATCOLORINGJP);CHKERRQ(ierr);
55234a97f8cSStefano Zampini   ierr = MatColoringApply(coloring_obj,&coloring_parallel_problems);CHKERRQ(ierr);
55334a97f8cSStefano Zampini   ierr = ISColoringGetIS(coloring_parallel_problems,&deluxe_ctx->par_colors,&par_is_colors);CHKERRQ(ierr);
55434a97f8cSStefano Zampini   if (pcbddc->dbg_flag) {
55534a97f8cSStefano Zampini     ierr = PetscViewerFlush(pcbddc->dbg_viewer);CHKERRQ(ierr);
55634a97f8cSStefano Zampini     ierr = PetscViewerASCIIPrintf(pcbddc->dbg_viewer,"Number of colors %d for parallel part of deluxe\n",deluxe_ctx->par_colors);CHKERRQ(ierr);
55734a97f8cSStefano Zampini   }
55834a97f8cSStefano Zampini 
55934a97f8cSStefano Zampini   /* all procs should know the color distribution */
56034a97f8cSStefano Zampini   ierr = PetscMalloc1(deluxe_ctx->par_colors,&is_colors);CHKERRQ(ierr);
56134a97f8cSStefano Zampini   for (i=0;i<deluxe_ctx->par_colors;i++) {
56234a97f8cSStefano Zampini     if (pcbddc->dbg_flag) {
56334a97f8cSStefano Zampini       ierr = PetscViewerASCIIPrintf(pcbddc->dbg_viewer,"Global problem indexes for color %d\n",i);CHKERRQ(ierr);
56434a97f8cSStefano Zampini       ierr = ISView(par_is_colors[i],pcbddc->dbg_viewer);CHKERRQ(ierr);
56534a97f8cSStefano Zampini       ierr = PetscViewerFlush(pcbddc->dbg_viewer);CHKERRQ(ierr);
56634a97f8cSStefano Zampini     }
56734a97f8cSStefano Zampini     ierr = ISAllGather(par_is_colors[i],&is_colors[i]);CHKERRQ(ierr);
56834a97f8cSStefano Zampini   }
56934a97f8cSStefano Zampini 
57034a97f8cSStefano Zampini   /* free unneeded objects */
57134a97f8cSStefano Zampini   ierr = ISColoringRestoreIS(coloring_parallel_problems,&par_is_colors);CHKERRQ(ierr);
57234a97f8cSStefano Zampini   ierr = ISColoringDestroy(&coloring_parallel_problems);CHKERRQ(ierr);
57334a97f8cSStefano Zampini   ierr = MatColoringDestroy(&coloring_obj);CHKERRQ(ierr);
57434a97f8cSStefano Zampini   ierr = MatDestroy(&parallel_problems);CHKERRQ(ierr);
57534a97f8cSStefano Zampini 
57634a97f8cSStefano Zampini   /* allocate deluxe arrays for parallel problems */
57734a97f8cSStefano Zampini   ierr = PetscMalloc5(deluxe_ctx->par_colors,&deluxe_ctx->par_ksp,
57834a97f8cSStefano Zampini                       deluxe_ctx->par_colors,&deluxe_ctx->par_scctx_s,
57934a97f8cSStefano Zampini                       deluxe_ctx->par_colors,&deluxe_ctx->par_scctx_p,
58034a97f8cSStefano Zampini                       deluxe_ctx->par_colors,&deluxe_ctx->par_vec,
58134a97f8cSStefano Zampini                       deluxe_ctx->par_colors,&deluxe_ctx->par_col2sub);CHKERRQ(ierr);
58234a97f8cSStefano Zampini 
58334a97f8cSStefano Zampini   /* cycle on colors */
58434a97f8cSStefano Zampini   for (i=0;i<deluxe_ctx->par_colors;i++) {
58534a97f8cSStefano Zampini     PetscSubcomm    par_subcomm;
58634a97f8cSStefano Zampini     const PetscInt* idxs_subproblems;
58734a97f8cSStefano Zampini     PetscInt        color_size;
58834a97f8cSStefano Zampini     PetscMPIInt     rank,active_color;
58934a97f8cSStefano Zampini 
59034a97f8cSStefano Zampini     /* get local index of i-th parallel colored problem */
59134a97f8cSStefano Zampini     ierr = ISGetLocalSize(is_colors[i],&color_size);CHKERRQ(ierr);
59234a97f8cSStefano Zampini     ierr = ISGetIndices(is_colors[i],&idxs_subproblems);CHKERRQ(ierr);
59334a97f8cSStefano Zampini     /* split comm for computing parallel problems for this color */
59434a97f8cSStefano Zampini     /* Processes not partecipating at this stage will have color = color_size */
59534a97f8cSStefano Zampini     /* because PetscCommDuplicate does not handle MPI_COMM_NULL */
59634a97f8cSStefano Zampini     active_color = color_size;
59734a97f8cSStefano Zampini     deluxe_ctx->par_col2sub[i] = -1;
59834a97f8cSStefano Zampini     for (j=0;j<n_local_parallel_problems;j++) {
59934a97f8cSStefano Zampini       PetscInt local_idx;
60034a97f8cSStefano Zampini       ierr = PetscFindInt(global_parallel[j],color_size,idxs_subproblems,&local_idx);CHKERRQ(ierr);
60134a97f8cSStefano Zampini       if (local_idx > -1) {
60234a97f8cSStefano Zampini         ierr = PetscMPIIntCast(local_idx,&active_color);CHKERRQ(ierr);
60334a97f8cSStefano Zampini         deluxe_ctx->par_col2sub[i] = index_parallel[j];
60434a97f8cSStefano Zampini         break;
60534a97f8cSStefano Zampini       }
60634a97f8cSStefano Zampini     }
60734a97f8cSStefano Zampini     ierr = ISRestoreIndices(is_colors[i],&idxs_subproblems);CHKERRQ(ierr);
60834a97f8cSStefano Zampini     ierr = PetscSubcommCreate(PetscObjectComm((PetscObject)pc),&par_subcomm);CHKERRQ(ierr);
60934a97f8cSStefano Zampini     ierr = PetscSubcommSetNumber(par_subcomm,color_size+1);CHKERRQ(ierr);
61034a97f8cSStefano Zampini     ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)pc),&rank);CHKERRQ(ierr);
61134a97f8cSStefano Zampini     ierr = PetscSubcommSetTypeGeneral(par_subcomm,active_color,rank);CHKERRQ(ierr);
61234a97f8cSStefano Zampini     /* print debug info */
61334a97f8cSStefano Zampini     if (pcbddc->dbg_flag) {
61434a97f8cSStefano Zampini       PetscMPIInt crank,csize;
61534a97f8cSStefano Zampini       ierr = MPI_Comm_rank(par_subcomm->comm,&crank);CHKERRQ(ierr);
61634a97f8cSStefano Zampini       ierr = MPI_Comm_size(par_subcomm->comm,&csize);CHKERRQ(ierr);
61734a97f8cSStefano Zampini       ierr = PetscViewerASCIIPrintf(pcbddc->dbg_viewer,"Color %d: size %d, details follows.\n",i,color_size);CHKERRQ(ierr);
61834a97f8cSStefano Zampini       ierr = PetscViewerFlush(pcbddc->dbg_viewer);CHKERRQ(ierr);
61934a97f8cSStefano Zampini       ierr = PetscViewerASCIISynchronizedAllow(pcbddc->dbg_viewer,PETSC_TRUE);CHKERRQ(ierr);
62034a97f8cSStefano 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);
62134a97f8cSStefano Zampini       ierr = PetscViewerFlush(pcbddc->dbg_viewer);CHKERRQ(ierr);
62234a97f8cSStefano Zampini     }
62334a97f8cSStefano Zampini 
62434a97f8cSStefano Zampini     if (deluxe_ctx->par_col2sub[i] >= 0) {
6255e8657edSStefano Zampini       PC                     pctemp;
6265e8657edSStefano Zampini       PC_IS                  *pcis=(PC_IS*)pc->data;
62734a97f8cSStefano Zampini       Mat                    color_mat,color_mat_is,temp_mat;
62834a97f8cSStefano Zampini       ISLocalToGlobalMapping WtoNmap,l2gmap_subset;
62934a97f8cSStefano Zampini       IS                     is_local_numbering,isB_local,isW_local,isW;
630b1b3d7a2SStefano Zampini       PCBDDCSubSchurs        sub_schurs = pcbddc->sub_schurs[1];
63134a97f8cSStefano Zampini       PetscInt               subidx,n_local_dofs,n_global_dofs;
63234a97f8cSStefano Zampini       PetscInt               *global_numbering,*local_numbering;
63334a97f8cSStefano Zampini       char                   ksp_prefix[256];
63434a97f8cSStefano Zampini       size_t                 len;
63534a97f8cSStefano Zampini 
63634a97f8cSStefano Zampini       /* Local index for schur complement on subset */
63734a97f8cSStefano Zampini       subidx = deluxe_ctx->par_col2sub[i];
63834a97f8cSStefano Zampini 
63934a97f8cSStefano Zampini       /* Parallel numbering for dofs in colored subset */
64068270318SStefano Zampini       ierr = ISSum(sub_schurs->is_I_layer,sub_schurs->is_AEj_B[subidx],&is_local_numbering);CHKERRQ(ierr);
64134a97f8cSStefano Zampini       ierr = ISGetLocalSize(is_local_numbering,&n_local_dofs);CHKERRQ(ierr);
64234a97f8cSStefano Zampini       ierr = ISGetIndices(is_local_numbering,(const PetscInt **)&local_numbering);CHKERRQ(ierr);
64334a97f8cSStefano Zampini       ierr = PCBDDCSubsetNumbering(par_subcomm->comm,pcbddc->mat_graph->l2gmap,n_local_dofs,local_numbering,PETSC_NULL,&n_global_dofs,&global_numbering);CHKERRQ(ierr);
64434a97f8cSStefano Zampini       ierr = ISRestoreIndices(is_local_numbering,(const PetscInt **)&local_numbering);CHKERRQ(ierr);
64534a97f8cSStefano Zampini 
64634a97f8cSStefano Zampini       /* L2Gmap from relevant dofs to local dofs */
64734a97f8cSStefano Zampini       ierr = ISLocalToGlobalMappingCreateIS(is_local_numbering,&WtoNmap);CHKERRQ(ierr);
64834a97f8cSStefano Zampini 
64934a97f8cSStefano Zampini       /* L2Gmap from local to global dofs */
65034a97f8cSStefano Zampini       ierr = ISLocalToGlobalMappingCreate(par_subcomm->comm,1,n_local_dofs,global_numbering,PETSC_COPY_VALUES,&l2gmap_subset);CHKERRQ(ierr);
65134a97f8cSStefano Zampini 
65234a97f8cSStefano Zampini       /* compute parallel matrix (extended dirichlet problem on subset) */
65334a97f8cSStefano Zampini       ierr = MatCreateIS(par_subcomm->comm,1,PETSC_DECIDE,PETSC_DECIDE,n_global_dofs,n_global_dofs,l2gmap_subset,&color_mat_is);CHKERRQ(ierr);
65434a97f8cSStefano Zampini       ierr = MatGetSubMatrix(pcbddc->local_mat,is_local_numbering,is_local_numbering,MAT_INITIAL_MATRIX,&temp_mat);CHKERRQ(ierr);
65534a97f8cSStefano Zampini       ierr = MatISSetLocalMat(color_mat_is,temp_mat);CHKERRQ(ierr);
65634a97f8cSStefano Zampini       ierr = MatDestroy(&temp_mat);CHKERRQ(ierr);
65734a97f8cSStefano Zampini       ierr = MatISGetMPIXAIJ(color_mat_is,MAT_INITIAL_MATRIX,&color_mat);CHKERRQ(ierr);
65834a97f8cSStefano Zampini       ierr = MatDestroy(&color_mat_is);CHKERRQ(ierr);
65934a97f8cSStefano Zampini 
66034a97f8cSStefano Zampini       /* work vector for (parallel) extended dirichlet problem */
6618a26ef87SStefano Zampini       ierr = MatCreateVecs(color_mat,&deluxe_ctx->par_vec[i],NULL);CHKERRQ(ierr);
66234a97f8cSStefano Zampini 
66334a97f8cSStefano Zampini       /* compute scatters */
66434a97f8cSStefano Zampini       /* deluxe_ctx->par_scctx_p[i] extension from local subset to extended dirichlet problem
66534a97f8cSStefano Zampini          deluxe_ctx->par_scctx_s[i] restriction from local boundary to subset -> simple copy of selected values */
6665e8657edSStefano Zampini       ierr = ISGlobalToLocalMappingApplyIS(pcis->BtoNmap,IS_GTOLM_DROP,sub_schurs->is_AEj_B[subidx],&isB_local);CHKERRQ(ierr);
66734a97f8cSStefano Zampini       ierr = VecScatterCreate(pcbddc->work_scaling,isB_local,sub_schurs->work1[subidx],NULL,&deluxe_ctx->par_scctx_s[i]);CHKERRQ(ierr);
66834a97f8cSStefano Zampini       ierr = ISGlobalToLocalMappingApplyIS(WtoNmap,IS_GTOLM_DROP,sub_schurs->is_AEj_B[subidx],&isW_local);CHKERRQ(ierr);
66934a97f8cSStefano Zampini       ierr = ISLocalToGlobalMappingApplyIS(l2gmap_subset,isW_local,&isW);CHKERRQ(ierr);
67034a97f8cSStefano Zampini       ierr = VecScatterCreate(sub_schurs->work1[subidx],NULL,deluxe_ctx->par_vec[i],isW,&deluxe_ctx->par_scctx_p[i]);CHKERRQ(ierr);
67134a97f8cSStefano Zampini 
67234a97f8cSStefano Zampini       /* free objects no longer neeeded */
67334a97f8cSStefano Zampini       ierr = ISDestroy(&isW);CHKERRQ(ierr);
67434a97f8cSStefano Zampini       ierr = ISDestroy(&isW_local);CHKERRQ(ierr);
67534a97f8cSStefano Zampini       ierr = ISDestroy(&isB_local);CHKERRQ(ierr);
67634a97f8cSStefano Zampini       ierr = ISLocalToGlobalMappingDestroy(&WtoNmap);CHKERRQ(ierr);
67734a97f8cSStefano Zampini       ierr = ISLocalToGlobalMappingDestroy(&l2gmap_subset);CHKERRQ(ierr);
67834a97f8cSStefano Zampini       ierr = ISDestroy(&is_local_numbering);CHKERRQ(ierr);
67934a97f8cSStefano Zampini       ierr = PetscFree(global_numbering);CHKERRQ(ierr);
68034a97f8cSStefano Zampini 
68134a97f8cSStefano Zampini       /* KSP for extended dirichlet problem */
68234a97f8cSStefano Zampini       ierr = KSPCreate(par_subcomm->comm,&deluxe_ctx->par_ksp[i]);CHKERRQ(ierr);
68334a97f8cSStefano Zampini       ierr = KSPSetOperators(deluxe_ctx->par_ksp[i],color_mat,color_mat);CHKERRQ(ierr);
68434a97f8cSStefano Zampini       ierr = KSPSetTolerances(deluxe_ctx->par_ksp[i],1.e-12,1.e-12,1.e10,10000);CHKERRQ(ierr);
68534a97f8cSStefano Zampini       ierr = KSPSetType(deluxe_ctx->par_ksp[i],KSPPREONLY);CHKERRQ(ierr);
6865e8657edSStefano Zampini       ierr = KSPGetPC(deluxe_ctx->par_ksp[i],&pctemp);CHKERRQ(ierr);
6875e8657edSStefano Zampini       ierr = PCSetType(pctemp,PCREDUNDANT);CHKERRQ(ierr);
68834a97f8cSStefano Zampini       ierr = PetscStrlen(((PetscObject)(pcbddc->ksp_D))->prefix,&len);CHKERRQ(ierr);
6898856534dSStefano Zampini       len -= 10; /* remove "dirichlet_" */
6908856534dSStefano Zampini       ierr = PetscStrncpy(ksp_prefix,((PetscObject)(pcbddc->ksp_D))->prefix,len+1);CHKERRQ(ierr); /* PetscStrncpy puts a terminating char at the end */
69134a97f8cSStefano Zampini       ierr = PetscStrcat(ksp_prefix,"deluxe_par_");CHKERRQ(ierr);
69234a97f8cSStefano Zampini       ierr = KSPSetOptionsPrefix(deluxe_ctx->par_ksp[i],ksp_prefix);CHKERRQ(ierr);
69334a97f8cSStefano Zampini       ierr = KSPSetFromOptions(deluxe_ctx->par_ksp[i]);CHKERRQ(ierr);
69434a97f8cSStefano Zampini       ierr = KSPSetUp(deluxe_ctx->par_ksp[i]);CHKERRQ(ierr);
69534a97f8cSStefano Zampini       ierr = MatDestroy(&color_mat);CHKERRQ(ierr);
69634a97f8cSStefano Zampini     } else { /* not partecipating in color */
69734a97f8cSStefano Zampini       deluxe_ctx->par_ksp[i] = 0;
69834a97f8cSStefano Zampini       deluxe_ctx->par_vec[i] = 0;
69934a97f8cSStefano Zampini       deluxe_ctx->par_scctx_p[i] = 0;
70034a97f8cSStefano Zampini       deluxe_ctx->par_scctx_s[i] = 0;
70134a97f8cSStefano Zampini     }
70234a97f8cSStefano Zampini     ierr = PetscSubcommDestroy(&par_subcomm);CHKERRQ(ierr);
70334a97f8cSStefano Zampini   }
70434a97f8cSStefano Zampini   for (i=0;i<deluxe_ctx->par_colors;i++) {
70534a97f8cSStefano Zampini     ierr = ISDestroy(&is_colors[i]);CHKERRQ(ierr);
70634a97f8cSStefano Zampini   }
70734a97f8cSStefano Zampini   ierr = PetscFree(is_colors);CHKERRQ(ierr);
70834a97f8cSStefano Zampini 
70934a97f8cSStefano Zampini   if (pcbddc->dbg_flag) {
71034a97f8cSStefano Zampini     Vec test_vec;
71134a97f8cSStefano Zampini     PetscReal error;
712b1b3d7a2SStefano Zampini     PCBDDCSubSchurs sub_schurs = pcbddc->sub_schurs[1];
71334a97f8cSStefano Zampini     /* test partition of unity of coloured schur complements  */
71434a97f8cSStefano Zampini     for (i=0;i<deluxe_ctx->par_colors;i++) {
71534a97f8cSStefano Zampini       PetscInt  subidx = deluxe_ctx->par_col2sub[i];
71634a97f8cSStefano Zampini       PetscBool error_found = PETSC_FALSE;
71734a97f8cSStefano Zampini       ierr = PetscViewerASCIISynchronizedAllow(pcbddc->dbg_viewer,PETSC_TRUE);CHKERRQ(ierr);
71834a97f8cSStefano Zampini 
71934a97f8cSStefano Zampini       if (deluxe_ctx->par_ksp[i]) {
72034a97f8cSStefano Zampini         /* create random test vec being zero on internal nodes of the extende dirichlet problem */
72134a97f8cSStefano Zampini         ierr = VecDuplicate(deluxe_ctx->par_vec[i],&test_vec);CHKERRQ(ierr);
72234a97f8cSStefano Zampini         ierr = VecSetRandom(sub_schurs->work1[subidx],PETSC_NULL);CHKERRQ(ierr);
72334a97f8cSStefano Zampini         ierr = VecSet(test_vec,0.0);CHKERRQ(ierr);
72434a97f8cSStefano Zampini         ierr = VecScatterBegin(deluxe_ctx->par_scctx_p[i],sub_schurs->work1[subidx],test_vec,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
72534a97f8cSStefano Zampini         ierr = VecScatterEnd(deluxe_ctx->par_scctx_p[i],sub_schurs->work1[subidx],test_vec,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
72634a97f8cSStefano Zampini         /* w_j */
72734a97f8cSStefano Zampini         ierr = VecScatterBegin(deluxe_ctx->par_scctx_p[i],test_vec,sub_schurs->work1[subidx],INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
72834a97f8cSStefano Zampini         ierr = VecScatterEnd(deluxe_ctx->par_scctx_p[i],test_vec,sub_schurs->work1[subidx],INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
72934a97f8cSStefano Zampini         /* S_j*w_j */
73034a97f8cSStefano Zampini         ierr = MatMult(sub_schurs->S_Ej[subidx],sub_schurs->work1[subidx],sub_schurs->work2[subidx]);CHKERRQ(ierr);
73134a97f8cSStefano Zampini         /* \sum_j S_j*w_j */
73234a97f8cSStefano Zampini         ierr = VecSet(deluxe_ctx->par_vec[i],0.0);CHKERRQ(ierr);
73334a97f8cSStefano Zampini         ierr = VecScatterBegin(deluxe_ctx->par_scctx_p[i],sub_schurs->work2[subidx],deluxe_ctx->par_vec[i],ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
73434a97f8cSStefano Zampini         ierr = VecScatterEnd(deluxe_ctx->par_scctx_p[i],sub_schurs->work2[subidx],deluxe_ctx->par_vec[i],ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
73534a97f8cSStefano Zampini         /* (\sum_j S_j)^(-1)(\sum_j S_j*w_j) */
73634a97f8cSStefano Zampini         ierr = KSPSolve(deluxe_ctx->par_ksp[i],deluxe_ctx->par_vec[i],deluxe_ctx->par_vec[i]);CHKERRQ(ierr);
73734a97f8cSStefano Zampini         ierr = VecScatterBegin(deluxe_ctx->par_scctx_p[i],deluxe_ctx->par_vec[i],sub_schurs->work1[subidx],INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
73834a97f8cSStefano Zampini         ierr = VecScatterEnd(deluxe_ctx->par_scctx_p[i],deluxe_ctx->par_vec[i],sub_schurs->work1[subidx],INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
73934a97f8cSStefano Zampini         ierr = VecSet(deluxe_ctx->par_vec[i],0.0);CHKERRQ(ierr);
74034a97f8cSStefano Zampini         ierr = VecScatterBegin(deluxe_ctx->par_scctx_p[i],sub_schurs->work1[subidx],deluxe_ctx->par_vec[i],INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
74134a97f8cSStefano Zampini         ierr = VecScatterEnd(deluxe_ctx->par_scctx_p[i],sub_schurs->work1[subidx],deluxe_ctx->par_vec[i],INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
74234a97f8cSStefano Zampini         /* test partition of unity */
74334a97f8cSStefano Zampini         ierr = VecAXPY(test_vec,-1.0,deluxe_ctx->par_vec[i]);CHKERRQ(ierr);
74434a97f8cSStefano Zampini         ierr = VecNorm(test_vec,NORM_INFINITY,&error);CHKERRQ(ierr);
745c63f45b2SStefano Zampini         if (PetscAbsReal(error) > 1.e-2) {
74634a97f8cSStefano Zampini           /* ierr = VecView(test_vec,0);CHKERRQ(ierr); */
74734a97f8cSStefano Zampini           error_found = PETSC_TRUE;
74834a97f8cSStefano Zampini         }
74934a97f8cSStefano Zampini         ierr = VecDestroy(&test_vec);CHKERRQ(ierr);
75034a97f8cSStefano Zampini       }
75134a97f8cSStefano Zampini       if (error_found) {
75234a97f8cSStefano Zampini         ierr = PetscViewerASCIISynchronizedPrintf(pcbddc->dbg_viewer,"Error testing local schur for color %d and subdomain %d\n",i,PetscGlobalRank);CHKERRQ(ierr);
75334a97f8cSStefano Zampini       }
75434a97f8cSStefano Zampini       ierr = PetscViewerFlush(pcbddc->dbg_viewer);CHKERRQ(ierr);
75534a97f8cSStefano Zampini     }
75634a97f8cSStefano Zampini   }
75734a97f8cSStefano Zampini   PetscFunctionReturn(0);
75834a97f8cSStefano Zampini }
75934a97f8cSStefano Zampini 
7605db18549SStefano Zampini #undef __FUNCT__
761883469d8SStefano Zampini #define __FUNCT__ "PCBDDCScalingSetUp_Deluxe_Seq"
762883469d8SStefano Zampini static PetscErrorCode PCBDDCScalingSetUp_Deluxe_Seq(PC pc)
7635db18549SStefano Zampini {
7645db18549SStefano Zampini   PC_BDDC                *pcbddc=(PC_BDDC*)pc->data;
7655db18549SStefano Zampini   PCBDDCDeluxeScaling    deluxe_ctx=pcbddc->deluxe_ctx;
7665db18549SStefano Zampini   PCBDDCSubSchurs        sub_schurs = pcbddc->sub_schurs[1];
7675db18549SStefano Zampini   PC                     pc_temp;
7685db18549SStefano Zampini   MatSolverPackage       solver=NULL;
7695db18549SStefano Zampini   char                   ksp_prefix[256];
7705db18549SStefano Zampini   size_t                 len;
7715db18549SStefano Zampini   PetscInt               local_size;
7725db18549SStefano Zampini   PetscErrorCode         ierr;
7735db18549SStefano Zampini 
7745db18549SStefano Zampini   PetscFunctionBegin;
7759221af80SStefano Zampini   if (!sub_schurs->n_subs_seq_g) {
7769221af80SStefano Zampini     PetscFunctionReturn(0);
7779221af80SStefano Zampini   }
7789221af80SStefano Zampini 
7795db18549SStefano Zampini   /* Create work vectors for sequential part of deluxe */
780*65d8bf0aSStefano Zampini   ierr = MatCreateVecs(sub_schurs->sum_S_Ej_all,&deluxe_ctx->seq_work1,&deluxe_ctx->seq_work2);CHKERRQ(ierr);
7815db18549SStefano Zampini 
7825db18549SStefano Zampini   /* Compute deluxe sequential scatter */
7835db18549SStefano Zampini   ierr = VecScatterCreate(pcbddc->work_scaling,sub_schurs->is_Ej_all,deluxe_ctx->seq_work1,NULL,&deluxe_ctx->seq_scctx);CHKERRQ(ierr);
7845db18549SStefano Zampini 
7855db18549SStefano Zampini   /* Create KSP object for sequential part of deluxe scaling */
7865db18549SStefano Zampini   ierr = KSPCreate(PETSC_COMM_SELF,&deluxe_ctx->seq_ksp);CHKERRQ(ierr);
7875db18549SStefano Zampini   ierr = KSPSetOperators(deluxe_ctx->seq_ksp,sub_schurs->sum_S_Ej_all,sub_schurs->sum_S_Ej_all);CHKERRQ(ierr);
7885db18549SStefano Zampini   ierr = KSPSetType(deluxe_ctx->seq_ksp,KSPPREONLY);CHKERRQ(ierr);
7895db18549SStefano Zampini   ierr = KSPGetPC(deluxe_ctx->seq_ksp,&pc_temp);CHKERRQ(ierr);
7905db18549SStefano Zampini   ierr = PCSetType(pc_temp,PCLU);CHKERRQ(ierr);
7915db18549SStefano Zampini   ierr = KSPGetPC(pcbddc->ksp_D,&pc_temp);CHKERRQ(ierr);
7925db18549SStefano Zampini   ierr = PCFactorGetMatSolverPackage(pc_temp,(const MatSolverPackage*)&solver);CHKERRQ(ierr);
7935db18549SStefano Zampini   ierr = MatGetSize(sub_schurs->sum_S_Ej_all,&local_size,NULL);CHKERRQ(ierr);
7945db18549SStefano Zampini   if (solver && local_size) { /* if local_size is null, some external packages will report errors */
7955db18549SStefano Zampini     PC     new_pc;
7965db18549SStefano Zampini     PCType type;
7975db18549SStefano Zampini     ierr = PCGetType(pc_temp,&type);CHKERRQ(ierr);
7985db18549SStefano Zampini     ierr = KSPGetPC(deluxe_ctx->seq_ksp,&new_pc);CHKERRQ(ierr);
7995db18549SStefano Zampini     ierr = PCSetType(new_pc,type);CHKERRQ(ierr);
8005db18549SStefano Zampini     ierr = PCFactorSetMatSolverPackage(new_pc,solver);CHKERRQ(ierr);
8015db18549SStefano Zampini   }
8025db18549SStefano Zampini   ierr = PetscStrlen(((PetscObject)(pcbddc->ksp_D))->prefix,&len);CHKERRQ(ierr);
8035db18549SStefano Zampini   len -= 10; /* remove "dirichlet_" */
8045db18549SStefano Zampini   ierr = PetscStrncpy(ksp_prefix,((PetscObject)(pcbddc->ksp_D))->prefix,len+1);CHKERRQ(ierr);
8055db18549SStefano Zampini   ierr = PetscStrcat(ksp_prefix,"deluxe_seq_");CHKERRQ(ierr);
8065db18549SStefano Zampini   ierr = KSPSetOptionsPrefix(deluxe_ctx->seq_ksp,ksp_prefix);CHKERRQ(ierr);
8075db18549SStefano Zampini   if (local_size) {
8085db18549SStefano Zampini     ierr = KSPSetFromOptions(deluxe_ctx->seq_ksp);CHKERRQ(ierr);
8095db18549SStefano Zampini   }
8105db18549SStefano Zampini   ierr = KSPSetUp(deluxe_ctx->seq_ksp);CHKERRQ(ierr);
8115db18549SStefano Zampini   PetscFunctionReturn(0);
8125db18549SStefano Zampini }
813