<div dir="ltr"><div dir="ltr">Abdullah,<div><br></div><div>The "Neumann" problems Jed is referring to result from assembling your problem on each subdomain ( = MPI process) separately.</div><div>Assuming you are using FEM, these problems have been historically  named "Neumann" as they correspond to a problem with natural boundary conditions (Neumann bc for Poisson).</div><div>Note that in PETSc the subdomain decomposition is associated with the mesh decomposition.</div><div><br></div><div>When converting from an assembled AIJ matrix to a MATIS format, such "Neumann" information is lost.</div><div>You can disassemble an AIJ matrix, in the sense that you can find local matrices A_j such that A = \sum_j R^T_j A_j R_j (as it is done in ex72.c), but you cannot guarantee (unless if you solve an optimization problem) that the disassembling will produce subdomain Neumann problems that are consistent with your FEM problem.</div><div><br></div><div>I have added such disassembling code a few months ago, just to have another alternative for preconditioning AIJ matrices in PETSc; there are few tweaks one can do to improve the quality of the disassembling, but I discourage its usage unless you don't have access to the FEM assembly code.</div><div><br></div><div>With that said, what problem are you trying to solve? Are you using DMDA or DMPlex? What are the results you obtained with using the automatic disassembling?</div></div></div><br><div class="gmail_quote"><div dir="ltr">Il giorno mer 24 ott 2018 alle ore 08:14 Abdullah Ali Sivas <<a href="mailto:abdullahasivas@gmail.com">abdullahasivas@gmail.com</a>> ha scritto:<br></div><blockquote class="gmail_quote" style="margin:0 0 0 .8ex;border-left:1px #ccc solid;padding-left:1ex"><div dir="ltr"><div>Hi Jed,</div><div><br></div><div>Thanks for your reply. The assembled matrix I have corresponds to the full problem on the full mesh. There are no "Neumann" problems (or any sort of domain decomposition) defined in the code generates the matrix. However, I think assembling the full problem is equivalent to implicitly assembling the "Neumann" problems, since the system can be partitioned as;</div><div><br></div><div>[A_{LL} | A_{LI}]  [u_L]     [F]</div><div>-----------|------------ -------- = -----<br></div><div>[A_{IL}  |A_{II} ]   [u_I]      [G]<br></div><div><br></div><div>and G should correspond to the Neumann problem. I might be thinking wrong (or maybe I completely misunderstood the idea), if so please correct me. But I think that the problem is that I am not explicitly telling PCBDDC which dofs are interface dofs.</div><div><br></div><div>Regards,</div><div>Abdullah Ali Sivas<br></div></div><br><div class="gmail_quote"><div dir="ltr">On Tue, 23 Oct 2018 at 23:16, Jed Brown <<a href="mailto:jed@jedbrown.org" target="_blank">jed@jedbrown.org</a>> wrote:<br></div><blockquote class="gmail_quote" style="margin:0 0 0 .8ex;border-left:1px #ccc solid;padding-left:1ex">Did you assemble "Neumann" problems that are compatible with your<br>
definition of interior/interface degrees of freedom?<br>
<br>
Abdullah Ali Sivas <<a href="mailto:abdullahasivas@gmail.com" target="_blank">abdullahasivas@gmail.com</a>> writes:<br>
<br>
> Dear all,<br>
><br>
> I have a series of linear systems coming from a PDE for which BDDC is an<br>
> optimal preconditioner. These linear systems are assembled and I read them<br>
> from a file, then convert into MATIS as required (as in<br>
> <a href="https://www.mcs.anl.gov/petsc/petsc-current/src/ksp/ksp/examples/tutorials/ex72.c.html" rel="noreferrer" target="_blank">https://www.mcs.anl.gov/petsc/petsc-current/src/ksp/ksp/examples/tutorials/ex72.c.html</a><br>
> ). I expect each of the systems converge to the solution in almost same<br>
> number of iterations but I don't observe it. I think it is because I do not<br>
> provide enough information to the preconditioner. I can get a list of inner<br>
> dofs and interface dofs. However, I do not know how to use them. Has anyone<br>
> have any insights about it or done something similar?<br>
><br>
> Best wishes,<br>
> Abdullah Ali Sivas<br>
</blockquote></div>
</blockquote></div><br clear="all"><div><br></div>-- <br><div dir="ltr" class="gmail_signature" data-smartmail="gmail_signature">Stefano</div>