xref: /petsc/src/ksp/pc/impls/bddc/bddcscalingbasic.c (revision a4e17c6708a50153b2912a92c0eeebc47f7a1da7)
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[]);
934a97f8cSStefano Zampini static PetscErrorCode PCBDDCScalingSetUp_Deluxe_Seq(PC,PetscInt,PetscInt,PetscInt[],PetscInt[]);
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;
3634a97f8cSStefano Zampini   PCBDDCSubSchurs     sub_schurs = deluxe_ctx->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 */
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 */
57674ae819SStefano Zampini   if (deluxe_ctx->seq_mat) {
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);
60674ae819SStefano Zampini     ierr = MatMult(deluxe_ctx->seq_mat,deluxe_ctx->seq_work1,deluxe_ctx->seq_work2);CHKERRQ(ierr);
61674ae819SStefano Zampini     ierr = KSPSolve(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;
13834a97f8cSStefano Zampini   PCBDDCSubSchurs     sub_schurs = deluxe_ctx->sub_schurs;
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 */
158674ae819SStefano Zampini   if (deluxe_ctx->seq_mat) {
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);
161674ae819SStefano Zampini     ierr = KSPSolveTranspose(deluxe_ctx->seq_ksp,deluxe_ctx->seq_work1,deluxe_ctx->seq_work2);CHKERRQ(ierr);
162674ae819SStefano Zampini     ierr = MatMultTranspose(deluxe_ctx->seq_mat,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   ierr = PCBDDCSubSchursCreate(&deluxe_ctx->sub_schurs);CHKERRQ(ierr);
32034a97f8cSStefano Zampini   pcbddc->deluxe_ctx = deluxe_ctx;
32134a97f8cSStefano Zampini   PetscFunctionReturn(0);
32234a97f8cSStefano Zampini }
32334a97f8cSStefano Zampini 
32434a97f8cSStefano Zampini #undef __FUNCT__
32534a97f8cSStefano Zampini #define __FUNCT__ "PCBDDCScalingDestroy_Deluxe"
32634a97f8cSStefano Zampini static PetscErrorCode PCBDDCScalingDestroy_Deluxe(PC pc)
32734a97f8cSStefano Zampini {
32834a97f8cSStefano Zampini   PC_BDDC*            pcbddc=(PC_BDDC*)pc->data;
32934a97f8cSStefano Zampini   PetscErrorCode      ierr;
33034a97f8cSStefano Zampini 
33134a97f8cSStefano Zampini   PetscFunctionBegin;
33234a97f8cSStefano Zampini   ierr = PCBDDCScalingReset_Deluxe_Solvers(pcbddc->deluxe_ctx);CHKERRQ(ierr);
33334a97f8cSStefano Zampini   ierr = PCBDDCSubSchursDestroy(&(pcbddc->deluxe_ctx->sub_schurs));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;
34734a97f8cSStefano Zampini   if (deluxe_ctx->seq_mat) {
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 = MatDestroy(&deluxe_ctx->seq_mat);CHKERRQ(ierr);
35234a97f8cSStefano Zampini     ierr = KSPDestroy(&deluxe_ctx->seq_ksp);CHKERRQ(ierr);
35334a97f8cSStefano Zampini   }
35434a97f8cSStefano Zampini   if (deluxe_ctx->par_colors) {
35534a97f8cSStefano Zampini     PetscInt i;
35634a97f8cSStefano Zampini     for (i=0;i<deluxe_ctx->par_colors;i++) {
35734a97f8cSStefano Zampini       ierr = VecScatterDestroy(&deluxe_ctx->par_scctx_s[i]);CHKERRQ(ierr);
35834a97f8cSStefano Zampini       ierr = VecScatterDestroy(&deluxe_ctx->par_scctx_p[i]);CHKERRQ(ierr);
35934a97f8cSStefano Zampini       ierr = VecDestroy(&deluxe_ctx->par_vec[i]);CHKERRQ(ierr);
36034a97f8cSStefano Zampini       ierr = KSPDestroy(&deluxe_ctx->par_ksp[i]);CHKERRQ(ierr);
36134a97f8cSStefano Zampini     }
36234a97f8cSStefano Zampini     ierr = PetscFree5(deluxe_ctx->par_ksp,
36334a97f8cSStefano Zampini                       deluxe_ctx->par_scctx_s,
36434a97f8cSStefano Zampini                       deluxe_ctx->par_scctx_p,
36534a97f8cSStefano Zampini                       deluxe_ctx->par_vec,
36634a97f8cSStefano Zampini                       deluxe_ctx->par_col2sub);CHKERRQ(ierr);
36734a97f8cSStefano Zampini   }
36834a97f8cSStefano Zampini   deluxe_ctx->par_colors = 0;
36934a97f8cSStefano Zampini   PetscFunctionReturn(0);
37034a97f8cSStefano Zampini }
37134a97f8cSStefano Zampini 
37234a97f8cSStefano Zampini #undef __FUNCT__
37334a97f8cSStefano Zampini #define __FUNCT__ "PCBDDCScalingSetUp_Deluxe"
37434a97f8cSStefano Zampini static PetscErrorCode PCBDDCScalingSetUp_Deluxe(PC pc)
37534a97f8cSStefano Zampini {
37634a97f8cSStefano Zampini   PC_IS               *pcis=(PC_IS*)pc->data;
37734a97f8cSStefano Zampini   PC_BDDC             *pcbddc=(PC_BDDC*)pc->data;
37834a97f8cSStefano Zampini   PCBDDCDeluxeScaling deluxe_ctx=pcbddc->deluxe_ctx;
37934a97f8cSStefano Zampini   PCBDDCSubSchurs     sub_schurs=deluxe_ctx->sub_schurs;
38034a97f8cSStefano Zampini   PCBDDCGraph         graph;
3817821e9e7SStefano Zampini   IS                  *faces,*edges,*all_cc;
3827821e9e7SStefano Zampini   PetscBT             bitmask;
38334a97f8cSStefano Zampini   PetscInt            *index_sequential,*index_parallel;
38434a97f8cSStefano Zampini   PetscInt            *auxlocal_sequential,*auxlocal_parallel;
38534a97f8cSStefano Zampini   PetscInt            *auxglobal_sequential,*auxglobal_parallel;
3867821e9e7SStefano Zampini   PetscInt            *auxmapping,*idxs;
3871cef56d8SStefano Zampini   PetscInt            i,max_subset_size;
38834a97f8cSStefano Zampini   PetscInt            n_sequential_problems,n_local_sequential_problems,n_parallel_problems,n_local_parallel_problems;
3897821e9e7SStefano Zampini   PetscInt            n_faces,n_edges,n_all_cc;
39034a97f8cSStefano Zampini   PetscErrorCode      ierr;
39134a97f8cSStefano Zampini 
39234a97f8cSStefano Zampini   PetscFunctionBegin;
39334a97f8cSStefano Zampini   /* throw away the solvers */
39434a97f8cSStefano Zampini   ierr = PCBDDCScalingReset_Deluxe_Solvers(deluxe_ctx);CHKERRQ(ierr);
39534a97f8cSStefano Zampini 
39634a97f8cSStefano Zampini   /* attach interface graph for determining subsets */
39734a97f8cSStefano Zampini   if (pcbddc->deluxe_rebuild) { /* in case rebuild has been requested, it uses a graph generated only by the neighbouring information */
39834a97f8cSStefano Zampini     PetscInt *idx_V_N;
39934a97f8cSStefano Zampini     IS       verticesIS;
40034a97f8cSStefano Zampini     ierr = PCBDDCGetPrimalVerticesLocalIdx(pc,&i,&idx_V_N);CHKERRQ(ierr);
40134a97f8cSStefano Zampini     ierr = ISCreateGeneral(PETSC_COMM_SELF,i,idx_V_N,PETSC_OWN_POINTER,&verticesIS);CHKERRQ(ierr);
40234a97f8cSStefano Zampini     ierr = PCBDDCGraphCreate(&graph);CHKERRQ(ierr);
40334a97f8cSStefano Zampini     ierr = PCBDDCGraphInit(graph,pcbddc->mat_graph->l2gmap);CHKERRQ(ierr);
4040a95f60dSStefano Zampini     ierr = PCBDDCGraphSetUp(graph,0,NULL,pcbddc->DirichletBoundariesLocal,0,NULL,verticesIS);CHKERRQ(ierr);
40534a97f8cSStefano Zampini     ierr = PCBDDCGraphComputeConnectedComponents(graph);CHKERRQ(ierr);
40634a97f8cSStefano Zampini     ierr = ISDestroy(&verticesIS);CHKERRQ(ierr);
40734a97f8cSStefano Zampini /*
40834a97f8cSStefano Zampini     if (pcbddc->dbg_flag) {
40934a97f8cSStefano Zampini       ierr = PCBDDCGraphASCIIView(graph,pcbddc->dbg_flag,pcbddc->dbg_viewer);CHKERRQ(ierr);
41034a97f8cSStefano Zampini     }
41134a97f8cSStefano Zampini */
41234a97f8cSStefano Zampini   } else {
41334a97f8cSStefano Zampini     graph = pcbddc->mat_graph;
41434a97f8cSStefano Zampini   }
41534a97f8cSStefano Zampini 
4167821e9e7SStefano Zampini   /* get index sets for faces and edges */
4177821e9e7SStefano Zampini   ierr = PCBDDCGraphGetCandidatesIS(graph,&n_faces,&faces,&n_edges,&edges,NULL);CHKERRQ(ierr);
4187821e9e7SStefano Zampini   n_all_cc = n_faces+n_edges;
4197821e9e7SStefano Zampini   ierr = PetscMalloc1(n_all_cc,&all_cc);CHKERRQ(ierr);
4207821e9e7SStefano Zampini   for (i=0;i<n_faces;i++) {
4217821e9e7SStefano Zampini     all_cc[i] = faces[i];
4227821e9e7SStefano Zampini   }
4237821e9e7SStefano Zampini   for (i=0;i<n_edges;i++) {
4247821e9e7SStefano Zampini     all_cc[n_faces+i] = edges[i];
4257821e9e7SStefano Zampini   }
4267821e9e7SStefano Zampini 
42734a97f8cSStefano Zampini   /* map interface's subsets */
4281cef56d8SStefano Zampini   max_subset_size = 0;
4297821e9e7SStefano Zampini   for (i=0;i<n_all_cc;i++) {
4307821e9e7SStefano Zampini     PetscInt subset_size;
4317821e9e7SStefano Zampini     ierr = ISGetLocalSize(all_cc[i],&subset_size);CHKERRQ(ierr);
4327821e9e7SStefano Zampini     max_subset_size = PetscMax(max_subset_size,subset_size);
43334a97f8cSStefano Zampini   }
4341cef56d8SStefano Zampini   ierr = PetscMalloc5(max_subset_size,&auxmapping,
43534a97f8cSStefano Zampini                       graph->ncc,&auxlocal_sequential,
43634a97f8cSStefano Zampini                       graph->ncc,&auxlocal_parallel,
43734a97f8cSStefano Zampini                       graph->ncc,&index_sequential,
43834a97f8cSStefano Zampini                       graph->ncc,&index_parallel);CHKERRQ(ierr);
43934a97f8cSStefano Zampini 
4401cef56d8SStefano Zampini   /* if threshold is negative, uses all sequential problems */
4411cef56d8SStefano Zampini   if (pcbddc->deluxe_threshold < 0) pcbddc->deluxe_threshold = max_subset_size;
4421cef56d8SStefano Zampini 
4437821e9e7SStefano Zampini   /* workspace */
4447821e9e7SStefano Zampini   ierr = PetscBTCreate(pcis->n,&bitmask);CHKERRQ(ierr);
4457821e9e7SStefano Zampini   ierr = ISGetIndices(pcis->is_I_local,(const PetscInt**)&idxs);CHKERRQ(ierr);
4467821e9e7SStefano Zampini   for (i=0;i<pcis->n-pcis->n_B;i++) {
4477821e9e7SStefano Zampini     ierr = PetscBTSet(bitmask,idxs[i]);CHKERRQ(ierr);
4487821e9e7SStefano Zampini   }
4497821e9e7SStefano Zampini   ierr = ISRestoreIndices(pcis->is_I_local,(const PetscInt**)&idxs);CHKERRQ(ierr);
4507821e9e7SStefano Zampini 
4517821e9e7SStefano Zampini   /* determine which problem has to be solved in parallel or sequentially */
45234a97f8cSStefano Zampini   n_local_sequential_problems = 0;
45334a97f8cSStefano Zampini   n_local_parallel_problems = 0;
4547821e9e7SStefano Zampini   for (i=0;i<n_all_cc;i++) {
4557821e9e7SStefano Zampini     PetscInt subset_size,j,min_loc = 0;
4567821e9e7SStefano Zampini 
4577821e9e7SStefano Zampini     ierr = ISGetLocalSize(all_cc[i],&subset_size);CHKERRQ(ierr);
4587821e9e7SStefano Zampini     ierr = ISGetIndices(all_cc[i],(const PetscInt**)&idxs);CHKERRQ(ierr);
4597821e9e7SStefano Zampini     for (j=0;j<subset_size;j++) {
4607821e9e7SStefano Zampini       ierr = PetscBTSet(bitmask,idxs[j]);CHKERRQ(ierr);
4617821e9e7SStefano Zampini     }
4627821e9e7SStefano Zampini     ierr = ISLocalToGlobalMappingApply(graph->l2gmap,subset_size,idxs,auxmapping);CHKERRQ(ierr);
46334a97f8cSStefano Zampini     for (j=1;j<subset_size;j++) {
46434a97f8cSStefano Zampini       if (auxmapping[j]<auxmapping[min_loc]) {
46534a97f8cSStefano Zampini         min_loc = j;
46634a97f8cSStefano Zampini       }
46734a97f8cSStefano Zampini     }
46834a97f8cSStefano Zampini     if (subset_size > pcbddc->deluxe_threshold) {
46934a97f8cSStefano Zampini       index_parallel[n_local_parallel_problems] = i;
4707821e9e7SStefano Zampini       auxlocal_parallel[n_local_parallel_problems] = idxs[min_loc];
47134a97f8cSStefano Zampini       n_local_parallel_problems++;
47234a97f8cSStefano Zampini     } else {
47334a97f8cSStefano Zampini       index_sequential[n_local_sequential_problems] = i;
4747821e9e7SStefano Zampini       auxlocal_sequential[n_local_sequential_problems] = idxs[min_loc];
47534a97f8cSStefano Zampini       n_local_sequential_problems++;
47634a97f8cSStefano Zampini     }
4777821e9e7SStefano Zampini     ierr = ISRestoreIndices(all_cc[i],(const PetscInt**)&idxs);CHKERRQ(ierr);
47834a97f8cSStefano Zampini   }
47934a97f8cSStefano Zampini 
4807821e9e7SStefano Zampini   /* diagonal scaling on interface dofs not contained in cc */
4817821e9e7SStefano Zampini   deluxe_ctx->n_simple = 0;
4827821e9e7SStefano Zampini   for (i=0;i<pcis->n;i++) {
4837821e9e7SStefano Zampini     if (!PetscBTLookup(bitmask,i)) {
48434a97f8cSStefano Zampini       deluxe_ctx->n_simple++;
48534a97f8cSStefano Zampini     }
48634a97f8cSStefano Zampini   }
48734a97f8cSStefano Zampini   ierr = PetscMalloc1(deluxe_ctx->n_simple,&deluxe_ctx->idx_simple_B);CHKERRQ(ierr);
48834a97f8cSStefano Zampini   deluxe_ctx->n_simple = 0;
4897821e9e7SStefano Zampini   for (i=0;i<pcis->n;i++) {
4907821e9e7SStefano Zampini     if (!PetscBTLookup(bitmask,i)) {
49134a97f8cSStefano Zampini       deluxe_ctx->idx_simple_B[deluxe_ctx->n_simple++] = i;
49234a97f8cSStefano Zampini     }
49334a97f8cSStefano Zampini   }
49434a97f8cSStefano Zampini   ierr = ISGlobalToLocalMappingApply(pcbddc->BtoNmap,IS_GTOLM_DROP,deluxe_ctx->n_simple,deluxe_ctx->idx_simple_B,&i,deluxe_ctx->idx_simple_B);CHKERRQ(ierr);
49534a97f8cSStefano Zampini   if (i != deluxe_ctx->n_simple) {
49634a97f8cSStefano Zampini     SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Error when mapping simple scaling dofs! %d != %d",i,deluxe_ctx->n_simple);
49734a97f8cSStefano Zampini   }
4987821e9e7SStefano Zampini   ierr = PetscBTDestroy(&bitmask);CHKERRQ(ierr);
49934a97f8cSStefano Zampini 
50034a97f8cSStefano Zampini   /* SetUp local schur complements on subsets TODO better reuse procedure */
50134a97f8cSStefano Zampini   if (!sub_schurs->n_subs) {
50234a97f8cSStefano Zampini     Mat       S_j;
50334a97f8cSStefano Zampini     PetscBool free_used_adj;
50434a97f8cSStefano Zampini     PetscInt  *used_xadj,*used_adjncy;
50534a97f8cSStefano Zampini 
50634a97f8cSStefano Zampini     /* decide the adjacency to be used for determining internal problems for local schur on subsets */
50734a97f8cSStefano Zampini     free_used_adj = PETSC_FALSE;
50834a97f8cSStefano Zampini     if (pcbddc->deluxe_layers == -1) {
50934a97f8cSStefano Zampini       used_xadj = NULL;
51034a97f8cSStefano Zampini       used_adjncy = NULL;
51134a97f8cSStefano Zampini     } else {
51234a97f8cSStefano Zampini       if ((pcbddc->deluxe_use_useradj && pcbddc->mat_graph->xadj) || !pcbddc->deluxe_compute_rowadj) {
51334a97f8cSStefano Zampini         used_xadj = pcbddc->mat_graph->xadj;
51434a97f8cSStefano Zampini         used_adjncy = pcbddc->mat_graph->adjncy;
51534a97f8cSStefano Zampini       } else {
51634a97f8cSStefano Zampini         Mat            mat_adj;
51734a97f8cSStefano Zampini         PetscBool      flg_row=PETSC_TRUE;
51834a97f8cSStefano Zampini         const PetscInt *xadj,*adjncy;
51934a97f8cSStefano Zampini         PetscInt       nvtxs;
52034a97f8cSStefano Zampini 
52134a97f8cSStefano Zampini         ierr = MatConvert(pcbddc->local_mat,MATMPIADJ,MAT_INITIAL_MATRIX,&mat_adj);CHKERRQ(ierr);
52234a97f8cSStefano Zampini         ierr = MatGetRowIJ(mat_adj,0,PETSC_TRUE,PETSC_FALSE,&nvtxs,&xadj,&adjncy,&flg_row);CHKERRQ(ierr);
52334a97f8cSStefano Zampini         if (!flg_row) {
52434a97f8cSStefano Zampini           SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Error in MatGetRowIJ called in %s\n",__FUNCT__);
52534a97f8cSStefano Zampini         }
52634a97f8cSStefano Zampini         ierr = PetscMalloc2(nvtxs+1,&used_xadj,xadj[nvtxs],&used_adjncy);CHKERRQ(ierr);
52734a97f8cSStefano Zampini         ierr = PetscMemcpy(used_xadj,xadj,(nvtxs+1)*sizeof(*xadj));CHKERRQ(ierr);
52834a97f8cSStefano Zampini         ierr = PetscMemcpy(used_adjncy,adjncy,(xadj[nvtxs])*sizeof(*adjncy));CHKERRQ(ierr);
52934a97f8cSStefano Zampini         ierr = MatRestoreRowIJ(mat_adj,0,PETSC_TRUE,PETSC_FALSE,&nvtxs,&xadj,&adjncy,&flg_row);CHKERRQ(ierr);
53034a97f8cSStefano Zampini         if (!flg_row) {
53134a97f8cSStefano Zampini           SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Error in MatRestoreRowIJ called in %s\n",__FUNCT__);
53234a97f8cSStefano Zampini         }
53334a97f8cSStefano Zampini         ierr = MatDestroy(&mat_adj);CHKERRQ(ierr);
53434a97f8cSStefano Zampini         free_used_adj = PETSC_TRUE;
53534a97f8cSStefano Zampini       }
53634a97f8cSStefano Zampini     }
53734a97f8cSStefano Zampini 
53834a97f8cSStefano Zampini     /* Create Schur complement matrix */
53934a97f8cSStefano Zampini     ierr = MatCreateSchurComplement(pcis->A_II,pcis->A_II,pcis->A_IB,pcis->A_BI,pcis->A_BB,&S_j);CHKERRQ(ierr);
54034a97f8cSStefano Zampini     ierr = MatSchurComplementSetKSP(S_j,pcbddc->ksp_D);CHKERRQ(ierr);
5417821e9e7SStefano Zampini 
54234a97f8cSStefano Zampini     /* setup Schur complements on subsets */
5437821e9e7SStefano Zampini     ierr = PCBDDCSubSchursSetUp(sub_schurs,S_j,pcis->is_I_local,pcis->is_B_local,n_all_cc,all_cc,used_xadj,used_adjncy,pcbddc->deluxe_layers);CHKERRQ(ierr);
54434a97f8cSStefano Zampini     ierr = MatDestroy(&S_j);CHKERRQ(ierr);
54534a97f8cSStefano Zampini     /* free adjacency */
54634a97f8cSStefano Zampini     if (free_used_adj) {
54734a97f8cSStefano Zampini       ierr = PetscFree2(used_xadj,used_adjncy);CHKERRQ(ierr);
54834a97f8cSStefano Zampini     }
54934a97f8cSStefano Zampini   }
5507821e9e7SStefano Zampini   for (i=0;i<n_all_cc;i++) {
5517821e9e7SStefano Zampini     ierr = ISDestroy(&all_cc[i]);CHKERRQ(ierr);
5527821e9e7SStefano Zampini   }
5537821e9e7SStefano Zampini   ierr = PetscFree(all_cc);CHKERRQ(ierr);
55434a97f8cSStefano Zampini 
55534a97f8cSStefano Zampini   /* Number parallel problems */
55634a97f8cSStefano Zampini   auxglobal_parallel = 0;
55734a97f8cSStefano Zampini   ierr = PCBDDCSubsetNumbering(PetscObjectComm((PetscObject)pc),graph->l2gmap,n_local_parallel_problems,auxlocal_parallel,PETSC_NULL,&n_parallel_problems,&auxglobal_parallel);CHKERRQ(ierr);
55834a97f8cSStefano Zampini   if (pcbddc->dbg_flag) {
55934a97f8cSStefano Zampini     ierr = PetscViewerASCIIPrintf(pcbddc->dbg_viewer,"Deluxe global number of parallel subproblems: %d\n",n_parallel_problems);
56034a97f8cSStefano Zampini   }
56134a97f8cSStefano Zampini 
56234a97f8cSStefano Zampini   /* Compute data structures to solve parallel problems */
56334a97f8cSStefano Zampini   ierr = PCBDDCScalingSetUp_Deluxe_Par(pc,n_local_parallel_problems,n_parallel_problems,auxglobal_parallel,index_parallel);CHKERRQ(ierr);
56434a97f8cSStefano Zampini   ierr = PetscFree(auxglobal_parallel);CHKERRQ(ierr);
56534a97f8cSStefano Zampini 
56634a97f8cSStefano Zampini 
56734a97f8cSStefano Zampini   /* Number sequential problems */
56834a97f8cSStefano Zampini   auxglobal_sequential = 0;
56934a97f8cSStefano Zampini   ierr = PCBDDCSubsetNumbering(PetscObjectComm((PetscObject)pc),graph->l2gmap,n_local_sequential_problems,auxlocal_sequential,PETSC_NULL,&n_sequential_problems,&auxglobal_sequential);CHKERRQ(ierr);
57034a97f8cSStefano Zampini   if (pcbddc->dbg_flag) {
57134a97f8cSStefano Zampini     ierr = PetscViewerASCIIPrintf(pcbddc->dbg_viewer,"Deluxe global number of sequential subproblems: %d\n",n_sequential_problems);
57234a97f8cSStefano Zampini   }
57334a97f8cSStefano Zampini 
57434a97f8cSStefano Zampini   /* Compute data structures to solve sequential problems */
57534a97f8cSStefano Zampini   ierr = PCBDDCScalingSetUp_Deluxe_Seq(pc,n_local_sequential_problems,n_sequential_problems,auxglobal_sequential,index_sequential);CHKERRQ(ierr);
57634a97f8cSStefano Zampini   ierr = PetscFree(auxglobal_sequential);CHKERRQ(ierr);
57734a97f8cSStefano Zampini 
57834a97f8cSStefano Zampini   /* free workspace */
57934a97f8cSStefano Zampini   ierr = PetscFree5(auxmapping,auxlocal_sequential,auxlocal_parallel,index_sequential,index_parallel);CHKERRQ(ierr);
58034a97f8cSStefano Zampini 
58134a97f8cSStefano Zampini   /* free graph struct */
58234a97f8cSStefano Zampini   if (pcbddc->deluxe_rebuild) {
58334a97f8cSStefano Zampini     ierr = PCBDDCGraphDestroy(&graph);CHKERRQ(ierr);
58434a97f8cSStefano Zampini   }
58534a97f8cSStefano Zampini   PetscFunctionReturn(0);
58634a97f8cSStefano Zampini }
58734a97f8cSStefano Zampini 
58834a97f8cSStefano Zampini #undef __FUNCT__
58934a97f8cSStefano Zampini #define __FUNCT__ "PCBDDCScalingSetUp_Deluxe_Par"
59034a97f8cSStefano Zampini static PetscErrorCode PCBDDCScalingSetUp_Deluxe_Par(PC pc, PetscInt n_local_parallel_problems,PetscInt n_parallel_problems,PetscInt global_parallel[],PetscInt index_parallel[])
59134a97f8cSStefano Zampini {
59234a97f8cSStefano Zampini   PC_BDDC                *pcbddc=(PC_BDDC*)pc->data;
59334a97f8cSStefano Zampini   PCBDDCDeluxeScaling    deluxe_ctx=pcbddc->deluxe_ctx;
59434a97f8cSStefano Zampini   /* coloring */
59534a97f8cSStefano Zampini   Mat                    parallel_problems;
59634a97f8cSStefano Zampini   MatColoring            coloring_obj;
59734a97f8cSStefano Zampini   ISColoring             coloring_parallel_problems;
59834a97f8cSStefano Zampini   IS                     *par_is_colors,*is_colors;
59934a97f8cSStefano Zampini   /* working stuff */
60034a97f8cSStefano Zampini   PetscInt               i,j;
60134a97f8cSStefano Zampini   PetscErrorCode         ierr;
60234a97f8cSStefano Zampini 
60334a97f8cSStefano Zampini   PetscFunctionBegin;
60434a97f8cSStefano Zampini   if (!n_parallel_problems) {
60534a97f8cSStefano Zampini     PetscFunctionReturn(0);
60634a97f8cSStefano Zampini   }
60734a97f8cSStefano Zampini   /* Color parallel subproblems */
60834a97f8cSStefano Zampini   ierr = MatCreate(PetscObjectComm((PetscObject)pc),&parallel_problems);CHKERRQ(ierr);
60934a97f8cSStefano Zampini   ierr = MatSetSizes(parallel_problems,PETSC_DECIDE,PETSC_DECIDE,n_parallel_problems,n_parallel_problems);CHKERRQ(ierr);
61034a97f8cSStefano Zampini   ierr = MatSetType(parallel_problems,MATAIJ);CHKERRQ(ierr);
61134a97f8cSStefano Zampini   ierr = MatSetUp(parallel_problems);CHKERRQ(ierr);
61234a97f8cSStefano Zampini   ierr = MatSetOption(parallel_problems,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_FALSE);CHKERRQ(ierr);
61334a97f8cSStefano Zampini   ierr = MatSetOption(parallel_problems,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_FALSE);CHKERRQ(ierr);
61434a97f8cSStefano Zampini   for (i=0;i<n_local_parallel_problems;i++) {
61534a97f8cSStefano Zampini     PetscInt row = global_parallel[i];
61634a97f8cSStefano Zampini     for (j=0;j<n_local_parallel_problems;j++) {
61734a97f8cSStefano Zampini       PetscInt col = global_parallel[j];
61834a97f8cSStefano Zampini       if (row != col) {
61934a97f8cSStefano Zampini         ierr = MatSetValue(parallel_problems,row,col,1.0,INSERT_VALUES);CHKERRQ(ierr);
62034a97f8cSStefano Zampini       }
62134a97f8cSStefano Zampini     }
62234a97f8cSStefano Zampini   }
62334a97f8cSStefano Zampini   ierr = MatAssemblyBegin(parallel_problems,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
62434a97f8cSStefano Zampini   ierr = MatAssemblyEnd(parallel_problems,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
62534a97f8cSStefano Zampini   if (pcbddc->dbg_flag > 1) {
62634a97f8cSStefano Zampini     ierr = PetscViewerFlush(pcbddc->dbg_viewer);CHKERRQ(ierr);
62734a97f8cSStefano Zampini     ierr = PetscViewerASCIIPrintf(pcbddc->dbg_viewer,"Adj matrix for deluxe parallel problems\n");CHKERRQ(ierr);
62834a97f8cSStefano Zampini     ierr = MatView(parallel_problems,pcbddc->dbg_viewer);CHKERRQ(ierr);
62934a97f8cSStefano Zampini   }
63034a97f8cSStefano Zampini   ierr = MatColoringCreate(parallel_problems,&coloring_obj);CHKERRQ(ierr);
63134a97f8cSStefano Zampini   ierr = MatColoringSetDistance(coloring_obj,1);CHKERRQ(ierr);
63234a97f8cSStefano Zampini   ierr = MatColoringSetType(coloring_obj,MATCOLORINGJP);CHKERRQ(ierr);
63334a97f8cSStefano Zampini   ierr = MatColoringApply(coloring_obj,&coloring_parallel_problems);CHKERRQ(ierr);
63434a97f8cSStefano Zampini   ierr = ISColoringGetIS(coloring_parallel_problems,&deluxe_ctx->par_colors,&par_is_colors);CHKERRQ(ierr);
63534a97f8cSStefano Zampini   if (pcbddc->dbg_flag) {
63634a97f8cSStefano Zampini     ierr = PetscViewerFlush(pcbddc->dbg_viewer);CHKERRQ(ierr);
63734a97f8cSStefano Zampini     ierr = PetscViewerASCIIPrintf(pcbddc->dbg_viewer,"Number of colors %d for parallel part of deluxe\n",deluxe_ctx->par_colors);CHKERRQ(ierr);
63834a97f8cSStefano Zampini   }
63934a97f8cSStefano Zampini 
64034a97f8cSStefano Zampini   /* all procs should know the color distribution */
64134a97f8cSStefano Zampini   ierr = PetscMalloc1(deluxe_ctx->par_colors,&is_colors);CHKERRQ(ierr);
64234a97f8cSStefano Zampini   for (i=0;i<deluxe_ctx->par_colors;i++) {
64334a97f8cSStefano Zampini     if (pcbddc->dbg_flag) {
64434a97f8cSStefano Zampini       ierr = PetscViewerASCIIPrintf(pcbddc->dbg_viewer,"Global problem indexes for color %d\n",i);CHKERRQ(ierr);
64534a97f8cSStefano Zampini       ierr = ISView(par_is_colors[i],pcbddc->dbg_viewer);CHKERRQ(ierr);
64634a97f8cSStefano Zampini       ierr = PetscViewerFlush(pcbddc->dbg_viewer);CHKERRQ(ierr);
64734a97f8cSStefano Zampini     }
64834a97f8cSStefano Zampini     ierr = ISAllGather(par_is_colors[i],&is_colors[i]);CHKERRQ(ierr);
64934a97f8cSStefano Zampini   }
65034a97f8cSStefano Zampini 
65134a97f8cSStefano Zampini   /* free unneeded objects */
65234a97f8cSStefano Zampini   ierr = ISColoringRestoreIS(coloring_parallel_problems,&par_is_colors);CHKERRQ(ierr);
65334a97f8cSStefano Zampini   ierr = ISColoringDestroy(&coloring_parallel_problems);CHKERRQ(ierr);
65434a97f8cSStefano Zampini   ierr = MatColoringDestroy(&coloring_obj);CHKERRQ(ierr);
65534a97f8cSStefano Zampini   ierr = MatDestroy(&parallel_problems);CHKERRQ(ierr);
65634a97f8cSStefano Zampini 
65734a97f8cSStefano Zampini   /* allocate deluxe arrays for parallel problems */
65834a97f8cSStefano Zampini   ierr = PetscMalloc5(deluxe_ctx->par_colors,&deluxe_ctx->par_ksp,
65934a97f8cSStefano Zampini                       deluxe_ctx->par_colors,&deluxe_ctx->par_scctx_s,
66034a97f8cSStefano Zampini                       deluxe_ctx->par_colors,&deluxe_ctx->par_scctx_p,
66134a97f8cSStefano Zampini                       deluxe_ctx->par_colors,&deluxe_ctx->par_vec,
66234a97f8cSStefano Zampini                       deluxe_ctx->par_colors,&deluxe_ctx->par_col2sub);CHKERRQ(ierr);
66334a97f8cSStefano Zampini 
66434a97f8cSStefano Zampini   /* cycle on colors */
66534a97f8cSStefano Zampini   for (i=0;i<deluxe_ctx->par_colors;i++) {
66634a97f8cSStefano Zampini     PetscSubcomm    par_subcomm;
66734a97f8cSStefano Zampini     const PetscInt* idxs_subproblems;
66834a97f8cSStefano Zampini     PetscInt        color_size;
66934a97f8cSStefano Zampini     PetscMPIInt     rank,active_color;
67034a97f8cSStefano Zampini 
67134a97f8cSStefano Zampini     /* get local index of i-th parallel colored problem */
67234a97f8cSStefano Zampini     ierr = ISGetLocalSize(is_colors[i],&color_size);CHKERRQ(ierr);
67334a97f8cSStefano Zampini     ierr = ISGetIndices(is_colors[i],&idxs_subproblems);CHKERRQ(ierr);
67434a97f8cSStefano Zampini     /* split comm for computing parallel problems for this color */
67534a97f8cSStefano Zampini     /* Processes not partecipating at this stage will have color = color_size */
67634a97f8cSStefano Zampini     /* because PetscCommDuplicate does not handle MPI_COMM_NULL */
67734a97f8cSStefano Zampini     active_color = color_size;
67834a97f8cSStefano Zampini     deluxe_ctx->par_col2sub[i] = -1;
67934a97f8cSStefano Zampini     for (j=0;j<n_local_parallel_problems;j++) {
68034a97f8cSStefano Zampini       PetscInt local_idx;
68134a97f8cSStefano Zampini       ierr = PetscFindInt(global_parallel[j],color_size,idxs_subproblems,&local_idx);CHKERRQ(ierr);
68234a97f8cSStefano Zampini       if (local_idx > -1) {
68334a97f8cSStefano Zampini         ierr = PetscMPIIntCast(local_idx,&active_color);CHKERRQ(ierr);
68434a97f8cSStefano Zampini         deluxe_ctx->par_col2sub[i] = index_parallel[j];
68534a97f8cSStefano Zampini         break;
68634a97f8cSStefano Zampini       }
68734a97f8cSStefano Zampini     }
68834a97f8cSStefano Zampini     ierr = ISRestoreIndices(is_colors[i],&idxs_subproblems);CHKERRQ(ierr);
68934a97f8cSStefano Zampini     ierr = PetscSubcommCreate(PetscObjectComm((PetscObject)pc),&par_subcomm);CHKERRQ(ierr);
69034a97f8cSStefano Zampini     ierr = PetscSubcommSetNumber(par_subcomm,color_size+1);CHKERRQ(ierr);
69134a97f8cSStefano Zampini     ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)pc),&rank);CHKERRQ(ierr);
69234a97f8cSStefano Zampini     ierr = PetscSubcommSetTypeGeneral(par_subcomm,active_color,rank);CHKERRQ(ierr);
69334a97f8cSStefano Zampini     /* print debug info */
69434a97f8cSStefano Zampini     if (pcbddc->dbg_flag) {
69534a97f8cSStefano Zampini       PetscMPIInt crank,csize;
69634a97f8cSStefano Zampini       ierr = MPI_Comm_rank(par_subcomm->comm,&crank);CHKERRQ(ierr);
69734a97f8cSStefano Zampini       ierr = MPI_Comm_size(par_subcomm->comm,&csize);CHKERRQ(ierr);
69834a97f8cSStefano Zampini       ierr = PetscViewerASCIIPrintf(pcbddc->dbg_viewer,"Color %d: size %d, details follows.\n",i,color_size);CHKERRQ(ierr);
69934a97f8cSStefano Zampini       ierr = PetscViewerFlush(pcbddc->dbg_viewer);CHKERRQ(ierr);
70034a97f8cSStefano Zampini       ierr = PetscViewerASCIISynchronizedAllow(pcbddc->dbg_viewer,PETSC_TRUE);CHKERRQ(ierr);
70134a97f8cSStefano 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);
70234a97f8cSStefano Zampini       ierr = PetscViewerFlush(pcbddc->dbg_viewer);CHKERRQ(ierr);
70334a97f8cSStefano Zampini     }
70434a97f8cSStefano Zampini 
70534a97f8cSStefano Zampini     if (deluxe_ctx->par_col2sub[i] >= 0) {
70634a97f8cSStefano Zampini       PC                     pc;
70734a97f8cSStefano Zampini       Mat                    color_mat,color_mat_is,temp_mat;
70834a97f8cSStefano Zampini       ISLocalToGlobalMapping WtoNmap,l2gmap_subset;
70934a97f8cSStefano Zampini       IS                     is_local_numbering,isB_local,isW_local,isW;
71034a97f8cSStefano Zampini       PCBDDCSubSchurs        sub_schurs = deluxe_ctx->sub_schurs;
71134a97f8cSStefano Zampini       PetscInt               subidx,n_local_dofs,n_global_dofs;
71234a97f8cSStefano Zampini       PetscInt               *global_numbering,*local_numbering;
71334a97f8cSStefano Zampini       char                   ksp_prefix[256];
71434a97f8cSStefano Zampini       size_t                 len;
71534a97f8cSStefano Zampini 
71634a97f8cSStefano Zampini       /* Local index for schur complement on subset */
71734a97f8cSStefano Zampini       subidx = deluxe_ctx->par_col2sub[i];
71834a97f8cSStefano Zampini 
71934a97f8cSStefano Zampini       /* Parallel numbering for dofs in colored subset */
72034a97f8cSStefano Zampini       ierr = ISSum(sub_schurs->is_AEj_I[subidx],sub_schurs->is_AEj_B[subidx],&is_local_numbering);CHKERRQ(ierr);
72134a97f8cSStefano Zampini       ierr = ISGetLocalSize(is_local_numbering,&n_local_dofs);CHKERRQ(ierr);
72234a97f8cSStefano Zampini       ierr = ISGetIndices(is_local_numbering,(const PetscInt **)&local_numbering);CHKERRQ(ierr);
72334a97f8cSStefano Zampini       ierr = PCBDDCSubsetNumbering(par_subcomm->comm,pcbddc->mat_graph->l2gmap,n_local_dofs,local_numbering,PETSC_NULL,&n_global_dofs,&global_numbering);CHKERRQ(ierr);
72434a97f8cSStefano Zampini       ierr = ISRestoreIndices(is_local_numbering,(const PetscInt **)&local_numbering);CHKERRQ(ierr);
72534a97f8cSStefano Zampini 
72634a97f8cSStefano Zampini       /* L2Gmap from relevant dofs to local dofs */
72734a97f8cSStefano Zampini       ierr = ISLocalToGlobalMappingCreateIS(is_local_numbering,&WtoNmap);CHKERRQ(ierr);
72834a97f8cSStefano Zampini 
72934a97f8cSStefano Zampini       /* L2Gmap from local to global dofs */
73034a97f8cSStefano Zampini       ierr = ISLocalToGlobalMappingCreate(par_subcomm->comm,1,n_local_dofs,global_numbering,PETSC_COPY_VALUES,&l2gmap_subset);CHKERRQ(ierr);
73134a97f8cSStefano Zampini 
73234a97f8cSStefano Zampini       /* compute parallel matrix (extended dirichlet problem on subset) */
73334a97f8cSStefano Zampini       ierr = MatCreateIS(par_subcomm->comm,1,PETSC_DECIDE,PETSC_DECIDE,n_global_dofs,n_global_dofs,l2gmap_subset,&color_mat_is);CHKERRQ(ierr);
73434a97f8cSStefano Zampini       ierr = MatGetSubMatrix(pcbddc->local_mat,is_local_numbering,is_local_numbering,MAT_INITIAL_MATRIX,&temp_mat);CHKERRQ(ierr);
73534a97f8cSStefano Zampini       ierr = MatISSetLocalMat(color_mat_is,temp_mat);CHKERRQ(ierr);
73634a97f8cSStefano Zampini       ierr = MatDestroy(&temp_mat);CHKERRQ(ierr);
73734a97f8cSStefano Zampini       ierr = MatISGetMPIXAIJ(color_mat_is,MAT_INITIAL_MATRIX,&color_mat);CHKERRQ(ierr);
73834a97f8cSStefano Zampini       ierr = MatDestroy(&color_mat_is);CHKERRQ(ierr);
73934a97f8cSStefano Zampini 
74034a97f8cSStefano Zampini       /* work vector for (parallel) extended dirichlet problem */
7418a26ef87SStefano Zampini       ierr = MatCreateVecs(color_mat,&deluxe_ctx->par_vec[i],NULL);CHKERRQ(ierr);
74234a97f8cSStefano Zampini 
74334a97f8cSStefano Zampini       /* compute scatters */
74434a97f8cSStefano Zampini       /* deluxe_ctx->par_scctx_p[i] extension from local subset to extended dirichlet problem
74534a97f8cSStefano Zampini          deluxe_ctx->par_scctx_s[i] restriction from local boundary to subset -> simple copy of selected values */
74634a97f8cSStefano Zampini       ierr = ISGlobalToLocalMappingApplyIS(pcbddc->BtoNmap,IS_GTOLM_DROP,sub_schurs->is_AEj_B[subidx],&isB_local);CHKERRQ(ierr);
74734a97f8cSStefano Zampini       ierr = VecScatterCreate(pcbddc->work_scaling,isB_local,sub_schurs->work1[subidx],NULL,&deluxe_ctx->par_scctx_s[i]);CHKERRQ(ierr);
74834a97f8cSStefano Zampini       ierr = ISGlobalToLocalMappingApplyIS(WtoNmap,IS_GTOLM_DROP,sub_schurs->is_AEj_B[subidx],&isW_local);CHKERRQ(ierr);
74934a97f8cSStefano Zampini       ierr = ISLocalToGlobalMappingApplyIS(l2gmap_subset,isW_local,&isW);CHKERRQ(ierr);
75034a97f8cSStefano Zampini       ierr = VecScatterCreate(sub_schurs->work1[subidx],NULL,deluxe_ctx->par_vec[i],isW,&deluxe_ctx->par_scctx_p[i]);CHKERRQ(ierr);
75134a97f8cSStefano Zampini 
75234a97f8cSStefano Zampini       /* free objects no longer neeeded */
75334a97f8cSStefano Zampini       ierr = ISDestroy(&isW);CHKERRQ(ierr);
75434a97f8cSStefano Zampini       ierr = ISDestroy(&isW_local);CHKERRQ(ierr);
75534a97f8cSStefano Zampini       ierr = ISDestroy(&isB_local);CHKERRQ(ierr);
75634a97f8cSStefano Zampini       ierr = ISLocalToGlobalMappingDestroy(&WtoNmap);CHKERRQ(ierr);
75734a97f8cSStefano Zampini       ierr = ISLocalToGlobalMappingDestroy(&l2gmap_subset);CHKERRQ(ierr);
75834a97f8cSStefano Zampini       ierr = ISDestroy(&is_local_numbering);CHKERRQ(ierr);
75934a97f8cSStefano Zampini       ierr = PetscFree(global_numbering);CHKERRQ(ierr);
76034a97f8cSStefano Zampini 
76134a97f8cSStefano Zampini       /* KSP for extended dirichlet problem */
76234a97f8cSStefano Zampini       ierr = KSPCreate(par_subcomm->comm,&deluxe_ctx->par_ksp[i]);CHKERRQ(ierr);
76334a97f8cSStefano Zampini       ierr = KSPSetOperators(deluxe_ctx->par_ksp[i],color_mat,color_mat);CHKERRQ(ierr);
76434a97f8cSStefano Zampini       ierr = KSPSetTolerances(deluxe_ctx->par_ksp[i],1.e-12,1.e-12,1.e10,10000);CHKERRQ(ierr);
76534a97f8cSStefano Zampini       ierr = KSPSetType(deluxe_ctx->par_ksp[i],KSPPREONLY);CHKERRQ(ierr);
76634a97f8cSStefano Zampini       ierr = KSPGetPC(deluxe_ctx->par_ksp[i],&pc);CHKERRQ(ierr);
76734a97f8cSStefano Zampini       ierr = PCSetType(pc,PCREDUNDANT);CHKERRQ(ierr);
76834a97f8cSStefano Zampini       ierr = PetscStrlen(((PetscObject)(pcbddc->ksp_D))->prefix,&len);CHKERRQ(ierr);
7698856534dSStefano Zampini       len -= 10; /* remove "dirichlet_" */
7708856534dSStefano Zampini       ierr = PetscStrncpy(ksp_prefix,((PetscObject)(pcbddc->ksp_D))->prefix,len+1);CHKERRQ(ierr); /* PetscStrncpy puts a terminating char at the end */
77134a97f8cSStefano Zampini       ierr = PetscStrcat(ksp_prefix,"deluxe_par_");CHKERRQ(ierr);
77234a97f8cSStefano Zampini       ierr = KSPSetOptionsPrefix(deluxe_ctx->par_ksp[i],ksp_prefix);CHKERRQ(ierr);
77334a97f8cSStefano Zampini       ierr = KSPSetFromOptions(deluxe_ctx->par_ksp[i]);CHKERRQ(ierr);
77434a97f8cSStefano Zampini       ierr = KSPSetUp(deluxe_ctx->par_ksp[i]);CHKERRQ(ierr);
77534a97f8cSStefano Zampini       ierr = MatDestroy(&color_mat);CHKERRQ(ierr);
77634a97f8cSStefano Zampini     } else { /* not partecipating in color */
77734a97f8cSStefano Zampini       deluxe_ctx->par_ksp[i] = 0;
77834a97f8cSStefano Zampini       deluxe_ctx->par_vec[i] = 0;
77934a97f8cSStefano Zampini       deluxe_ctx->par_scctx_p[i] = 0;
78034a97f8cSStefano Zampini       deluxe_ctx->par_scctx_s[i] = 0;
78134a97f8cSStefano Zampini     }
78234a97f8cSStefano Zampini     ierr = PetscSubcommDestroy(&par_subcomm);CHKERRQ(ierr);
78334a97f8cSStefano Zampini   }
78434a97f8cSStefano Zampini   for (i=0;i<deluxe_ctx->par_colors;i++) {
78534a97f8cSStefano Zampini     ierr = ISDestroy(&is_colors[i]);CHKERRQ(ierr);
78634a97f8cSStefano Zampini   }
78734a97f8cSStefano Zampini   ierr = PetscFree(is_colors);CHKERRQ(ierr);
78834a97f8cSStefano Zampini 
78934a97f8cSStefano Zampini   if (pcbddc->dbg_flag) {
79034a97f8cSStefano Zampini     Vec test_vec;
79134a97f8cSStefano Zampini     PetscReal error;
79234a97f8cSStefano Zampini     PCBDDCSubSchurs sub_schurs = deluxe_ctx->sub_schurs;
79334a97f8cSStefano Zampini     /* test partition of unity of coloured schur complements  */
79434a97f8cSStefano Zampini     for (i=0;i<deluxe_ctx->par_colors;i++) {
79534a97f8cSStefano Zampini       PetscInt  subidx = deluxe_ctx->par_col2sub[i];
79634a97f8cSStefano Zampini       PetscBool error_found = PETSC_FALSE;
79734a97f8cSStefano Zampini       ierr = PetscViewerASCIISynchronizedAllow(pcbddc->dbg_viewer,PETSC_TRUE);CHKERRQ(ierr);
79834a97f8cSStefano Zampini 
79934a97f8cSStefano Zampini       if (deluxe_ctx->par_ksp[i]) {
80034a97f8cSStefano Zampini         /* create random test vec being zero on internal nodes of the extende dirichlet problem */
80134a97f8cSStefano Zampini         ierr = VecDuplicate(deluxe_ctx->par_vec[i],&test_vec);CHKERRQ(ierr);
80234a97f8cSStefano Zampini         ierr = VecSetRandom(sub_schurs->work1[subidx],PETSC_NULL);CHKERRQ(ierr);
80334a97f8cSStefano Zampini         ierr = VecSet(test_vec,0.0);CHKERRQ(ierr);
80434a97f8cSStefano Zampini         ierr = VecScatterBegin(deluxe_ctx->par_scctx_p[i],sub_schurs->work1[subidx],test_vec,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
80534a97f8cSStefano Zampini         ierr = VecScatterEnd(deluxe_ctx->par_scctx_p[i],sub_schurs->work1[subidx],test_vec,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
80634a97f8cSStefano Zampini         /* w_j */
80734a97f8cSStefano Zampini         ierr = VecScatterBegin(deluxe_ctx->par_scctx_p[i],test_vec,sub_schurs->work1[subidx],INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
80834a97f8cSStefano Zampini         ierr = VecScatterEnd(deluxe_ctx->par_scctx_p[i],test_vec,sub_schurs->work1[subidx],INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
80934a97f8cSStefano Zampini         /* S_j*w_j */
81034a97f8cSStefano Zampini         ierr = MatMult(sub_schurs->S_Ej[subidx],sub_schurs->work1[subidx],sub_schurs->work2[subidx]);CHKERRQ(ierr);
81134a97f8cSStefano Zampini         /* \sum_j S_j*w_j */
81234a97f8cSStefano Zampini         ierr = VecSet(deluxe_ctx->par_vec[i],0.0);CHKERRQ(ierr);
81334a97f8cSStefano Zampini         ierr = VecScatterBegin(deluxe_ctx->par_scctx_p[i],sub_schurs->work2[subidx],deluxe_ctx->par_vec[i],ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
81434a97f8cSStefano Zampini         ierr = VecScatterEnd(deluxe_ctx->par_scctx_p[i],sub_schurs->work2[subidx],deluxe_ctx->par_vec[i],ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
81534a97f8cSStefano Zampini         /* (\sum_j S_j)^(-1)(\sum_j S_j*w_j) */
81634a97f8cSStefano Zampini         ierr = KSPSolve(deluxe_ctx->par_ksp[i],deluxe_ctx->par_vec[i],deluxe_ctx->par_vec[i]);CHKERRQ(ierr);
81734a97f8cSStefano Zampini         ierr = VecScatterBegin(deluxe_ctx->par_scctx_p[i],deluxe_ctx->par_vec[i],sub_schurs->work1[subidx],INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
81834a97f8cSStefano Zampini         ierr = VecScatterEnd(deluxe_ctx->par_scctx_p[i],deluxe_ctx->par_vec[i],sub_schurs->work1[subidx],INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
81934a97f8cSStefano Zampini         ierr = VecSet(deluxe_ctx->par_vec[i],0.0);CHKERRQ(ierr);
82034a97f8cSStefano Zampini         ierr = VecScatterBegin(deluxe_ctx->par_scctx_p[i],sub_schurs->work1[subidx],deluxe_ctx->par_vec[i],INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
82134a97f8cSStefano Zampini         ierr = VecScatterEnd(deluxe_ctx->par_scctx_p[i],sub_schurs->work1[subidx],deluxe_ctx->par_vec[i],INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
82234a97f8cSStefano Zampini         /* test partition of unity */
82334a97f8cSStefano Zampini         ierr = VecAXPY(test_vec,-1.0,deluxe_ctx->par_vec[i]);CHKERRQ(ierr);
82434a97f8cSStefano Zampini         ierr = VecNorm(test_vec,NORM_INFINITY,&error);CHKERRQ(ierr);
825c63f45b2SStefano Zampini         if (PetscAbsReal(error) > 1.e-2) {
82634a97f8cSStefano Zampini           /* ierr = VecView(test_vec,0);CHKERRQ(ierr); */
82734a97f8cSStefano Zampini           error_found = PETSC_TRUE;
82834a97f8cSStefano Zampini         }
82934a97f8cSStefano Zampini         ierr = VecDestroy(&test_vec);CHKERRQ(ierr);
83034a97f8cSStefano Zampini       }
83134a97f8cSStefano Zampini       if (error_found) {
83234a97f8cSStefano Zampini         ierr = PetscViewerASCIISynchronizedPrintf(pcbddc->dbg_viewer,"Error testing local schur for color %d and subdomain %d\n",i,PetscGlobalRank);CHKERRQ(ierr);
83334a97f8cSStefano Zampini       }
83434a97f8cSStefano Zampini       ierr = PetscViewerFlush(pcbddc->dbg_viewer);CHKERRQ(ierr);
83534a97f8cSStefano Zampini     }
83634a97f8cSStefano Zampini   }
83734a97f8cSStefano Zampini   PetscFunctionReturn(0);
83834a97f8cSStefano Zampini }
83934a97f8cSStefano Zampini 
84034a97f8cSStefano Zampini 
84134a97f8cSStefano Zampini #undef __FUNCT__
84234a97f8cSStefano Zampini #define __FUNCT__ "PCBDDCScalingSetUp_Deluxe_Seq"
84334a97f8cSStefano Zampini static PetscErrorCode PCBDDCScalingSetUp_Deluxe_Seq(PC pc,PetscInt n_local_sequential_problems,PetscInt n_sequential_problems,PetscInt global_sequential[],PetscInt local_sequential[])
84434a97f8cSStefano Zampini {
84534a97f8cSStefano Zampini   PC_BDDC             *pcbddc=(PC_BDDC*)pc->data;
84634a97f8cSStefano Zampini   PCBDDCDeluxeScaling deluxe_ctx=pcbddc->deluxe_ctx;
84734a97f8cSStefano Zampini   PCBDDCSubSchurs     sub_schurs = deluxe_ctx->sub_schurs;
84834a97f8cSStefano Zampini   Mat                 global_schur_subsets,*submat_global_schur_subsets,work_mat;
84934a97f8cSStefano Zampini   IS                  is_to,is_from;
85034a97f8cSStefano Zampini   PetscScalar         *array,*fill_vals;
85134a97f8cSStefano Zampini   PetscInt            *all_local_idx_G,*all_local_idx_B,*all_local_idx_N,*all_permutation_G,*dummy_idx;
85234a97f8cSStefano Zampini   PetscInt            i,j,k,local_problem_index;
85334a97f8cSStefano Zampini   PetscInt            subset_size,max_subset_size,max_subset_size_red;
85434a97f8cSStefano Zampini   PetscInt            local_size,global_size;
855ccdc0be0SStefano Zampini   PC                  pc_temp;
856ccdc0be0SStefano Zampini   MatSolverPackage    solver=NULL;
85734a97f8cSStefano Zampini   char                ksp_prefix[256];
85834a97f8cSStefano Zampini   size_t              len;
85934a97f8cSStefano Zampini   PetscErrorCode      ierr;
86034a97f8cSStefano Zampini 
86134a97f8cSStefano Zampini   PetscFunctionBegin;
86234a97f8cSStefano Zampini   if (!n_sequential_problems) {
86334a97f8cSStefano Zampini     PetscFunctionReturn(0);
86434a97f8cSStefano Zampini   }
86534a97f8cSStefano Zampini   /* Get info on subset sizes and sum of all subsets sizes */
86634a97f8cSStefano Zampini   max_subset_size = 0;
86734a97f8cSStefano Zampini   local_size = 0;
86834a97f8cSStefano Zampini   for (i=0;i<n_local_sequential_problems;i++) {
86934a97f8cSStefano Zampini     local_problem_index = local_sequential[i];
87034a97f8cSStefano Zampini     ierr = ISGetLocalSize(sub_schurs->is_AEj_B[local_problem_index],&subset_size);CHKERRQ(ierr);
87134a97f8cSStefano Zampini     max_subset_size = PetscMax(subset_size,max_subset_size);
87234a97f8cSStefano Zampini     local_size += subset_size;
87334a97f8cSStefano Zampini   }
87434a97f8cSStefano Zampini 
87534a97f8cSStefano Zampini   /* Work arrays for local indices */
87634a97f8cSStefano Zampini   ierr = PetscMalloc1(local_size,&all_local_idx_B);CHKERRQ(ierr);
87734a97f8cSStefano Zampini   ierr = PetscMalloc1(local_size,&all_local_idx_N);CHKERRQ(ierr);
87834a97f8cSStefano Zampini 
87934a97f8cSStefano Zampini   /* Get local indices in local whole numbering and local boundary numbering */
88034a97f8cSStefano Zampini   local_size = 0;
88134a97f8cSStefano Zampini   for (i=0;i<n_local_sequential_problems;i++) {
88234a97f8cSStefano Zampini     PetscInt *idxs;
88334a97f8cSStefano Zampini     /* get info on local problem */
88434a97f8cSStefano Zampini     local_problem_index = local_sequential[i];
88534a97f8cSStefano Zampini     ierr = ISGetLocalSize(sub_schurs->is_AEj_B[local_problem_index],&subset_size);CHKERRQ(ierr);
88634a97f8cSStefano Zampini     ierr = ISGetIndices(sub_schurs->is_AEj_B[local_problem_index],(const PetscInt**)&idxs);CHKERRQ(ierr);
88734a97f8cSStefano Zampini     /* subset indices in local numbering */
88834a97f8cSStefano Zampini     ierr = PetscMemcpy(all_local_idx_N+local_size,idxs,subset_size*sizeof(PetscInt));CHKERRQ(ierr);
88934a97f8cSStefano Zampini     /* subset indices in local boundary numbering */
89034a97f8cSStefano Zampini     ierr = ISGlobalToLocalMappingApply(pcbddc->BtoNmap,IS_GTOLM_DROP,subset_size,idxs,&j,&all_local_idx_B[local_size]);CHKERRQ(ierr);
89134a97f8cSStefano Zampini     ierr = ISRestoreIndices(sub_schurs->is_AEj_B[local_problem_index],(const PetscInt**)&idxs);CHKERRQ(ierr);
89234a97f8cSStefano Zampini     if (j != subset_size) {
89334a97f8cSStefano Zampini       SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Error in BDDC deluxe serial %d (BtoNmap)! %d != %d\n",local_problem_index,subset_size,j);
89434a97f8cSStefano Zampini     }
89534a97f8cSStefano Zampini     local_size += subset_size;
89634a97f8cSStefano Zampini   }
89734a97f8cSStefano Zampini 
89834a97f8cSStefano Zampini   /* Number dofs on all subsets (parallel) and sort numbering */
89934a97f8cSStefano Zampini   ierr = PCBDDCSubsetNumbering(PetscObjectComm((PetscObject)pc),pcbddc->mat_graph->l2gmap,local_size,all_local_idx_N,PETSC_NULL,&global_size,&all_local_idx_G);CHKERRQ(ierr);
90034a97f8cSStefano Zampini   ierr = PetscMalloc1(local_size,&all_permutation_G);CHKERRQ(ierr);
90134a97f8cSStefano Zampini   for (i=0;i<local_size;i++) {
90234a97f8cSStefano Zampini     all_permutation_G[i]=i;
90334a97f8cSStefano Zampini   }
90434a97f8cSStefano Zampini   ierr = PetscSortIntWithPermutation(local_size,all_local_idx_G,all_permutation_G);CHKERRQ(ierr);
90534a97f8cSStefano Zampini 
90634a97f8cSStefano Zampini   /* Local matrix of all local Schur on subsets */
90734a97f8cSStefano Zampini   ierr = MatCreate(PETSC_COMM_SELF,&deluxe_ctx->seq_mat);CHKERRQ(ierr);
90834a97f8cSStefano Zampini   ierr = MatSetSizes(deluxe_ctx->seq_mat,PETSC_DECIDE,PETSC_DECIDE,local_size,local_size);CHKERRQ(ierr);
90934a97f8cSStefano Zampini   ierr = MatSetType(deluxe_ctx->seq_mat,MATAIJ);CHKERRQ(ierr);
91034a97f8cSStefano Zampini   ierr = MatSeqAIJSetPreallocation(deluxe_ctx->seq_mat,max_subset_size,PETSC_NULL);CHKERRQ(ierr);
91134a97f8cSStefano Zampini 
91234a97f8cSStefano Zampini   /* Global matrix of all assembled Schur on subsets */
91334a97f8cSStefano Zampini   ierr = MatCreate(PetscObjectComm((PetscObject)pc),&global_schur_subsets);CHKERRQ(ierr);
91434a97f8cSStefano Zampini   ierr = MatSetSizes(global_schur_subsets,PETSC_DECIDE,PETSC_DECIDE,global_size,global_size);CHKERRQ(ierr);
91534a97f8cSStefano Zampini   ierr = MatSetType(global_schur_subsets,MATAIJ);CHKERRQ(ierr);
91634a97f8cSStefano Zampini   ierr = MPI_Allreduce(&max_subset_size,&max_subset_size_red,1,MPIU_INT,MPI_MAX,PetscObjectComm((PetscObject)pc));CHKERRQ(ierr);
91734a97f8cSStefano Zampini   ierr = MatMPIAIJSetPreallocation(global_schur_subsets,max_subset_size_red,PETSC_NULL,max_subset_size_red,PETSC_NULL);CHKERRQ(ierr);
91834a97f8cSStefano Zampini 
91934a97f8cSStefano Zampini   /* Work arrays */
92034a97f8cSStefano Zampini   ierr = PetscMalloc2(max_subset_size,&dummy_idx,max_subset_size*max_subset_size,&fill_vals);CHKERRQ(ierr);
92134a97f8cSStefano Zampini 
92234a97f8cSStefano Zampini   /* Loop on local problems to compute Schur complements explicitly */
92334a97f8cSStefano Zampini   local_size = 0;
92434a97f8cSStefano Zampini   for (i=0;i<n_local_sequential_problems;i++) {
92534a97f8cSStefano Zampini     /* get info on local problem */
92634a97f8cSStefano Zampini     local_problem_index = local_sequential[i];
92734a97f8cSStefano Zampini     ierr = ISGetLocalSize(sub_schurs->is_AEj_B[local_problem_index],&subset_size);CHKERRQ(ierr);
92834a97f8cSStefano Zampini     /* local Schur */
92934a97f8cSStefano Zampini     for (j=0;j<subset_size;j++) {
93034a97f8cSStefano Zampini       ierr = VecSet(sub_schurs->work1[local_problem_index],0.0);CHKERRQ(ierr);
93134a97f8cSStefano Zampini       ierr = VecSetValue(sub_schurs->work1[local_problem_index],j,1.0,INSERT_VALUES);CHKERRQ(ierr);
93234a97f8cSStefano Zampini       ierr = MatMult(sub_schurs->S_Ej[local_problem_index],sub_schurs->work1[local_problem_index],sub_schurs->work2[local_problem_index]);CHKERRQ(ierr);
93334a97f8cSStefano Zampini       /* store vals */
93434a97f8cSStefano Zampini       ierr = VecGetArray(sub_schurs->work2[local_problem_index],&array);CHKERRQ(ierr);
93534a97f8cSStefano Zampini       for (k=0;k<subset_size;k++) {
93634a97f8cSStefano Zampini         fill_vals[k*subset_size+j] = array[k];
93734a97f8cSStefano Zampini       }
93834a97f8cSStefano Zampini       ierr = VecRestoreArray(sub_schurs->work2[local_problem_index],&array);CHKERRQ(ierr);
93934a97f8cSStefano Zampini     }
94034a97f8cSStefano Zampini     for (j=0;j<subset_size;j++) {
94134a97f8cSStefano Zampini       dummy_idx[j]=local_size+j;
94234a97f8cSStefano Zampini     }
94334a97f8cSStefano Zampini     ierr = MatSetValues(deluxe_ctx->seq_mat,subset_size,dummy_idx,subset_size,dummy_idx,fill_vals,INSERT_VALUES);CHKERRQ(ierr);
94434a97f8cSStefano Zampini     ierr = MatSetValues(global_schur_subsets,subset_size,&all_local_idx_G[local_size],subset_size,&all_local_idx_G[local_size],fill_vals,ADD_VALUES);CHKERRQ(ierr);
94534a97f8cSStefano Zampini     local_size += subset_size;
94634a97f8cSStefano Zampini   }
94734a97f8cSStefano Zampini   ierr = MatAssemblyBegin(deluxe_ctx->seq_mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
94834a97f8cSStefano Zampini   ierr = MatAssemblyEnd(deluxe_ctx->seq_mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
94934a97f8cSStefano Zampini   ierr = MatAssemblyBegin(global_schur_subsets,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
95034a97f8cSStefano Zampini   ierr = MatAssemblyEnd(global_schur_subsets,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
95134a97f8cSStefano Zampini   ierr = PetscFree2(dummy_idx,fill_vals);CHKERRQ(ierr);
95234a97f8cSStefano Zampini 
95334a97f8cSStefano Zampini   /* Create work vectors for sequential part of deluxe */
9548a26ef87SStefano Zampini   ierr = MatCreateVecs(deluxe_ctx->seq_mat,&deluxe_ctx->seq_work1,&deluxe_ctx->seq_work2);CHKERRQ(ierr);
95534a97f8cSStefano Zampini 
95634a97f8cSStefano Zampini   /* Compute deluxe sequential scatter */
95734a97f8cSStefano Zampini   ierr = ISCreateGeneral(PETSC_COMM_SELF,local_size,all_local_idx_B,PETSC_OWN_POINTER,&is_from);CHKERRQ(ierr);
95834a97f8cSStefano Zampini   ierr = VecScatterCreate(pcbddc->work_scaling,is_from,deluxe_ctx->seq_work1,NULL,&deluxe_ctx->seq_scctx);CHKERRQ(ierr);
95934a97f8cSStefano Zampini   ierr = ISDestroy(&is_from);CHKERRQ(ierr);
96034a97f8cSStefano Zampini 
96134a97f8cSStefano Zampini   /* Get local part of (\sum_j S_Ej) */
96234a97f8cSStefano Zampini   for (i=0;i<local_size;i++) {
96334a97f8cSStefano Zampini     all_local_idx_N[i] = all_local_idx_G[all_permutation_G[i]];
96434a97f8cSStefano Zampini   }
96534a97f8cSStefano Zampini   ierr = ISCreateGeneral(PetscObjectComm((PetscObject)pc),local_size,all_local_idx_N,PETSC_OWN_POINTER,&is_to);CHKERRQ(ierr);
96634a97f8cSStefano Zampini   ierr = MatGetSubMatrices(global_schur_subsets,1,&is_to,&is_to,MAT_INITIAL_MATRIX,&submat_global_schur_subsets);CHKERRQ(ierr);
96734a97f8cSStefano Zampini   ierr = MatDestroy(&global_schur_subsets);CHKERRQ(ierr);
96834a97f8cSStefano Zampini   ierr = ISDestroy(&is_to);CHKERRQ(ierr);
96934a97f8cSStefano Zampini   for (i=0;i<local_size;i++) {
97034a97f8cSStefano Zampini     all_local_idx_G[all_permutation_G[i]] = i;
97134a97f8cSStefano Zampini   }
97234a97f8cSStefano Zampini   ierr = ISCreateGeneral(PETSC_COMM_SELF,local_size,all_local_idx_G,PETSC_OWN_POINTER,&is_from);CHKERRQ(ierr);
97334a97f8cSStefano Zampini   ierr = ISSetPermutation(is_from);CHKERRQ(ierr);
97434a97f8cSStefano Zampini   ierr = MatPermute(submat_global_schur_subsets[0],is_from,is_from,&work_mat);CHKERRQ(ierr);
97534a97f8cSStefano Zampini   ierr = MatDestroyMatrices(1,&submat_global_schur_subsets);CHKERRQ(ierr);
97634a97f8cSStefano Zampini   ierr = ISDestroy(&is_from);CHKERRQ(ierr);
97734a97f8cSStefano Zampini   ierr = PetscFree(all_permutation_G);CHKERRQ(ierr);
97834a97f8cSStefano Zampini 
97934a97f8cSStefano Zampini   /* Create KSP object for sequential part of deluxe scaling */
98034a97f8cSStefano Zampini   ierr = KSPCreate(PETSC_COMM_SELF,&deluxe_ctx->seq_ksp);CHKERRQ(ierr);
98134a97f8cSStefano Zampini   ierr = KSPSetOperators(deluxe_ctx->seq_ksp,work_mat,work_mat);CHKERRQ(ierr);
98234a97f8cSStefano Zampini   ierr = KSPSetType(deluxe_ctx->seq_ksp,KSPPREONLY);CHKERRQ(ierr);
983ccdc0be0SStefano Zampini   ierr = KSPGetPC(deluxe_ctx->seq_ksp,&pc_temp);CHKERRQ(ierr);
984ccdc0be0SStefano Zampini   ierr = PCSetType(pc_temp,PCLU);CHKERRQ(ierr);
985ccdc0be0SStefano Zampini   ierr = KSPGetPC(pcbddc->ksp_D,&pc_temp);CHKERRQ(ierr);
986ccdc0be0SStefano Zampini   ierr = PCFactorGetMatSolverPackage(pc_temp,(const MatSolverPackage*)&solver);CHKERRQ(ierr);
987ccdc0be0SStefano Zampini   if (solver) {
988*a4e17c67SStefano Zampini     PC     new_pc;
989*a4e17c67SStefano Zampini     PCType type;
990*a4e17c67SStefano Zampini     ierr = PCGetType(pc_temp,&type);CHKERRQ(ierr);
991*a4e17c67SStefano Zampini     ierr = KSPGetPC(deluxe_ctx->seq_ksp,&new_pc);CHKERRQ(ierr);
992*a4e17c67SStefano Zampini     ierr = PCSetType(new_pc,type);CHKERRQ(ierr);
993*a4e17c67SStefano Zampini     ierr = PCFactorSetMatSolverPackage(new_pc,solver);CHKERRQ(ierr);
994ccdc0be0SStefano Zampini   }
99534a97f8cSStefano Zampini   ierr = PetscStrlen(((PetscObject)(pcbddc->ksp_D))->prefix,&len);CHKERRQ(ierr);
99634a97f8cSStefano Zampini   len -= 10; /* remove "dirichlet_" */
9978856534dSStefano Zampini   ierr = PetscStrncpy(ksp_prefix,((PetscObject)(pcbddc->ksp_D))->prefix,len+1);CHKERRQ(ierr);
99834a97f8cSStefano Zampini   ierr = PetscStrcat(ksp_prefix,"deluxe_seq_");CHKERRQ(ierr);
99934a97f8cSStefano Zampini   ierr = KSPSetOptionsPrefix(deluxe_ctx->seq_ksp,ksp_prefix);CHKERRQ(ierr);
100034a97f8cSStefano Zampini   ierr = KSPSetFromOptions(deluxe_ctx->seq_ksp);CHKERRQ(ierr);
100134a97f8cSStefano Zampini   ierr = KSPSetUp(deluxe_ctx->seq_ksp);CHKERRQ(ierr);
100234a97f8cSStefano Zampini   ierr = MatDestroy(&work_mat);CHKERRQ(ierr);
100334a97f8cSStefano Zampini   PetscFunctionReturn(0);
100434a97f8cSStefano Zampini }
1005