<html><head><meta http-equiv="Content-Type" content="text/html; charset=utf-8"></head><body style="word-wrap: break-word; -webkit-nbsp-mode: space; line-break: after-white-space;" class="">Manav,<div class=""><br class=""></div><div class="">The conversion MATAIJ -> MATIS is suboptimal if you can generate the matrices (e.g. with libMesh).</div><div class=""><br class=""></div><div class="">For best performance of these methods, you really need to assemble the subdomain Neumann problems.</div><div class="">This has been done for other libraries (e.g. FEniCS, MFEM, FreeFem++, I also did it for deal.II at some point), and it can be done from libMesh too. I’m not aware of the low level details of the library, but I can tell you how the MATIS object should be properly setup.</div><div class=""><br class=""></div><div class="">This is the sequence of calls (they don’t need to be in the same function call, but the order should be respected)</div><div class=""><br class=""></div><div class="">1) MatCreate</div><div class="">2) MatSetSizes</div><div class="">3) MatSetType<br class=""><div>4) MatSetLocalToGlobalMapping</div><div>5) MatISSetPreallocation</div><div><br class=""></div><div>Then, either MatSetValues or MatSetValuesLocal can be called and the Neumann problems will be automatically assembled.</div><div>Most of the matrix operations are supported (e.g. MatMult, MatZeroRows, etc)</div><div><br class=""></div><div>The crucial part is 4); all distributed memory FEM  libraries have this information hidden somewhere, as it is needed to properly assemble the matrix anyway in a standard local loop on the elements owned by the MPI Process.</div><div>It is a list of indices in the global matrix that are owned by the local MPI process (e.g. corresponding to the dofs of the elements that are in the interior of the local part of the mesh) or shared (e.g. does of the elements lying on the boundary of the local mesh).</div><div><div class="">Probably, libMesh experts are reading this message and they can tell you where to find the information you need.</div><div class=""><br class=""></div></div><div>Some notes:</div><div><br class=""></div><div>- 2) and 3) can be swapped</div><div>- 5) MatXAIJSetPreallocation can be called instead</div><div><br class=""></div><div>Let me know if you need further clarifications</div><div>Stefano</div><div><br class=""><blockquote type="cite" class=""><div class="">On Dec 28, 2018, at 9:05 PM, Manav Bhatia via petsc-users <<a href="mailto:petsc-users@mcs.anl.gov" class="">petsc-users@mcs.anl.gov</a>> wrote:</div><br class="Apple-interchange-newline"><div class=""><meta http-equiv="Content-Type" content="text/html; charset=us-ascii" class=""><div style="word-wrap: break-word; -webkit-nbsp-mode: space; line-break: after-white-space;" class=""><div class="">Hi, </div><div class=""><br class=""></div><div class="">   I am currently accessing PETSc through libMesh where the matrices of type MPIAIJ are being initialized using the MatSetSizes, MatSetFromOptions, and MatMPIAISetPreallocation calls. </div><div class=""><br class=""></div><div class="">   I am interested in accessing the BDDC solver and have looked through the PDF manual that talks about BDDC setup through MATIS. However, it is not clear what would be the most efficient route to accessing BDDC for my unstructured FE code. </div><div class=""><br class=""></div><div class="">   Example 72 (<a href="https://www.mcs.anl.gov/petsc/petsc-3.10/src/ksp/ksp/examples/tutorials/ex72.c.html" class="">https://www.mcs.anl.gov/petsc/petsc-3.10/src/ksp/ksp/examples/tutorials/ex72.c.html</a>) seems to be calling MatConvert to make the conversion to MATIS. However, the documentation of MATIS says that <a href="https://www.mcs.anl.gov/petsc/petsc-3.10/docs/manualpages/Mat/MatSetLocalToGlobalMapping.html#MatSetLocalToGlobalMapping" style="font-family: -webkit-standard;" class="">MatSetLocalToGlobalMapping</a> should be called before initializing MATIS, which is not done in Example 72. I am also not sure if I need to setup any DM/DA constructs to aid in this process. </div><div class=""><br class=""></div><div class="">   I would appreciate if someone could provide an outline of steps required to make the conversion of MatMPIAIJ in order to use BDDC/FETI preconditioners. </div><div class=""><br class=""></div><div class="">Thanks,</div><div class="">Manav</div><div class=""><br class=""></div></div></div></blockquote></div><br class=""></div></body></html>