<div dir="ltr"><div dir="ltr"><div dir="ltr"><div dir="ltr"><div dir="ltr"><div dir="ltr"><div dir="ltr"><div dir="ltr"><div dir="ltr"><div dir="ltr"><div dir="ltr">Thanks a lot, Stefano. <div>Yes. I am using 2nd-order Nedelec elements. <span style="color:rgb(0,0,0);font-family:Menlo;font-size:21px;font-variant-ligatures:no-common-ligatures">-pc_</span><span style="font-family:Menlo;font-size:21px;font-variant-ligatures:no-common-ligatures;color:rgb(180,36,25)"><b>bddc_use_lo</b></span><span style="color:rgb(0,0,0);font-family:Menlo;font-size:21px;font-variant-ligatures:no-common-ligatures">cal_mat_graph 0 </span><span style="color:rgb(0,0,0);font-family:Menlo;font-variant-ligatures:no-common-ligatures"><font size="4">can make the code run. I am testing more cpu #. </font></span></div><div><span style="color:rgb(0,0,0);font-family:Menlo;font-variant-ligatures:no-common-ligatures"><font size="4">I am testing my code using, </font></span></div><div><div dir="ltr">petsc-3.21.1/petsc/arch-linux-c-opt/bin/mpirun -n 4 ./app -pc_type bddc -pc_bddc_coarse_redundant_pc_type svd  -pc_bddc_use_vertices -ksp_error_if_not_converged -mat_type is -ksp_monitor -ksp_rtol 1e-8 -ksp_gmres_restart 2000 -ksp_view -malloc_view -pc_bddc_use_local_mat_graph 0 -ksp_converged_reason -pc_bddc_neumann_pc_type gamg -pc_bddc_neumann_pc_gamg_esteig_ksp_max_it 10 -ksp_converged_reason -pc_bddc_neumann_approximate -pc_bddc_dirichlet_pc_type gamg -pc_bddc_dirichlet_pc_gamg_esteig_ksp_max_it 10 -ksp_converged_reason -pc_bddc_dirichlet_approximate<br></div><div class="gmail-HOEnZb"><div class="gmail-adm"><div id="gmail-q_9" class="gmail-ajR gmail-h4"></div></div></div></div><div><span style="color:rgb(0,0,0);font-family:Menlo;font-variant-ligatures:no-common-ligatures"><font size="4"><br></font></span></div><div><font color="#000000" face="Menlo" size="4"><span style="font-variant-ligatures:no-common-ligatures">The residual dropped to 6e-5 very fast and then continued to reduce very slowly.  Do you have any suggestions to improve this ?</span></font></div><div><font color="#000000" face="Menlo" size="4"><span style="font-variant-ligatures:no-common-ligatures"><br></span></font></div><div><span style="color:rgb(0,0,0);font-family:Menlo;font-variant-ligatures:no-common-ligatures">Will it be necessary to </span>change the basis for BDDC in order to accelerate the convergence ? In addition, I tried -pc_bddc_use_deluxe_scaling, but it showed some errors. It seems deluxe scaling obviously requires a much larger size (<b>Global size overflow 3051678564</b>) than my problem. </div><div><br></div><div>Thanks,</div><div>Xiaodong </div><div><div dir="ltr"><br></div><div dir="ltr">[0]PETSC ERROR: --------------------- Error Message --------------------------------------------------------------<br>[0]PETSC ERROR: Overflow in integer operation: <a href="https://urldefense.us/v3/__https://petsc.org/release/faq/*64-bit-indices__;Iw!!G_uCfscf7eWS!ZVyzxJb4s9N1kzsS2BV7raG-kJIn8X6skBNtfsvA8aHyjWPm8oYGfzk83j1n0PFstGE6nDCHpOIpMvkLFZcexA$" target="_blank">https://petsc.org/release/faq/#64-bit-indices</a><br>[0]PETSC ERROR: Global size overflow 3051678564. You may consider ./configure PETSc with --with-64-bit-indices for the case you are running<br>[0]PETSC ERROR: WARNING! There are unused option(s) set! Could be the program crashed before usage or a spelling mistake, etc!<br>[0]PETSC ERROR:   Option left: name:-ksp_converged_reason (no value) source: command line<br>[0]PETSC ERROR:   Option left: name:-pc_bddc_coarse_redundant_pc_type value: svd source: command line<br>[0]PETSC ERROR:   Option left: name:-pc_bddc_neumann_pc_gamg_esteig_ksp_max_it value: 10 source: command line<br>[0]PETSC ERROR:   Option left: name:-pc_bddc_neumann_pc_type value: gamg source: command line<br>[0]PETSC ERROR: See <a href="https://urldefense.us/v3/__https://petsc.org/release/faq/__;!!G_uCfscf7eWS!ZVyzxJb4s9N1kzsS2BV7raG-kJIn8X6skBNtfsvA8aHyjWPm8oYGfzk83j1n0PFstGE6nDCHpOIpMvlAB6PtAg$" target="_blank">https://petsc.org/release/faq/</a> for trouble shooting.<br>[0]PETSC ERROR: Petsc Release Version 3.21.1, unknown<br>[0]PETSC ERROR: Configure options --with-cc=gcc --with-fc=gfortran --with-cxx=g++ --download-fblaslapack --download-mpich --with-scalar-type=complex --download-triangle --with-debugging=no<br>[0]PETSC ERROR: #1 PetscSplitOwnership() at /Documents/petsc-3.21.1/petsc/src/sys/utils/psplit.c:86<br>[0]PETSC ERROR: #2 PetscLayoutSetUp() at  /Documents/petsc-3.21.1/petsc/src/vec/is/utils/pmap.c:244<br>[0]PETSC ERROR: #3 PetscLayoutCreateFromSizes() at /Documents/petsc-3.21.1/petsc/src/vec/is/utils/pmap.c:107<br>[0]PETSC ERROR: #4 ISGeneralSetIndices_General() at  /Documents/petsc-3.21.1/petsc/src/vec/is/is/impls/general/general.c:569<br>[0]PETSC ERROR: #5 ISGeneralSetIndices() at /Documents/petsc-3.21.1/petsc/src/vec/is/is/impls/general/general.c:559<br>[0]PETSC ERROR: #6 ISCreateGeneral() at /Documents/petsc-3.21.1/petsc/src/vec/is/is/impls/general/general.c:530<br>[0]PETSC ERROR: #7 ISRenumber() at /Documents/petsc-3.21.1/petsc/src/vec/is/is/interface/index.c:198<br>[0]PETSC ERROR: #8 PCBDDCSubSchursSetUp() at /Documents/petsc-3.21.1/petsc/src/ksp/pc/impls/bddc/bddcschurs.c:646<br>[0]PETSC ERROR: #9 PCBDDCSetUpSubSchurs() at /Documents/petsc-3.21.1/petsc/src/ksp/pc/impls/bddc/bddcprivate.c:9348<br>[0]PETSC ERROR: #10 PCSetUp_BDDC() at /Documents/petsc-3.21.1/petsc/src/ksp/pc/impls/bddc/bddc.c:1564<br>[0]PETSC ERROR: #11 PCSetUp() at /Documents/petsc-3.21.1/petsc/src/ksp/pc/interface/precon.c:1079<br>[0]PETSC ERROR: #12 KSPSetUp() at /Documents/petsc-3.21.1/petsc/src/ksp/ksp/interface/itfunc.c:415<br>[0]PETSC ERROR: #13 KSPSolve_Private() at Documents/petsc-3.21.1/petsc/src/ksp/ksp/interface/itfunc.c:831<br>[0]PETSC ERROR: #14 KSPSolve() at /Documents/petsc-3.21.1/petsc/src/ksp/ksp/interface/itfunc.c:1078<br></div><div><div><br></div><div class="gmail-HOEnZb"><div class="gmail-adm"><div id="gmail-q_9" class="gmail-ajR gmail-h4"></div></div></div></div><div class="gmail-HOEnZb"><div class="gmail-adm"><div id="gmail-q_9" class="gmail-ajR gmail-h4"></div></div></div></div></div></div></div></div></div></div></div></div></div></div></div><br><div class="gmail_quote"><div dir="ltr" class="gmail_attr">On Wed, Aug 14, 2024 at 11:54 AM Stefano Zampini <<a href="mailto:stefano.zampini@gmail.com">stefano.zampini@gmail.com</a>> wrote:<br></div><blockquote class="gmail_quote" style="margin:0px 0px 0px 0.8ex;border-left:1px solid rgb(204,204,204);padding-left:1ex"><div dir="ltr">Ok, the problem is that the default algorithm for detecting the connected components of the interface finds a lot of disconnected dofs.<div>What discretization is this? Nedelec elements? Can you try using <span style="color:rgb(0,0,0);font-family:Menlo;font-size:21px;font-variant-ligatures:no-common-ligatures">-pc_</span><span style="font-family:Menlo;font-size:21px;font-variant-ligatures:no-common-ligatures;color:rgb(180,36,25)"><b>bddc_use_lo</b></span><span style="color:rgb(0,0,0);font-family:Menlo;font-size:21px;font-variant-ligatures:no-common-ligatures">cal_mat_graph 0?</span></div><div>Also, you are using -pc_bddc_monolithic, but you only have one field. That flag aggregates different fields, but you only have one. </div><div>Note that with Nedelec elements, you need a special change of basis for BDDC to work, see e.g. <a href="https://urldefense.us/v3/__https://www.osti.gov/servlets/purl/1377770__;!!G_uCfscf7eWS!ZVyzxJb4s9N1kzsS2BV7raG-kJIn8X6skBNtfsvA8aHyjWPm8oYGfzk83j1n0PFstGE6nDCHpOIpMvlI_QH81A$" target="_blank">https://www.osti.gov/servlets/purl/1377770</a></div>





</div><br><div class="gmail_quote"><div dir="ltr" class="gmail_attr">Il giorno mer 14 ago 2024 alle ore 05:15 neil liu <<a href="mailto:liufield@gmail.com" target="_blank">liufield@gmail.com</a>> ha scritto:<br></div><blockquote class="gmail_quote" style="margin:0px 0px 0px 0.8ex;border-left:1px solid rgb(204,204,204);padding-left:1ex"><div dir="ltr">Hi, Stefano, <div><br></div><div>Please see the attached for the smaller case(successful with BDDC).</div><div>and the Error_largerMesh shows the error with the large mesh using petsc debug mode. </div><div><br></div><div>Thanks a lot,</div><div><br></div><div>Xiaodong</div><div><br></div></div><br><div class="gmail_quote"><div dir="ltr" class="gmail_attr">On Tue, Aug 13, 2024 at 5:47 PM Stefano Zampini <<a href="mailto:stefano.zampini@gmail.com" target="_blank">stefano.zampini@gmail.com</a>> wrote:<br></div><blockquote class="gmail_quote" style="margin:0px 0px 0px 0.8ex;border-left:1px solid rgb(204,204,204);padding-left:1ex"><div dir="ltr">can you run the same options and add "-ksp_view -pc_bddc_check_level 1" for the smaller case? Also, can you send the full stack trace of the out-of-memory error using a debug version of PETSc?<div>A note aside: you should not need pc_bddc_use_vertices (which is on by default)</div></div><br><div class="gmail_quote"><div dir="ltr" class="gmail_attr">Il giorno mar 13 ago 2024 alle ore 23:17 neil liu <<a href="mailto:liufield@gmail.com" target="_blank">liufield@gmail.com</a>> ha scritto:<br></div><blockquote class="gmail_quote" style="margin:0px 0px 0px 0.8ex;border-left:1px solid rgb(204,204,204);padding-left:1ex"><div dir="ltr"><div dir="ltr"><div dir="ltr"><div dir="ltr"><div dir="ltr"><div dir="ltr"><div dir="ltr">Dear Petsc developers,<div><br></div><div>I am testing PCBDDC for my vector based FEM solver(complex system). It can work well on a coarse mesh(tetrahedra cell #: 6,108; dof # : 39,596). Then I tried a finer mesh (tetrahedra cell #: 32,036; dof # : 206,362). It seems ASM can work well with </div><div><br></div><div><div>petsc-3.21.1/petsc/arch-linux-c-opt/bin/mpirun -n 4 ./app -pc_type asm  -ksp_converged_reason -ksp_monitor  -ksp_gmres_restart 100 -ksp_rtol 1e-4 -pc_asm_overalp 4 -sub_pc_type ilu -malloc_view<br></div><div><br></div></div><div>while  PCBDDC eats up the memory (61 GB) when I tried </div><div><br></div><div>petsc-3.21.1/petsc/arch-linux-c-opt/bin/mpirun -n 4 ./app -pc_type bddc -pc_bddc_coarse_redundant_pc_type ilu  -pc_bddc_use_vertices -ksp_error_if_not_converged -mat_type is -ksp_monitor -ksp_rtol 1e-8 -ksp_gmres_restart 30 -ksp_view -malloc_view -pc_bddc_monolithic -pc_bddc_neumann_pc_type ilu -pc_bddc_dirichlet_pc_type ilu <br></div><div><br></div><div>The following errors with BDDC came out. The memory usage for PCBDDC (different from PCASM) is also listed (I am assuming the unit is Bytes, right?). <b>Although the BDDC requires more memory, it still seems normal, right? </b></div><div><br></div><div><div>[0]PETSC ERROR: --------------------- Error Message --------------------------------------------------------------</div><div>[0]PETSC ERROR: Out of memory. This could be due to allocating</div><div>[0]PETSC ERROR: too large an object or bleeding by not properly</div><div>[0]PETSC ERROR: destroying unneeded objects.</div><div>[0] Maximum memory PetscMalloc()ed 30829727808 maximum size of entire process 16899194880</div><div>[0] Memory usage sorted by function</div></div><div>....</div><div><div><b>[0] 1 240 PCBDDCGraphCreate()</b></div><div><b>[0] 1 3551136 PCBDDCGraphInit()</b></div><div><b>[0] 2045 32720 PCBDDCGraphSetUp()</b></div><div><b>[0] 2 8345696 PCBDDCSetLocalAdjacencyGraph_BDDC()</b></div><div><b>[0] 1 784 PCCreate()</b></div><div><b>[0] 1 1216 PCCreate_BDDC()</b></div></div><div>....</div><div><br></div><div>Thanks for your help. </div><div><br></div><div>Xiaodong </div><div><div><br></div><div><br></div></div></div></div></div></div></div></div></div>
</blockquote></div><br clear="all"><div><br></div><span class="gmail_signature_prefix">-- </span><br><div dir="ltr" class="gmail_signature">Stefano</div>
</blockquote></div>
</blockquote></div><br clear="all"><div><br></div><span class="gmail_signature_prefix">-- </span><br><div dir="ltr" class="gmail_signature">Stefano</div>
</blockquote></div>