xref: /petsc/src/ksp/pc/impls/mg/gdsw.c (revision 32fe681d04f1d7625d6fc236e17224268fcf5f36)
1*32fe681dSStefano Zampini #include <petsc/private/pcmgimpl.h>
2*32fe681dSStefano Zampini #include <petsc/private/pcbddcimpl.h>
3*32fe681dSStefano Zampini #include <petsc/private/pcbddcprivateimpl.h>
4*32fe681dSStefano Zampini 
5*32fe681dSStefano Zampini static PetscErrorCode PCMGGDSWSetUp(PC pc, PetscInt l, DM dm, KSP smooth, PetscInt Nc, Mat A, PetscInt *ns, Mat **sA_IG_n, KSP **sksp_n, IS **sI_n, IS **sG_n, Mat **sGf_n, IS **sGi_n, IS **sGiM_n)
6*32fe681dSStefano Zampini {
7*32fe681dSStefano Zampini   KSP                    *sksp;
8*32fe681dSStefano Zampini   PC                     pcbddc = NULL, smoothpc;
9*32fe681dSStefano Zampini   PC_BDDC                *ipcbddc;
10*32fe681dSStefano Zampini   PC_IS                  *ipcis;
11*32fe681dSStefano Zampini   Mat                    *sA_IG, *sGf, cmat, lA;
12*32fe681dSStefano Zampini   ISLocalToGlobalMapping l2g;
13*32fe681dSStefano Zampini   IS                     *sI, *sG, *sGi, *sGiM, cref;
14*32fe681dSStefano Zampini   PCBDDCSubSchurs        sub_schurs = NULL;
15*32fe681dSStefano Zampini   PCBDDCGraph            graph;
16*32fe681dSStefano Zampini   const char             *prefix;
17*32fe681dSStefano Zampini   const PetscScalar      *tdata;
18*32fe681dSStefano Zampini   PetscScalar            *data, *cdata;
19*32fe681dSStefano Zampini   PetscReal              tol = 0.0, otol;
20*32fe681dSStefano Zampini   const PetscInt         *ia, *ja;
21*32fe681dSStefano Zampini   PetscInt               *ccii, *cridx;
22*32fe681dSStefano Zampini   PetscInt               i, j, ngct, ng, dbg = 0, odbg, minmax[2] = {0,PETSC_MAX_INT}, ominmax[2], vsize;
23*32fe681dSStefano Zampini   PetscBool              flg, userdefined = PETSC_TRUE, reuse_solver = PETSC_TRUE, reduced = PETSC_FALSE;
24*32fe681dSStefano Zampini 
25*32fe681dSStefano Zampini   PetscFunctionBegin;
26*32fe681dSStefano Zampini   PetscCall(MatGetBlockSize(A,&vsize));
27*32fe681dSStefano Zampini   PetscCall(KSPGetOptionsPrefix(smooth,&prefix));
28*32fe681dSStefano Zampini   PetscOptionsBegin(PetscObjectComm((PetscObject)smooth),prefix,"GDSW options","PC");
29*32fe681dSStefano Zampini   PetscCall(PetscOptionsReal("-gdsw_tolerance","Tolerance for eigenvalue problem",NULL,tol,&tol,NULL));
30*32fe681dSStefano Zampini   PetscCall(PetscOptionsBool("-gdsw_userdefined","Use user-defined functions in addition to those adaptively generated",NULL,userdefined,&userdefined,NULL));
31*32fe681dSStefano Zampini   PetscCall(PetscOptionsIntArray("-gdsw_minmax","Minimum and maximum number of basis functions per connected component for adaptive GDSW",NULL,minmax,(i=2,&i),NULL));
32*32fe681dSStefano Zampini   PetscCall(PetscOptionsInt("-gdsw_vertex_size","Connected components smaller or equal to vertex size will be considered as vertices",NULL,vsize,&vsize,NULL));
33*32fe681dSStefano Zampini   PetscCall(PetscOptionsBool("-gdsw_reuse","Reuse interior solver from Schur complement computations",NULL,reuse_solver,&reuse_solver,NULL));
34*32fe681dSStefano Zampini   PetscCall(PetscOptionsBool("-gdsw_reduced","Reduced GDSW",NULL,reduced,&reduced,NULL));
35*32fe681dSStefano Zampini   PetscCall(PetscOptionsInt("-gdsw_debug","Debug output",NULL,dbg,&dbg,NULL));
36*32fe681dSStefano Zampini   PetscOptionsEnd();
37*32fe681dSStefano Zampini 
38*32fe681dSStefano Zampini   PetscCall(PetscObjectTypeCompare((PetscObject)A,MATIS,&flg));
39*32fe681dSStefano Zampini   if (!flg) {
40*32fe681dSStefano Zampini     MatNullSpace nnsp;
41*32fe681dSStefano Zampini 
42*32fe681dSStefano Zampini     PetscCall(MatGetNearNullSpace(A,&nnsp));
43*32fe681dSStefano Zampini     PetscObjectReference((PetscObject)nnsp);
44*32fe681dSStefano Zampini     PetscCall(MatConvert(A,MATIS,MAT_INITIAL_MATRIX,&A));
45*32fe681dSStefano Zampini     PetscCall(MatSetNearNullSpace(A,nnsp));
46*32fe681dSStefano Zampini     PetscCall(MatNullSpaceDestroy(&nnsp));
47*32fe681dSStefano Zampini   } else PetscCall(PetscObjectReference((PetscObject)A));
48*32fe681dSStefano Zampini 
49*32fe681dSStefano Zampini   /* TODO Multi sub */
50*32fe681dSStefano Zampini   *ns = 1;
51*32fe681dSStefano Zampini   PetscCall(PetscMalloc1(*ns,&sA_IG));
52*32fe681dSStefano Zampini   PetscCall(PetscMalloc1(*ns,&sksp));
53*32fe681dSStefano Zampini   PetscCall(PetscMalloc1(*ns,&sI));
54*32fe681dSStefano Zampini   PetscCall(PetscMalloc1(*ns,&sG));
55*32fe681dSStefano Zampini   PetscCall(PetscMalloc1(*ns,&sGf));
56*32fe681dSStefano Zampini   PetscCall(PetscMalloc1(*ns,&sGi));
57*32fe681dSStefano Zampini   PetscCall(PetscMalloc1(*ns,&sGiM));
58*32fe681dSStefano Zampini   *sA_IG_n = sA_IG;
59*32fe681dSStefano Zampini   *sksp_n  = sksp;
60*32fe681dSStefano Zampini   *sI_n    = sI;
61*32fe681dSStefano Zampini   *sG_n    = sG;
62*32fe681dSStefano Zampini   *sGf_n   = sGf;
63*32fe681dSStefano Zampini   *sGi_n   = sGi;
64*32fe681dSStefano Zampini   *sGiM_n  = sGiM;
65*32fe681dSStefano Zampini 
66*32fe681dSStefano Zampini   /* submatrices and solvers */
67*32fe681dSStefano Zampini   PetscCall(KSPGetPC(smooth,&smoothpc));
68*32fe681dSStefano Zampini   PetscCall(PetscObjectTypeCompareAny((PetscObject)smoothpc,&flg,PCBDDC,""));
69*32fe681dSStefano Zampini   if (!flg) {
70*32fe681dSStefano Zampini     Mat smoothA;
71*32fe681dSStefano Zampini 
72*32fe681dSStefano Zampini     PetscCall(PCGetOperators(smoothpc,&smoothA,NULL));
73*32fe681dSStefano Zampini     PetscCall(PCCreate(PetscObjectComm((PetscObject)A),&pcbddc));
74*32fe681dSStefano Zampini     PetscCall(PetscLogObjectParent((PetscObject)pc,(PetscObject)pcbddc));
75*32fe681dSStefano Zampini     PetscCall(PCSetType(pcbddc,PCBDDC));
76*32fe681dSStefano Zampini     PetscCall(PCSetOperators(pcbddc,smoothA,A));
77*32fe681dSStefano Zampini     PetscCall(PCISSetUp(pcbddc,PETSC_TRUE,PETSC_FALSE));
78*32fe681dSStefano Zampini   } else {
79*32fe681dSStefano Zampini     PetscCall(PetscObjectReference((PetscObject)smoothpc));
80*32fe681dSStefano Zampini     pcbddc = smoothpc;
81*32fe681dSStefano Zampini   }
82*32fe681dSStefano Zampini   ipcis = (PC_IS*)pcbddc->data;
83*32fe681dSStefano Zampini   ipcbddc = (PC_BDDC*)pcbddc->data;
84*32fe681dSStefano Zampini   PetscCall(PetscObjectReference((PetscObject)ipcis->A_IB));
85*32fe681dSStefano Zampini   PetscCall(PetscObjectReference((PetscObject)ipcis->is_I_global));
86*32fe681dSStefano Zampini   PetscCall(PetscObjectReference((PetscObject)ipcis->is_B_global));
87*32fe681dSStefano Zampini   sA_IG[0] = ipcis->A_IB;
88*32fe681dSStefano Zampini   sI[0] = ipcis->is_I_global;
89*32fe681dSStefano Zampini   sG[0] = ipcis->is_B_global;
90*32fe681dSStefano Zampini 
91*32fe681dSStefano Zampini   PetscCall(KSPCreate(PetscObjectComm((PetscObject)ipcis->A_II),&sksp[0]));
92*32fe681dSStefano Zampini   PetscCall(KSPSetOperators(sksp[0],ipcis->A_II,ipcis->pA_II));
93*32fe681dSStefano Zampini   PetscCall(KSPSetOptionsPrefix(sksp[0],prefix));
94*32fe681dSStefano Zampini   PetscCall(KSPAppendOptionsPrefix(sksp[0],"gdsw_"));
95*32fe681dSStefano Zampini   PetscCall(KSPSetFromOptions(sksp[0]));
96*32fe681dSStefano Zampini 
97*32fe681dSStefano Zampini   /* analyze interface */
98*32fe681dSStefano Zampini   PetscCall(MatISGetLocalMat(A,&lA));
99*32fe681dSStefano Zampini   graph = ipcbddc->mat_graph;
100*32fe681dSStefano Zampini   if (!flg) {
101*32fe681dSStefano Zampini     PetscInt N;
102*32fe681dSStefano Zampini 
103*32fe681dSStefano Zampini     PetscCall(MatISGetLocalToGlobalMapping(A,&l2g,NULL));
104*32fe681dSStefano Zampini     PetscCall(MatGetSize(A,&N,NULL));
105*32fe681dSStefano Zampini     graph->commsizelimit = 0; /* don't use the COMM_SELF variant of the graph */
106*32fe681dSStefano Zampini     PetscCall(PCBDDCGraphInit(graph,l2g,N,PETSC_MAX_INT));
107*32fe681dSStefano Zampini     PetscCall(MatGetRowIJ(lA,0,PETSC_TRUE,PETSC_FALSE,&graph->nvtxs_csr,(const PetscInt**)&graph->xadj,(const PetscInt**)&graph->adjncy,&flg));
108*32fe681dSStefano Zampini     PetscCall(PCBDDCGraphSetUp(graph,vsize,NULL,NULL,0,NULL,NULL));
109*32fe681dSStefano Zampini     PetscCall(MatRestoreRowIJ(lA,0,PETSC_TRUE,PETSC_FALSE,&graph->nvtxs_csr,(const PetscInt**)&graph->xadj,(const PetscInt**)&graph->adjncy,&flg));
110*32fe681dSStefano Zampini     PetscCall(PCBDDCGraphComputeConnectedComponents(graph));
111*32fe681dSStefano Zampini   }
112*32fe681dSStefano Zampini   l2g = graph->l2gmap;
113*32fe681dSStefano Zampini   if (reduced) {
114*32fe681dSStefano Zampini     PetscContainer        gcand;
115*32fe681dSStefano Zampini     PCBDDCGraphCandidates cand;
116*32fe681dSStefano Zampini     PetscErrorCode        (*rgdsw)(DM,PetscInt*,IS**);
117*32fe681dSStefano Zampini 
118*32fe681dSStefano Zampini     PetscCall(PetscObjectQueryFunction((PetscObject)dm,"DMComputeLocalRGDSWSets",&rgdsw));
119*32fe681dSStefano Zampini     PetscCheck(rgdsw,PETSC_COMM_SELF,PETSC_ERR_SUP,"Not supported");
120*32fe681dSStefano Zampini     PetscCall(PetscNew(&cand));
121*32fe681dSStefano Zampini     PetscCall((*rgdsw)(dm,&cand->nfc,&cand->Faces));
122*32fe681dSStefano Zampini     /* filter interior (if any) and guarantee IS are ordered by global numbering */
123*32fe681dSStefano Zampini     for (i = 0; i < cand->nfc; i++) {
124*32fe681dSStefano Zampini       IS is,is2;
125*32fe681dSStefano Zampini 
126*32fe681dSStefano Zampini       PetscCall(ISLocalToGlobalMappingApplyIS(l2g,cand->Faces[i],&is));
127*32fe681dSStefano Zampini       PetscCall(ISDestroy(&cand->Faces[i]));
128*32fe681dSStefano Zampini       PetscCall(ISSort(is));
129*32fe681dSStefano Zampini       PetscCall(ISGlobalToLocalMappingApplyIS(l2g,IS_GTOLM_DROP,is,&is2));
130*32fe681dSStefano Zampini       PetscCall(ISDestroy(&is));
131*32fe681dSStefano Zampini       PetscCall(ISGlobalToLocalMappingApplyIS(ipcis->BtoNmap,IS_GTOLM_DROP,is2,&is));
132*32fe681dSStefano Zampini       PetscCall(ISDestroy(&is2));
133*32fe681dSStefano Zampini       PetscCall(ISLocalToGlobalMappingApplyIS(ipcis->BtoNmap,is,&cand->Faces[i]));
134*32fe681dSStefano Zampini       PetscCall(ISDestroy(&is));
135*32fe681dSStefano Zampini     }
136*32fe681dSStefano Zampini     PetscCall(PetscContainerCreate(PETSC_COMM_SELF,&gcand));
137*32fe681dSStefano Zampini     PetscCall(PetscContainerSetPointer(gcand,cand));
138*32fe681dSStefano Zampini     PetscCall(PetscContainerSetUserDestroy(gcand,PCBDDCDestroyGraphCandidatesIS));
139*32fe681dSStefano Zampini     PetscCall(PetscObjectCompose((PetscObject)l2g,"_PCBDDCGraphCandidatesIS",(PetscObject)gcand));
140*32fe681dSStefano Zampini     PetscCall(PetscContainerDestroy(&gcand));
141*32fe681dSStefano Zampini   }
142*32fe681dSStefano Zampini 
143*32fe681dSStefano Zampini   /* interface functions */
144*32fe681dSStefano Zampini   otol = ipcbddc->adaptive_threshold[1];
145*32fe681dSStefano Zampini   odbg = ipcbddc->dbg_flag;
146*32fe681dSStefano Zampini   ominmax[0] = ipcbddc->adaptive_nmin;
147*32fe681dSStefano Zampini   ominmax[1] = ipcbddc->adaptive_nmax;
148*32fe681dSStefano Zampini   ipcbddc->adaptive_threshold[1] = tol;
149*32fe681dSStefano Zampini   ipcbddc->dbg_flag = dbg;
150*32fe681dSStefano Zampini   ipcbddc->adaptive_nmin = minmax[0];
151*32fe681dSStefano Zampini   ipcbddc->adaptive_nmax = minmax[1];
152*32fe681dSStefano Zampini   if (tol != 0.0) { /* adaptive */
153*32fe681dSStefano Zampini     Mat lS;
154*32fe681dSStefano Zampini 
155*32fe681dSStefano Zampini     PetscCall(MatCreateSchurComplement(ipcis->A_II,ipcis->pA_II,ipcis->A_IB,ipcis->A_BI,ipcis->A_BB,&lS));
156*32fe681dSStefano Zampini     PetscCall(KSPGetOptionsPrefix(sksp[0],&prefix));
157*32fe681dSStefano Zampini     PetscCall(PCBDDCSubSchursCreate(&sub_schurs));
158*32fe681dSStefano Zampini     PetscCall(PCBDDCSubSchursInit(sub_schurs,prefix,ipcis->is_I_local,ipcis->is_B_local,graph,ipcis->BtoNmap,PETSC_FALSE,PETSC_TRUE));
159*32fe681dSStefano Zampini     if (userdefined) PetscCall(PCBDDCComputeFakeChange(pcbddc,PETSC_FALSE,graph,NULL,&cmat,&cref,NULL,&flg));
160*32fe681dSStefano Zampini     else { cmat = NULL; cref = NULL; }
161*32fe681dSStefano Zampini     PetscCall(PCBDDCSubSchursSetUp(sub_schurs,lA,lS,PETSC_TRUE,NULL,NULL,-1,NULL,PETSC_TRUE,reuse_solver,PETSC_FALSE,0,NULL,NULL,cmat,cref));
162*32fe681dSStefano Zampini     PetscCall(MatDestroy(&lS));
163*32fe681dSStefano Zampini     PetscCall(MatDestroy(&cmat));
164*32fe681dSStefano Zampini     PetscCall(ISDestroy(&cref));
165*32fe681dSStefano Zampini     if (sub_schurs->reuse_solver) {
166*32fe681dSStefano Zampini       PetscCall(KSPSetPC(sksp[0],sub_schurs->reuse_solver->interior_solver));
167*32fe681dSStefano Zampini       PetscCall(PCDestroy(&sub_schurs->reuse_solver->interior_solver));
168*32fe681dSStefano Zampini       sub_schurs->reuse_solver = NULL;
169*32fe681dSStefano Zampini     }
170*32fe681dSStefano Zampini   }
171*32fe681dSStefano Zampini   PetscCall(PCBDDCComputeFakeChange(pcbddc,PETSC_TRUE,graph,sub_schurs,&cmat,&cref,&sGiM[0],NULL));
172*32fe681dSStefano Zampini   PetscCall(PCBDDCSubSchursDestroy(&sub_schurs));
173*32fe681dSStefano Zampini   ipcbddc->adaptive_threshold[1] = otol;
174*32fe681dSStefano Zampini   ipcbddc->dbg_flag = odbg;
175*32fe681dSStefano Zampini   ipcbddc->adaptive_nmin = ominmax[0];
176*32fe681dSStefano Zampini   ipcbddc->adaptive_nmax = ominmax[1];
177*32fe681dSStefano Zampini 
178*32fe681dSStefano Zampini   PetscCall(ISLocalToGlobalMappingApplyIS(l2g,cref,&sGi[0]));
179*32fe681dSStefano Zampini   PetscCall(ISDestroy(&cref));
180*32fe681dSStefano Zampini 
181*32fe681dSStefano Zampini   PetscCall(MatSeqAIJGetArrayRead(cmat,&tdata));
182*32fe681dSStefano Zampini   PetscCall(MatGetRowIJ(cmat,0,PETSC_FALSE,PETSC_FALSE,&ngct,&ia,&ja,&flg));
183*32fe681dSStefano Zampini   PetscCheck(flg,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Error in MatGetRowIJ");
184*32fe681dSStefano Zampini 
185*32fe681dSStefano Zampini   PetscCall(PetscMalloc1(ngct+1,&ccii));
186*32fe681dSStefano Zampini   PetscCall(PetscMalloc1(ia[ngct],&cridx));
187*32fe681dSStefano Zampini   PetscCall(PetscMalloc1(ia[ngct],&cdata));
188*32fe681dSStefano Zampini 
189*32fe681dSStefano Zampini   PetscCall(PetscArraycpy(ccii,ia,ngct+1));
190*32fe681dSStefano Zampini   PetscCall(PetscArraycpy(cdata,tdata,ia[ngct]));
191*32fe681dSStefano Zampini   PetscCall(ISGlobalToLocalMappingApply(ipcis->BtoNmap,IS_GTOLM_DROP,ia[ngct],ja,&i,cridx));
192*32fe681dSStefano Zampini   PetscCheck(i == ia[ngct],PETSC_COMM_SELF,PETSC_ERR_PLIB,"Error in G2L");
193*32fe681dSStefano Zampini 
194*32fe681dSStefano Zampini   PetscCall(MatRestoreRowIJ(cmat,0,PETSC_FALSE,PETSC_FALSE,&i,&ia,&ja,&flg));
195*32fe681dSStefano Zampini   PetscCall(MatSeqAIJRestoreArrayRead(cmat,&tdata));
196*32fe681dSStefano Zampini   PetscCall(MatDestroy(&cmat));
197*32fe681dSStefano Zampini 
198*32fe681dSStefano Zampini   /* populate dense matrix */
199*32fe681dSStefano Zampini   PetscCall(ISGetLocalSize(sG[0],&ng));
200*32fe681dSStefano Zampini   PetscCall(MatCreateSeqDense(PETSC_COMM_SELF,ng,ngct,NULL,&sGf[0]));
201*32fe681dSStefano Zampini   PetscCall(MatDenseGetArrayWrite(sGf[0],&data));
202*32fe681dSStefano Zampini   for (i = 0; i < ngct; i++)
203*32fe681dSStefano Zampini     for (j = ccii[i]; j < ccii[i+1]; j++)
204*32fe681dSStefano Zampini       data[ng*i + cridx[j]] = cdata[j];
205*32fe681dSStefano Zampini   PetscCall(MatDenseRestoreArrayWrite(sGf[0],&data));
206*32fe681dSStefano Zampini 
207*32fe681dSStefano Zampini   PetscCall(PetscFree(cdata));
208*32fe681dSStefano Zampini   PetscCall(PetscFree(ccii));
209*32fe681dSStefano Zampini   PetscCall(PetscFree(cridx));
210*32fe681dSStefano Zampini   PetscCall(PCDestroy(&pcbddc));
211*32fe681dSStefano Zampini   PetscCall(MatDestroy(&A));
212*32fe681dSStefano Zampini   PetscFunctionReturn(0);
213*32fe681dSStefano Zampini }
214*32fe681dSStefano Zampini 
215*32fe681dSStefano Zampini PetscErrorCode PCMGGDSWCreateCoarseSpace_Private(PC pc, PetscInt l, DM dm, KSP smooth, PetscInt Nc, Mat guess, Mat* cspace)
216*32fe681dSStefano Zampini {
217*32fe681dSStefano Zampini   KSP            *sksp;
218*32fe681dSStefano Zampini   Mat            A, *sA_IG, *sGf, preallocator;
219*32fe681dSStefano Zampini   IS             Gidx, GidxMult, cG;
220*32fe681dSStefano Zampini   IS             *sI, *sG, *sGi, *sGiM;
221*32fe681dSStefano Zampini   const PetscInt *cidx;
222*32fe681dSStefano Zampini   PetscInt       NG, ns, n, i, c, rbs, cbs[2];
223*32fe681dSStefano Zampini   PetscBool      flg;
224*32fe681dSStefano Zampini   MatType        ptype;
225*32fe681dSStefano Zampini 
226*32fe681dSStefano Zampini   PetscFunctionBegin;
227*32fe681dSStefano Zampini   *cspace = NULL;
228*32fe681dSStefano Zampini   if (!l) PetscFunctionReturn(0);
229*32fe681dSStefano Zampini   if (pc->useAmat) {
230*32fe681dSStefano Zampini     PetscCall(KSPGetOperatorsSet(smooth,&flg,NULL));
231*32fe681dSStefano Zampini     PetscCheck(flg,PetscObjectComm((PetscObject)smooth),PETSC_ERR_ORDER,"Amat not set");
232*32fe681dSStefano Zampini     PetscCall(KSPGetOperators(smooth,&A,NULL));
233*32fe681dSStefano Zampini   } else {
234*32fe681dSStefano Zampini     PetscCall(KSPGetOperatorsSet(smooth,NULL,&flg));
235*32fe681dSStefano Zampini     PetscCheck(flg,PetscObjectComm((PetscObject)smooth),PETSC_ERR_ORDER,"Pmat not set");
236*32fe681dSStefano Zampini     PetscCall(KSPGetOperators(smooth,NULL,&A));
237*32fe681dSStefano Zampini   }
238*32fe681dSStefano Zampini 
239*32fe681dSStefano Zampini   /* Setup (also setup smoother here) */
240*32fe681dSStefano Zampini   if (!pc->setupcalled) PetscCall(KSPSetFromOptions(smooth));
241*32fe681dSStefano Zampini   PetscCall(KSPSetUp(smooth));
242*32fe681dSStefano Zampini   PetscCall(KSPSetUpOnBlocks(smooth));
243*32fe681dSStefano Zampini   PetscCall(PCMGGDSWSetUp(pc,l,dm,smooth,Nc,A,&ns,&sA_IG,&sksp,&sI,&sG,&sGf,&sGi,&sGiM));
244*32fe681dSStefano Zampini 
245*32fe681dSStefano Zampini   /* Number GDSW basis functions */
246*32fe681dSStefano Zampini   PetscCall(ISConcatenate(PetscObjectComm((PetscObject)A),ns,sGi,&Gidx));
247*32fe681dSStefano Zampini   PetscCall(ISConcatenate(PetscObjectComm((PetscObject)A),ns,sGiM,&GidxMult));
248*32fe681dSStefano Zampini   PetscCall(ISRenumber(Gidx,GidxMult,&NG,&cG));
249*32fe681dSStefano Zampini   PetscCall(ISDestroy(&Gidx));
250*32fe681dSStefano Zampini 
251*32fe681dSStefano Zampini   /* Detect column block size */
252*32fe681dSStefano Zampini   PetscCall(ISGetMinMax(GidxMult,&cbs[0],&cbs[1]));
253*32fe681dSStefano Zampini   PetscCall(PetscGlobalMinMaxInt(PetscObjectComm((PetscObject)A),cbs,cbs));
254*32fe681dSStefano Zampini   PetscCall(ISDestroy(&GidxMult));
255*32fe681dSStefano Zampini 
256*32fe681dSStefano Zampini   /* Construct global interpolation matrix */
257*32fe681dSStefano Zampini   PetscCall(MatGetLocalSize(A,NULL,&n));
258*32fe681dSStefano Zampini   PetscCall(MatCreate(PetscObjectComm((PetscObject)A),&preallocator));
259*32fe681dSStefano Zampini   PetscCall(MatSetSizes(preallocator,n,PETSC_DECIDE,PETSC_DECIDE,NG));
260*32fe681dSStefano Zampini   PetscCall(MatSetType(preallocator,MATPREALLOCATOR));
261*32fe681dSStefano Zampini   PetscCall(MatSetUp(preallocator));
262*32fe681dSStefano Zampini   PetscCall(ISGetIndices(cG,&cidx));
263*32fe681dSStefano Zampini   for (i = 0, c = 0; i < ns; i++) {
264*32fe681dSStefano Zampini     const PetscInt *ri, *rg;
265*32fe681dSStefano Zampini     PetscInt       nri, nrg, ncg;
266*32fe681dSStefano Zampini 
267*32fe681dSStefano Zampini     PetscCall(ISGetLocalSize(sI[i],&nri));
268*32fe681dSStefano Zampini     PetscCall(ISGetLocalSize(sG[i],&nrg));
269*32fe681dSStefano Zampini     PetscCall(ISGetIndices(sI[i],&ri));
270*32fe681dSStefano Zampini     PetscCall(ISGetIndices(sG[i],&rg));
271*32fe681dSStefano Zampini     PetscCall(MatGetSize(sGf[i],NULL,&ncg));
272*32fe681dSStefano Zampini     PetscCall(MatSetValues(preallocator,nri,ri,ncg,cidx+c,NULL,INSERT_VALUES));
273*32fe681dSStefano Zampini     PetscCall(MatSetValues(preallocator,nrg,rg,ncg,cidx+c,NULL,INSERT_VALUES));
274*32fe681dSStefano Zampini     PetscCall(ISRestoreIndices(sI[i],&ri));
275*32fe681dSStefano Zampini     PetscCall(ISRestoreIndices(sG[i],&rg));
276*32fe681dSStefano Zampini   }
277*32fe681dSStefano Zampini   PetscCall(MatAssemblyBegin(preallocator,MAT_FINAL_ASSEMBLY));
278*32fe681dSStefano Zampini   PetscCall(MatAssemblyEnd(preallocator,MAT_FINAL_ASSEMBLY));
279*32fe681dSStefano Zampini 
280*32fe681dSStefano Zampini   ptype = MATAIJ;
281*32fe681dSStefano Zampini   if (PetscDefined(HAVE_DEVICE)) {
282*32fe681dSStefano Zampini     PetscCall(MatBoundToCPU(A,&flg));
283*32fe681dSStefano Zampini     if (!flg) {
284*32fe681dSStefano Zampini       VecType vtype;
285*32fe681dSStefano Zampini       char    *found;
286*32fe681dSStefano Zampini 
287*32fe681dSStefano Zampini       PetscCall(MatGetVecType(A,&vtype));
288*32fe681dSStefano Zampini       PetscCall(PetscStrstr(vtype,"cuda",&found));
289*32fe681dSStefano Zampini       if (found) ptype = MATAIJCUSPARSE;
290*32fe681dSStefano Zampini     }
291*32fe681dSStefano Zampini   }
292*32fe681dSStefano Zampini   PetscCall(MatCreate(PetscObjectComm((PetscObject)A),cspace));
293*32fe681dSStefano Zampini   PetscCall(MatSetSizes(*cspace,n,PETSC_DECIDE,PETSC_DECIDE,NG));
294*32fe681dSStefano Zampini   PetscCall(MatSetType(*cspace,ptype));
295*32fe681dSStefano Zampini   PetscCall(MatGetBlockSizes(A,NULL,&rbs));
296*32fe681dSStefano Zampini   PetscCall(MatSetBlockSizes(*cspace,rbs,cbs[0] == cbs[1] ? cbs[0] : 1));
297*32fe681dSStefano Zampini   PetscCall(MatPreallocatorPreallocate(preallocator,PETSC_FALSE,*cspace));
298*32fe681dSStefano Zampini   PetscCall(MatDestroy(&preallocator));
299*32fe681dSStefano Zampini   PetscCall(MatSetOption(*cspace,MAT_ROW_ORIENTED,PETSC_FALSE));
300*32fe681dSStefano Zampini 
301*32fe681dSStefano Zampini   for (i = 0, c = 0; i < ns; i++) {
302*32fe681dSStefano Zampini     Mat               X, Y;
303*32fe681dSStefano Zampini     const PetscScalar *v;
304*32fe681dSStefano Zampini     const PetscInt    *ri, *rg;
305*32fe681dSStefano Zampini     PetscInt          nri, nrg, ncg;
306*32fe681dSStefano Zampini 
307*32fe681dSStefano Zampini     PetscCall(MatMatMult(sA_IG[i],sGf[i],MAT_INITIAL_MATRIX,PETSC_DEFAULT,&Y));
308*32fe681dSStefano Zampini     PetscCall(MatScale(Y,-1.0));
309*32fe681dSStefano Zampini     PetscCall(MatDuplicate(Y,MAT_DO_NOT_COPY_VALUES,&X));
310*32fe681dSStefano Zampini     PetscCall(KSPMatSolve(sksp[i],Y,X));
311*32fe681dSStefano Zampini 
312*32fe681dSStefano Zampini     PetscCall(ISGetLocalSize(sI[i],&nri));
313*32fe681dSStefano Zampini     PetscCall(ISGetLocalSize(sG[i],&nrg));
314*32fe681dSStefano Zampini     PetscCall(ISGetIndices(sI[i],&ri));
315*32fe681dSStefano Zampini     PetscCall(ISGetIndices(sG[i],&rg));
316*32fe681dSStefano Zampini     PetscCall(MatGetSize(sGf[i],NULL,&ncg));
317*32fe681dSStefano Zampini 
318*32fe681dSStefano Zampini     PetscCall(MatDenseGetArrayRead(X,&v));
319*32fe681dSStefano Zampini     PetscCall(MatSetValues(*cspace,nri,ri,ncg,cidx+c,v,INSERT_VALUES));
320*32fe681dSStefano Zampini     PetscCall(MatDenseRestoreArrayRead(X,&v));
321*32fe681dSStefano Zampini     PetscCall(MatDenseGetArrayRead(sGf[i],&v));
322*32fe681dSStefano Zampini     PetscCall(MatSetValues(*cspace,nrg,rg,ncg,cidx+c,v,INSERT_VALUES));
323*32fe681dSStefano Zampini     PetscCall(MatDenseRestoreArrayRead(sGf[i],&v));
324*32fe681dSStefano Zampini     PetscCall(ISRestoreIndices(sI[i],&ri));
325*32fe681dSStefano Zampini     PetscCall(ISRestoreIndices(sG[i],&rg));
326*32fe681dSStefano Zampini     PetscCall(MatDestroy(&Y));
327*32fe681dSStefano Zampini     PetscCall(MatDestroy(&X));
328*32fe681dSStefano Zampini   }
329*32fe681dSStefano Zampini   PetscCall(ISRestoreIndices(cG,&cidx));
330*32fe681dSStefano Zampini   PetscCall(ISDestroy(&cG));
331*32fe681dSStefano Zampini   PetscCall(MatAssemblyBegin(*cspace,MAT_FINAL_ASSEMBLY));
332*32fe681dSStefano Zampini   PetscCall(MatAssemblyEnd(*cspace,MAT_FINAL_ASSEMBLY));
333*32fe681dSStefano Zampini 
334*32fe681dSStefano Zampini   for (i = 0; i < ns; i++) {
335*32fe681dSStefano Zampini     PetscCall(KSPDestroy(&sksp[i]));
336*32fe681dSStefano Zampini     PetscCall(ISDestroy(&sI[i]));
337*32fe681dSStefano Zampini     PetscCall(ISDestroy(&sG[i]));
338*32fe681dSStefano Zampini     PetscCall(ISDestroy(&sGi[i]));
339*32fe681dSStefano Zampini     PetscCall(ISDestroy(&sGiM[i]));
340*32fe681dSStefano Zampini     PetscCall(MatDestroy(&sGf[i]));
341*32fe681dSStefano Zampini     PetscCall(MatDestroy(&sA_IG[i]));
342*32fe681dSStefano Zampini   }
343*32fe681dSStefano Zampini   PetscCall(PetscFree(sksp));
344*32fe681dSStefano Zampini   PetscCall(PetscFree(sI));
345*32fe681dSStefano Zampini   PetscCall(PetscFree(sG));
346*32fe681dSStefano Zampini   PetscCall(PetscFree(sGi));
347*32fe681dSStefano Zampini   PetscCall(PetscFree(sGiM));
348*32fe681dSStefano Zampini   PetscCall(PetscFree(sGf));
349*32fe681dSStefano Zampini   PetscCall(PetscFree(sA_IG));
350*32fe681dSStefano Zampini   PetscFunctionReturn(0);
351*32fe681dSStefano Zampini }
352