<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=""><br class=""><div><br class=""><blockquote type="cite" class=""><div class="">On Dec 30, 2018, at 6:30 PM, Manav Bhatia <<a href="mailto:bhatiamanav@gmail.com" class="">bhatiamanav@gmail.com</a>> wrote:</div><br class="Apple-interchange-newline"><div class=""><meta http-equiv="Content-Type" content="text/html; charset=utf-8" class=""><div style="word-wrap: break-word; -webkit-nbsp-mode: space; line-break: after-white-space;" class=""><div class="">Stefano, </div><div class=""><br class=""></div><div class="">   Thanks for the clarification.  </div><div class=""><br class=""></div><div class="">   If it is just about the index sets, then this may be relatively straightforward to setup since the information is already present in libMesh. </div><div class=""> </div><div class="">   A few more questions: </div><div class=""><br class=""></div><div class="">— If there are mixed types of elements (2D plates and 1D beams in a built-up wing structure) in the model, is it enough to define the index sets or does there need to be any information about types of elements? </div><div class=""><br class=""></div></div></div></blockquote><div><br class=""></div><div>When you have multiple fields, you can specify the index sets of the dots associated to each field via PCBDDCSetDofsSplitting. When this information is present, PCBDDC will define edge and faces for each field separately. Have you read this <a href="https://epubs.siam.org/doi/10.1137/15M1025785?" class="">https://epubs.siam.org/doi/10.1137/15M1025785?</a></div><br class=""><blockquote type="cite" class=""><div class=""><div style="word-wrap: break-word; -webkit-nbsp-mode: space; line-break: after-white-space;" class=""><div class="">— A structural analysis with shells typically has 6 variables: three displacements and three rotations. In this case, is any extra information required, such as  the block structure, or just the typical connectivity information suffice? </div></div></div></blockquote><div><br class=""></div>I don’t have experience with shells. I’m aware of a couple of publications about BDDC for shells. I suggest you to take a look at those.<br class=""><blockquote type="cite" class=""><div class=""><div style="word-wrap: break-word; -webkit-nbsp-mode: space; line-break: after-white-space;" class=""><div class=""><br class=""></div><div class="">— The MatSetLocalToGlobalMapping seems to primarily be concerned about index sets and the relation of local and global indexing (?). If so, then I am a bit confused by this comment in the PCBDDC manual page: “<span style="font-family: -webkit-standard; font-size: inherit; background-color: rgb(255, 255, 255);" class="">Boundary nodes are split in vertices, edges and faces classes using information from the local to global mapping of dofs and the local connectivity graph of nodes.</span>"</div><div class=""><br class=""></div><div class="">   Is it necessary to provide this information about edges/faces or does the solver figure it out from the indexing? Or is this optional? </div><div class=""><br class=""></div></div></div></blockquote><blockquote type="cite" class=""><div class=""><div style="word-wrap: break-word; -webkit-nbsp-mode: space; line-break: after-white-space;" class=""><div class="">— The previous point is in relation to PCBDDCSetLocalAdjacencyGraph. If this is optional then for what types of problems would you recommend setting this up? </div></div></div></blockquote><div><br class=""></div><div>Edges and faces are identified using connected component analysis. You can pass additional information to refine these sets, such as the fields specification (PCBDDCSetDofsSplitting), a custom connectivity graph between dots (PCBDDCSetLocalAdjacencyGraph), boundary dofs specification (PCBDDCSetDirichlet/NeumannBoundaries). You can also specify specific interface dofs as primal via PCBDDCSetPrimalVerticesIS. I suggest you to take a look at the paper about PCBDDC.</div><div class=""><br class=""></div><br class=""><blockquote type="cite" class=""><div class=""><div style="word-wrap: break-word; -webkit-nbsp-mode: space; line-break: after-white-space;" class=""><div class=""><br class=""></div><div class="">— I am particularly interested in applying the solver to large-scale models of built-up structures (for example, wings), which are made-up of different types of elements (1D/2D) and also involves non-manifold mesh due to internal spars/ribs connecting to upper and lower skins. Is there any guidance to setup the solver for these types of problems, or would the approach identified in your previous email be sufficient?</div></div></div></blockquote><div><br class=""></div>No guidance.</div><div><br class=""><blockquote type="cite" class=""><div class=""><div style="word-wrap: break-word; -webkit-nbsp-mode: space; line-break: after-white-space;" class=""><div class=""><br class=""></div></div></div></blockquote><blockquote type="cite" class=""><div class=""><div style="word-wrap: break-word; -webkit-nbsp-mode: space; line-break: after-white-space;" class=""><div class="">   Thank you for your time and guidance. </div><div class=""><br class=""></div><div class="">Regards,</div><div class="">Manav</div><br class=""><div class=""><br class=""><blockquote type="cite" class=""><div class="">On Dec 30, 2018, at 8:53 AM, Stefano Zampini <<a href="mailto:stefano.zampini@gmail.com" class="">stefano.zampini@gmail.com</a>> wrote:</div><br class="Apple-interchange-newline"><div class=""><meta http-equiv="Content-Type" content="text/html; charset=utf-8" class=""><div 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 class="">4) MatSetLocalToGlobalMapping</div><div class="">5) MatISSetPreallocation</div><div class=""><br class=""></div><div class="">Then, either MatSetValues or MatSetValuesLocal can be called and the Neumann problems will be automatically assembled.</div><div class="">Most of the matrix operations are supported (e.g. MatMult, MatZeroRows, etc)</div><div class=""><br class=""></div><div class="">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 class="">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 class=""><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 class="">Some notes:</div><div class=""><br class=""></div><div class="">- 2) and 3) can be swapped</div><div class="">- 5) MatXAIJSetPreallocation can be called instead</div><div class=""><br class=""></div><div class="">Let me know if you need further clarifications</div><div class="">Stefano</div><div class=""><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></div></div></blockquote></div><br class=""></div></div></blockquote></div><br class=""></body></html>