<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=""><div class=""><br class=""></div>  Sophie,<div class=""><br class=""></div><div class="">   Sorry, looks like an old bug in PETSc that was undetected due to lack of use. The code is trying to use the first of the two matrices to determine the preconditioner which won't work in your case since it is matrix-free. It should be using the second matrix.</div><div class=""><br class=""></div><div class="">  I hope the branch <span style="font-family: Menlo; font-size: 14px;" class="">barry/2020-09-01/fix-fieldsplit-mf </span>resolves this issue for you.</div><div class=""><br class=""></div><div class="">  Barry</div><div class=""><br class=""><div><br class=""><blockquote type="cite" class=""><div class="">On Sep 1, 2020, at 12:45 PM, Blondel, Sophie <<a href="mailto:sblondel@utk.edu" class="">sblondel@utk.edu</a>> wrote:</div><br class="Apple-interchange-newline"><div class=""><div style="font-style: normal; font-variant-caps: normal; font-weight: normal; letter-spacing: normal; text-align: start; text-indent: 0px; text-transform: none; white-space: normal; word-spacing: 0px; -webkit-text-stroke-width: 0px; text-decoration: none; font-family: Calibri, Arial, Helvetica, sans-serif; font-size: 12pt;" class="">Hi Barry,</div><div style="font-style: normal; font-variant-caps: normal; font-weight: normal; letter-spacing: normal; text-align: start; text-indent: 0px; text-transform: none; white-space: normal; word-spacing: 0px; -webkit-text-stroke-width: 0px; text-decoration: none; font-family: Calibri, Arial, Helvetica, sans-serif; font-size: 12pt;" class=""><br class=""></div><div style="font-style: normal; font-variant-caps: normal; font-weight: normal; letter-spacing: normal; text-align: start; text-indent: 0px; text-transform: none; white-space: normal; word-spacing: 0px; -webkit-text-stroke-width: 0px; text-decoration: none; font-family: Calibri, Arial, Helvetica, sans-serif; font-size: 12pt;" class="">I'm working through step 1) but I think I am doing something wrong. I'm using DMDASetBlockFillsSparse to set the non-zeros only for the diffusing clusters (small He clusters here, from size 1 to 7) and all the diagonal entries. Then I added a few lines in the code:<br class=""><div class="">Mat mat;</div><div class="">DMCreateMatrix(da, &mat);</div>MatSetOption(mat,MAT_NEW_NONZERO_LOCATIONS,PETSC_FALSE);<br class=""></div><div style="font-style: normal; font-variant-caps: normal; font-weight: normal; letter-spacing: normal; text-align: start; text-indent: 0px; text-transform: none; white-space: normal; word-spacing: 0px; -webkit-text-stroke-width: 0px; text-decoration: none; font-family: Calibri, Arial, Helvetica, sans-serif; font-size: 12pt;" class=""><br class=""></div><div style="font-style: normal; font-variant-caps: normal; font-weight: normal; letter-spacing: normal; text-align: start; text-indent: 0px; text-transform: none; white-space: normal; word-spacing: 0px; -webkit-text-stroke-width: 0px; text-decoration: none; font-family: Calibri, Arial, Helvetica, sans-serif; font-size: 12pt;" class="">When I try to run with the following options: -snes_mf_operator -ts_dt 1.0e-12 -ts_adapt_time_step_increase_delay 2 -snes_force_iteration -pc_fieldsplit_detect_coupling -pc_type fieldsplit -fieldsplit_0_pc_type jacobi -fieldsplit_1_pc_type redundant -ts_max_time 1000.0 -ts_adapt_dt_max 2.0e-3 -ts_adapt_wnormtype INFINITY -ts_exact_final_time stepover -ts_max_snes_failures -1 -ts_monitor -ts_max_steps 20</div><div style="font-style: normal; font-variant-caps: normal; font-weight: normal; letter-spacing: normal; text-align: start; text-indent: 0px; text-transform: none; white-space: normal; word-spacing: 0px; -webkit-text-stroke-width: 0px; text-decoration: none; font-family: Calibri, Arial, Helvetica, sans-serif; font-size: 12pt;" class=""><br class=""></div><div style="font-style: normal; font-variant-caps: normal; font-weight: normal; letter-spacing: normal; text-align: start; text-indent: 0px; text-transform: none; white-space: normal; word-spacing: 0px; -webkit-text-stroke-width: 0px; text-decoration: none; font-family: Calibri, Arial, Helvetica, sans-serif; font-size: 12pt;" class="">I get an error:</div><div style="font-style: normal; font-variant-caps: normal; font-weight: normal; letter-spacing: normal; text-align: start; text-indent: 0px; text-transform: none; white-space: normal; word-spacing: 0px; -webkit-text-stroke-width: 0px; text-decoration: none; font-family: Calibri, Arial, Helvetica, sans-serif; font-size: 12pt;" class="">[0]PETSC ERROR: --------------------- Error Message --------------------------------------------------------------<div class="">[0]PETSC ERROR: No support for this operation for this object type</div><div class="">[0]PETSC ERROR: Matrix type mffd does not have a find off block diagonal entries defined</div><div class="">[0]PETSC ERROR: See<span class="Apple-converted-space"> </span><a href="https://www.mcs.anl.gov/petsc/documentation/faq.html" class="">https://www.mcs.anl.gov/petsc/documentation/faq.html</a><span class="Apple-converted-space"> </span>for trouble shooting.</div><div class="">[0]PETSC ERROR: Petsc Development GIT revision: v3.13.4-851-gde18fec8da  GIT Date: 2020-08-28 16:47:50 +0000</div><div class="">[0]PETSC ERROR: Unknown Name on a 20200828 named sophie-Precision-5530 by sophie Tue Sep  1 10:58:44 2020</div><div class="">[0]PETSC ERROR: Configure options PETSC_DIR=/home/sophie/Code/petsc PETSC_ARCH=20200828 --with-cc=mpicc --with-cxx=mpicxx --with-fc=mpif77 --with-debugging=no --with-shared-libraries</div><div class="">[0]PETSC ERROR: #1 MatFindOffBlockDiagonalEntries() line 9847 in /home/sophie/Code/petsc/src/mat/interface/matrix.c</div><div class="">[0]PETSC ERROR: #2 PCFieldSplitSetDefaults() line 504 in /home/sophie/Code/petsc/src/ksp/pc/impls/fieldsplit/fieldsplit.c</div><div class="">[0]PETSC ERROR: #3 PCSetUp_FieldSplit() line 606 in /home/sophie/Code/petsc/src/ksp/pc/impls/fieldsplit/fieldsplit.c</div><div class="">[0]PETSC ERROR: #4 PCSetUp() line 1009 in /home/sophie/Code/petsc/src/ksp/pc/interface/precon.c</div><div class="">[0]PETSC ERROR: #5 KSPSetUp() line 406 in /home/sophie/Code/petsc/src/ksp/ksp/interface/itfunc.c</div><div class="">[0]PETSC ERROR: #6 KSPSolve_Private() line 658 in /home/sophie/Code/petsc/src/ksp/ksp/interface/itfunc.c</div><div class="">[0]PETSC ERROR: #7 KSPSolve() line 889 in /home/sophie/Code/petsc/src/ksp/ksp/interface/itfunc.c</div><div class="">[0]PETSC ERROR: #8 SNESSolve_NEWTONLS() line 225 in /home/sophie/Code/petsc/src/snes/impls/ls/ls.c</div><div class="">[0]PETSC ERROR: #9 SNESSolve() line 4524 in /home/sophie/Code/petsc/src/snes/interface/snes.c</div><div class="">[0]PETSC ERROR: #10 TSStep_ARKIMEX() line 811 in /home/sophie/Code/petsc/src/ts/impls/arkimex/arkimex.c</div><div class="">[0]PETSC ERROR: #11 TSStep() line 3731 in /home/sophie/Code/petsc/src/ts/interface/ts.c</div><div class="">[0]PETSC ERROR: #12 TSSolve() line 4128 in /home/sophie/Code/petsc/src/ts/interface/ts.c</div><div class="">PetscSolver::solve: TSSolve failed.</div><br class="">Cheers,</div><div style="font-style: normal; font-variant-caps: normal; font-weight: normal; letter-spacing: normal; text-align: start; text-indent: 0px; text-transform: none; white-space: normal; word-spacing: 0px; -webkit-text-stroke-width: 0px; text-decoration: none; font-family: Calibri, Arial, Helvetica, sans-serif; font-size: 12pt;" class=""><br class=""></div><div style="font-style: normal; font-variant-caps: normal; font-weight: normal; letter-spacing: normal; text-align: start; text-indent: 0px; text-transform: none; white-space: normal; word-spacing: 0px; -webkit-text-stroke-width: 0px; text-decoration: none; font-family: Calibri, Arial, Helvetica, sans-serif; font-size: 12pt;" class="">Sophie<br class=""></div><div id="appendonsend" style="caret-color: rgb(0, 0, 0); font-family: Helvetica; font-size: 18px; font-style: normal; font-variant-caps: normal; font-weight: normal; letter-spacing: normal; text-align: start; text-indent: 0px; text-transform: none; white-space: normal; word-spacing: 0px; -webkit-text-stroke-width: 0px; text-decoration: none;" class=""></div><hr tabindex="-1" style="caret-color: rgb(0, 0, 0); font-family: Helvetica; font-size: 18px; font-style: normal; font-variant-caps: normal; font-weight: normal; letter-spacing: normal; text-align: start; text-indent: 0px; text-transform: none; white-space: normal; word-spacing: 0px; -webkit-text-stroke-width: 0px; text-decoration: none; display: inline-block; width: 949.609375px;" class=""><span style="caret-color: rgb(0, 0, 0); font-family: Helvetica; font-size: 18px; font-style: normal; font-variant-caps: normal; font-weight: normal; letter-spacing: normal; text-align: start; text-indent: 0px; text-transform: none; white-space: normal; word-spacing: 0px; -webkit-text-stroke-width: 0px; text-decoration: none; float: none; display: inline !important;" class=""></span><div id="divRplyFwdMsg" dir="ltr" style="caret-color: rgb(0, 0, 0); font-family: Helvetica; font-size: 18px; font-style: normal; font-variant-caps: normal; font-weight: normal; letter-spacing: normal; text-align: start; text-indent: 0px; text-transform: none; white-space: normal; word-spacing: 0px; -webkit-text-stroke-width: 0px; text-decoration: none;" class=""><font face="Calibri, sans-serif" style="font-size: 11pt;" class=""><b class="">De :</b><span class="Apple-converted-space"> </span>Barry Smith <<a href="mailto:bsmith@petsc.dev" class="">bsmith@petsc.dev</a>><br class=""><b class="">Envoyé :</b><span class="Apple-converted-space"> </span>lundi 31 août 2020 14:50<br class=""><b class="">À :</b><span class="Apple-converted-space"> </span>Blondel, Sophie <<a href="mailto:sblondel@utk.edu" class="">sblondel@utk.edu</a>><br class=""><b class="">Cc :</b><span class="Apple-converted-space"> </span><a href="mailto:petsc-users@mcs.anl.gov" class="">petsc-users@mcs.anl.gov</a><span class="Apple-converted-space"> </span><<a href="mailto:petsc-users@mcs.anl.gov" class="">petsc-users@mcs.anl.gov</a>>;<span class="Apple-converted-space"> </span><a href="mailto:xolotl-psi-development@lists.sourceforge.net" class="">xolotl-psi-development@lists.sourceforge.net</a><span class="Apple-converted-space"> </span><<a href="mailto:xolotl-psi-development@lists.sourceforge.net" class="">xolotl-psi-development@lists.sourceforge.net</a>><br class=""><b class="">Objet :</b><span class="Apple-converted-space"> </span>Re: [petsc-users] Matrix Free Method questions</font><div class=""> </div></div><div class="" style="caret-color: rgb(0, 0, 0); font-family: Helvetica; font-size: 18px; font-style: normal; font-variant-caps: normal; font-weight: normal; letter-spacing: normal; text-align: start; text-indent: 0px; text-transform: none; white-space: normal; word-spacing: 0px; -webkit-text-stroke-width: 0px; text-decoration: none; word-wrap: break-word; line-break: after-white-space;"><div class=""><br class=""></div> Sophie,<div class=""><br class=""></div><div class="">   Thanks. </div><div class=""><br class=""></div><div class="">   The factor of 4 is lot, the 1.5 not so bad.</div><div class=""><br class=""></div><div class="">   You will definitely want to retain the full matrix assembly codes for speed and to verify a reduced matrix version.</div><div class=""><br class=""></div><div class="">   It is worth trying a "reduced matrix version" with matrix-free multiply based on these numbers. This reduced matrix Jacobian will only have the diagonals and all the terms connected to the cluster sizes that move. In other words you will be building just the part of the Jacobian needed for the new preconditioner (PC subtype for Jacobi) and doing the matrix-vector product matrix free. (SOR requires all the Jacobian entries).</div><div class=""><br class=""></div><div class="">   Fortunately this is hopefully pretty straightforward for this code. You will not have to change the structure of the main code at all.</div><div class=""><br class=""></div><div class="">  Step 1) create a new "sparse matrix" that will be passed to DMDASetBlockFillsSparse(). This new "sparse matrix" needs to retain all the diagonal entries and also all the entries that are associated with the variables that diffuse. If I remember correctly these are just the smallest cluster size, plain Helium?</div><div class=""><br class=""></div><div class="">  Call MatSetOptions(mat,MAT_NEW_NONZERO_LOCATIONS,PETSC_FALSE); </div><div class=""><br class=""></div><div class=""><div class="">Then you would run the code with -snes_mf_operator and the new PC subtype for Jacobi.</div><div class=""><br class=""></div><div class="">  A test that the new reduced Jacobian is correct will be that you get almost the same iterations as the runs you just make using the PC subtype of Jacobi. Hopefully not slower and using a great deal less memory. The iterations will not be identical because of the matrix-free multiple.</div></div><div class=""><br class=""></div><div class=""> Step 2) create a new version of the Jacobian computation routine. This routine should only compute the elements of the Jacobian needed for this reduced matrix Jacobian, so the diagonals and the diffusion/convection terms. </div><div class=""><br class=""></div><div class="">   Again run with with -snes_mf_operator and the new PC subtype for Jacobi and you should again get the same convergence history.</div><div class=""><br class=""></div><div class="">   I made two steps because it makes it easier to validate and debug to get the same results as before. The first step cheats in that it still computes the full Jacobian but ignores the entries that we don't need to store for the preconditioner. The second step is more efficient because it only computes the Jacobian entries needed for the preconditioner but it requires you going through the Jacobian code and making sure only the needed parts are computed.</div><div class=""><br class=""></div><div class="">  </div><div class="">  If you have any questions please let me know. </div><div class=""><br class=""></div><div class="">  Barry</div><div class=""><br class=""></div><div class=""><br class=""></div><div class=""><br class=""><div class=""><br class=""><blockquote type="cite" class=""><div class="">On Aug 31, 2020, at 1:13 PM, Blondel, Sophie <<a href="mailto:sblondel@utk.edu" class="">sblondel@utk.edu</a>> wrote:</div><br class="x_Apple-interchange-newline"><div class=""><div class="" style="font-style: normal; font-variant-caps: normal; font-weight: normal; letter-spacing: normal; text-align: start; text-indent: 0px; text-transform: none; white-space: normal; word-spacing: 0px; text-decoration: none; font-family: Calibri, Arial, Helvetica, sans-serif; font-size: 12pt;">Hi Barry,</div><div class="" style="font-style: normal; font-variant-caps: normal; font-weight: normal; letter-spacing: normal; text-align: start; text-indent: 0px; text-transform: none; white-space: normal; word-spacing: 0px; text-decoration: none; font-family: Calibri, Arial, Helvetica, sans-serif; font-size: 12pt;"><br class=""></div><div class="" style="font-style: normal; font-variant-caps: normal; font-weight: normal; letter-spacing: normal; text-align: start; text-indent: 0px; text-transform: none; white-space: normal; word-spacing: 0px; text-decoration: none; font-family: Calibri, Arial, Helvetica, sans-serif; font-size: 12pt;">I ran the 2 cases to look at the effect of the Jacobi pre-conditionner:</div><div class="" style="font-style: normal; font-variant-caps: normal; font-weight: normal; letter-spacing: normal; text-align: start; text-indent: 0px; text-transform: none; white-space: normal; word-spacing: 0px; text-decoration: none; font-family: Calibri, Arial, Helvetica, sans-serif; font-size: 12pt;"><ul class=""><li class="">1D with 200 grid points and 7759 DOF per grid point (for the PSI application), for 20 TS: the factor between SOR and Jacobi is ~4 (976 MatMult for SOR and 4162 MatMult for Jacobi)</li><li class="">2D with 63x63 grid points and 4124 DOF per grid point (for the NE application), for 20 TS: the factor is 1.5 (6657 for SOR, 10379 for Jacobi)</li></ul><div class="">Cheers,</div><div class=""><br class=""></div><div class="">Sophie<br class=""></div></div><div id="x_appendonsend" class="" style="font-family: Helvetica; font-size: 18px; font-style: normal; font-variant-caps: normal; font-weight: normal; letter-spacing: normal; text-align: start; text-indent: 0px; text-transform: none; white-space: normal; word-spacing: 0px; text-decoration: none;"></div><hr tabindex="-1" class="" style="font-family: Helvetica; font-size: 18px; font-style: normal; font-variant-caps: normal; font-weight: normal; letter-spacing: normal; text-align: start; text-indent: 0px; text-transform: none; white-space: normal; word-spacing: 0px; text-decoration: none; display: inline-block; width: 949.609375px;"><span class="" style="font-family: Helvetica; font-size: 18px; font-style: normal; font-variant-caps: normal; font-weight: normal; letter-spacing: normal; text-align: start; text-indent: 0px; text-transform: none; white-space: normal; word-spacing: 0px; text-decoration: none; float: none; display: inline !important;"></span><div id="x_divRplyFwdMsg" dir="ltr" class="" style="font-family: Helvetica; font-size: 18px; font-style: normal; font-variant-caps: normal; font-weight: normal; letter-spacing: normal; text-align: start; text-indent: 0px; text-transform: none; white-space: normal; word-spacing: 0px; text-decoration: none;"><font face="Calibri, sans-serif" class="" style="font-size: 11pt;"><b class="">De :</b><span class="x_Apple-converted-space"> </span>Barry Smith <<a href="mailto:bsmith@petsc.dev" class="">bsmith@petsc.dev</a>><br class=""><b class="">Envoyé :</b><span class="x_Apple-converted-space"> </span>vendredi 28 août 2020 18:31<br class=""><b class="">À :</b><span class="x_Apple-converted-space"> </span>Blondel, Sophie <<a href="mailto:sblondel@utk.edu" class="">sblondel@utk.edu</a>><br class=""><b class="">Cc :</b><span class="x_Apple-converted-space"> </span><a href="mailto:petsc-users@mcs.anl.gov" class="">petsc-users@mcs.anl.gov</a><span class="x_Apple-converted-space"> </span><<a href="mailto:petsc-users@mcs.anl.gov" class="">petsc-users@mcs.anl.gov</a>>;<span class="x_Apple-converted-space"> </span><a href="mailto:xolotl-psi-development@lists.sourceforge.net" class="">xolotl-psi-development@lists.sourceforge.net</a><span class="x_Apple-converted-space"> </span><<a href="mailto:xolotl-psi-development@lists.sourceforge.net" class="">xolotl-psi-development@lists.sourceforge.net</a>><br class=""><b class="">Objet :</b><span class="x_Apple-converted-space"> </span>Re: [petsc-users] Matrix Free Method questions</font><div class=""> </div></div><div class="" style="font-family: Helvetica; font-size: 18px; font-style: normal; font-variant-caps: normal; font-weight: normal; letter-spacing: normal; text-align: start; text-indent: 0px; text-transform: none; white-space: normal; word-spacing: 0px; text-decoration: none; word-wrap: break-word; line-break: after-white-space;"><br class=""><div class=""><br class=""><blockquote type="cite" class=""><div class="">On Aug 28, 2020, at 4:11 PM, Blondel, Sophie <<a href="mailto:sblondel@utk.edu" class="">sblondel@utk.edu</a>> wrote:</div><br class="x_x_Apple-interchange-newline"><div class=""><div class="" style="font-style: normal; font-variant-caps: normal; font-weight: normal; letter-spacing: normal; text-align: start; text-indent: 0px; text-transform: none; white-space: normal; word-spacing: 0px; text-decoration: none; font-family: Calibri, Arial, Helvetica, sans-serif; font-size: 12pt;">Thank you Jed and Barry,</div><div class="" style="font-style: normal; font-variant-caps: normal; font-weight: normal; letter-spacing: normal; text-align: start; text-indent: 0px; text-transform: none; white-space: normal; word-spacing: 0px; text-decoration: none; font-family: Calibri, Arial, Helvetica, sans-serif; font-size: 12pt;"><br class=""></div><div class="" style="font-style: normal; font-variant-caps: normal; font-weight: normal; letter-spacing: normal; text-align: start; text-indent: 0px; text-transform: none; white-space: normal; word-spacing: 0px; text-decoration: none; font-family: Calibri, Arial, Helvetica, sans-serif; font-size: 12pt;">First, attached are the logs from the benchmark runs I did without (log_std.txt) and with MF method (log_mf.txt). It took me some trouble to get the -log_view to work because I'm using push and pop for the options which means that PETSc is initialized with no argument so the command line argument was not taken into account, but I guess this is for a separate discussion.</div><div class="" style="font-style: normal; font-variant-caps: normal; font-weight: normal; letter-spacing: normal; text-align: start; text-indent: 0px; text-transform: none; white-space: normal; word-spacing: 0px; text-decoration: none; font-family: Calibri, Arial, Helvetica, sans-serif; font-size: 12pt;"><br class=""></div><div class="" style="font-style: normal; font-variant-caps: normal; font-weight: normal; letter-spacing: normal; text-align: start; text-indent: 0px; text-transform: none; white-space: normal; word-spacing: 0px; text-decoration: none; font-family: Calibri, Arial, Helvetica, sans-serif; font-size: 12pt;">To answer questions about the current per-conditioners:</div><div class="" style="font-style: normal; font-variant-caps: normal; font-weight: normal; letter-spacing: normal; text-align: start; text-indent: 0px; text-transform: none; white-space: normal; word-spacing: 0px; text-decoration: none; font-family: Calibri, Arial, Helvetica, sans-serif; font-size: 12pt;"><ul class=""><li class="">I used the same pre-conditioner options as listed in my previous email when I added the -snes_mf option; I did try to remove all the PC related options at one point with the MF method but didn't see a change in runtime so I put them back in</li><li class="">this benchmark is for a 1D DMDA using 20 grid points; when running in 2D or 3D I switch the PC options to: -pc_type fieldsplit -fieldsplit_0_pc_type sor -fieldsplit_1_pc_type gamg -fieldsplit_1_ksp_type gmres -ksp_type fgmres -fieldsplit_1_pc_gamg_threshold -1</li></ul><div class="">I haven't tried a Jacobi PC instead of SOR, I will run a set of more realistic runs (1D and 2D) without MF but with Jacobi and report on it next week. When you say "iterations" do you mean what is given by -ksp_monitor?</div></div></div></blockquote><div class=""><br class=""></div>  Yes, the number of MatMult is a good enough surrogate.</div><div class=""><br class=""></div><div class="">  So using matrix-free (which means no preconditioning) has </div><div class=""><br class=""></div><div class="">  35846/160<div class=""><br class=""></div><div class="">ans =</div><div class=""><br class=""></div><div class="">  224.0375</div><div class=""><br class=""></div><div class="">  or 224 as many iterations. So even for this modest 1d problem preconditioning is doing a great deal.</div><div class=""><br class=""></div><div class="">  Barry</div><div class=""><br class=""></div><div class=""><br class=""></div><div class=""><br class=""></div><blockquote type="cite" class=""><div class=""><div class="" style="font-style: normal; font-variant-caps: normal; font-weight: normal; letter-spacing: normal; text-align: start; text-indent: 0px; text-transform: none; white-space: normal; word-spacing: 0px; text-decoration: none; font-family: Calibri, Arial, Helvetica, sans-serif; font-size: 12pt;"><div class=""><br class=""></div><div class="">Cheers,</div><div class=""><br class=""></div><div class="">Sophie<br class=""></div></div><div id="x_x_appendonsend" class="" style="font-family: Helvetica; font-size: 18px; font-style: normal; font-variant-caps: normal; font-weight: normal; letter-spacing: normal; text-align: start; text-indent: 0px; text-transform: none; white-space: normal; word-spacing: 0px; text-decoration: none;"></div><hr tabindex="-1" class="" style="font-family: Helvetica; font-size: 18px; font-style: normal; font-variant-caps: normal; font-weight: normal; letter-spacing: normal; text-align: start; text-indent: 0px; text-transform: none; white-space: normal; word-spacing: 0px; text-decoration: none; display: inline-block; width: 949.609375px;"><span class="" style="font-family: Helvetica; font-size: 18px; font-style: normal; font-variant-caps: normal; font-weight: normal; letter-spacing: normal; text-align: start; text-indent: 0px; text-transform: none; white-space: normal; word-spacing: 0px; text-decoration: none; float: none; display: inline !important;"></span><div id="x_x_divRplyFwdMsg" dir="ltr" class="" style="font-family: Helvetica; font-size: 18px; font-style: normal; font-variant-caps: normal; font-weight: normal; letter-spacing: normal; text-align: start; text-indent: 0px; text-transform: none; white-space: normal; word-spacing: 0px; text-decoration: none;"><font face="Calibri, sans-serif" class="" style="font-size: 11pt;"><b class="">De :</b><span class="x_x_Apple-converted-space"> </span>Barry Smith <<a href="mailto:bsmith@petsc.dev" class="">bsmith@petsc.dev</a>><br class=""><b class="">Envoyé :</b><span class="x_x_Apple-converted-space"> </span>vendredi 28 août 2020 12:12<br class=""><b class="">À :</b><span class="x_x_Apple-converted-space"> </span>Blondel, Sophie <<a href="mailto:sblondel@utk.edu" class="">sblondel@utk.edu</a>><br class=""><b class="">Cc :</b><span class="x_x_Apple-converted-space"> </span><a href="mailto:petsc-users@mcs.anl.gov" class="">petsc-users@mcs.anl.gov</a><span class="x_x_Apple-converted-space"> </span><<a href="mailto:petsc-users@mcs.anl.gov" class="">petsc-users@mcs.anl.gov</a>>;<span class="x_x_Apple-converted-space"> </span><a href="mailto:xolotl-psi-development@lists.sourceforge.net" class="">xolotl-psi-development@lists.sourceforge.net</a><span class="x_x_Apple-converted-space"> </span><<a href="mailto:xolotl-psi-development@lists.sourceforge.net" class="">xolotl-psi-development@lists.sourceforge.net</a>><br class=""><b class="">Objet :</b><span class="x_x_Apple-converted-space"> </span>Re: [petsc-users] Matrix Free Method questions</font><div class=""> </div></div><div class="" style="font-family: Helvetica; font-size: 18px; font-style: normal; font-variant-caps: normal; font-weight: normal; letter-spacing: normal; text-align: start; text-indent: 0px; text-transform: none; white-space: normal; word-spacing: 0px; text-decoration: none; word-wrap: break-word; line-break: after-white-space;"><div class="" style="margin-top: 0px; margin-bottom: 0px;"><span class="" style="color: rgb(51, 51, 51); background-color: rgb(252, 235, 166); border: 10px solid rgb(252, 235, 166); font-family: Arial, sans-serif; margin-top: 12px;"><strong class="">[External Email]</strong></span></div><div class=""><div class=""><br class=""></div>  Sophie,<div class=""><br class=""></div><div class="">   This is exactly what i would expect. If you run with -ksp_monitor you will see the -snes_mf run takes many more iterations.</div><div class=""><br class=""></div><div class="">   I am puzzled that the argument <span class="" style="font-family: Calibri, Arial, Helvetica, sans-serif; font-size: 12pt;">-pc_type fieldsplit did not stop the run since this is under normal circumstances not a viable preconditioner with -snes_mf. Did you also remove the </span><span class="" style="font-family: Calibri, Arial, Helvetica, sans-serif; font-size: 12pt;">-pc_type fieldsplit</span><span class="" style="font-family: Calibri, Arial, Helvetica, sans-serif; font-size: 12pt;"> argument?</span></div><div class=""><span class="" style="font-family: Calibri, Arial, Helvetica, sans-serif; font-size: 12pt;"><br class=""></span></div><div class=""><span class="" style="font-family: Calibri, Arial, Helvetica, sans-serif; font-size: 12pt;">   In order to see how one can avoid forming the entire matrix and use matrix-free to do the matrix-vector but still have an effective preconditioner let's look at what the current preconditioner options do.</span></div><div class=""><span class="" style="font-family: Calibri, Arial, Helvetica, sans-serif; font-size: 12pt;"><br class=""></span></div><div class=""><blockquote type="cite" class=""><div class="" style="font-family: Calibri, Arial, Helvetica, sans-serif; font-size: 12pt;"> -pc_fieldsplit_detect_coupling </div></blockquote><div class=""><br class=""></div>creates two sub-preconditioners, the first for all the variables and the second for those that are coupled by the matrix to variables in neighboring cells Since only the smallest cluster sizes have diffusion/advection this second set contains only the cluster size one variables.</div><div class=""><br class=""></div><div class=""><blockquote type="cite" class=""><div class="" style="font-family: Calibri, Arial, Helvetica, sans-serif; font-size: 12pt;">-fieldsplit_0_pc_type sor </div></blockquote><br class=""></div><div class="">Runs SOR on all the variables; you can think of this as running SOR on the reactions, it is a pretty good preconditioner for the reactions since the reactions are local, per cell.</div><div class=""><br class=""></div><div class=""><blockquote type="cite" class=""><div class="" style="font-family: Calibri, Arial, Helvetica, sans-serif; font-size: 12pt;"> -fieldsplit_1_pc_type redundant</div></blockquote></div><div class=""><font face="Calibri, Arial, Helvetica, sans-serif" size="3" class="">  <br class=""></font><div class="">This runs the default preconditioner (ILU) on just the variables that diffuse, i.e. the elliptic part. For smallish problems this is fine, for larger problems and 2d and 3d presumably you have also -redundant_pc_type gamg to use algebraic multigrid for the diffusion.  This part of the matrix will always need to be formed and used in the preconditioner. It  is very important since the diffusion is what brings in most of the ill-conditioning for larger problems into the linear system. Note that it only needs the matrix entries for the cluster size of 1 so it is very small compared to the entire sparse matrix.</div><div class=""><br class=""></div><div class="">----</div><div class=""> The first preconditioner SOR requires ALL the matrix entries which are almost all (except for the diffusion terms) the coupling between different size clusters within a cell. Especially each cell has its own sparse matrix of the size of total number of clusters, it is sparse but not super sparse.</div><div class=""><br class=""></div><div class=""> So the to significantly lower memory usage we need to remove the SOR and the storing of all the matrix entries but still have an efficient preconditioner for the "reaction" terms. </div><div class=""><br class=""></div><div class=""> The simplest thing would be to use Jacobi instead of SOR for the first subpreconditioner since it only requires the diagonal entries in the matrix. But Jacobi is a worse preconditioner than SOR (since it totally ignores the matrix coupling) and sometimes can be much worse.</div><div class=""><br class=""></div><div class="">  Before anyone writes additional code we need to know if doing something along these lines does not ruin the convergence that.</div><div class=""><br class=""></div><div class=""> Have you used the same options as before but with  <span class="" style="font-family: Calibri, Arial, Helvetica, sans-serif; font-size: 12pt;">-fieldsplit_0_pc_type jacobi ? (Not using any matrix free). We need to get an idea of how many more linear iterations it requires (not time, comparing time won't be helpful for this exercise.) We also need this information for realistic size problems in 2 or 3 dimensions that you really want to run; for small problems this approach will work ok and give misleading information about what happens for large problems.</span></div><div class=""><span class="" style="font-family: Calibri, Arial, Helvetica, sans-serif; font-size: 12pt;"><br class=""></span></div><div class=""><font face="Calibri, Arial, Helvetica, sans-serif" size="3" class="">  I suspect the iteration counts will shot up. Can you run some cases and see how the iteration counts change?</font></div><div class=""><font face="Calibri, Arial, Helvetica, sans-serif" size="3" class=""><br class=""></font></div><div class=""><font face="Calibri, Arial, Helvetica, sans-serif" size="3" class="">  Based on that we can decide if we still retain "good convergence" by changing the SOR to Jacobi and then change the code to make this change efficient (basically by skipping the explicit computation of the reaction Jacobian terms and using matrix-free on the outside of the PCFIELDSPLIT.)</font></div><div class=""><font face="Calibri, Arial, Helvetica, sans-serif" size="3" class=""><br class=""></font></div><div class=""><font face="Calibri, Arial, Helvetica, sans-serif" size="3" class="">  Barry</font></div><div class=""><font face="Calibri, Arial, Helvetica, sans-serif" size="3" class=""><br class=""></font></div><div class=""><br class=""></div><div class="">  </div><div class=""><br class=""></div><div class="">  </div><div class=""><br class=""></div><div class=""><br class=""></div><div class=""><br class=""></div><div class=""><br class=""><blockquote type="cite" class=""><div class="">On Aug 28, 2020, at 9:49 AM, Blondel, Sophie via petsc-users <<a href="mailto:petsc-users@mcs.anl.gov" class="">petsc-users@mcs.anl.gov</a>> wrote:</div><br class="x_x_x_Apple-interchange-newline"><div class=""><div class="" style="font-style: normal; font-variant-caps: normal; font-weight: normal; letter-spacing: normal; text-align: start; text-indent: 0px; text-transform: none; white-space: normal; word-spacing: 0px; text-decoration: none; font-family: Calibri, Arial, Helvetica, sans-serif; font-size: 12pt;">Hi everyone,</div><div class="" style="font-style: normal; font-variant-caps: normal; font-weight: normal; letter-spacing: normal; text-align: start; text-indent: 0px; text-transform: none; white-space: normal; word-spacing: 0px; text-decoration: none; font-family: Calibri, Arial, Helvetica, sans-serif; font-size: 12pt;"><br class=""></div><div class="" style="font-style: normal; font-variant-caps: normal; font-weight: normal; letter-spacing: normal; text-align: start; text-indent: 0px; text-transform: none; white-space: normal; word-spacing: 0px; text-decoration: none; font-family: Calibri, Arial, Helvetica, sans-serif; font-size: 12pt;">I have been using PETSc for a few years with a fully implicit TS ARKIMEX method and am now exploring the matrix free method option. Here is the list of PETSc options I typically use: -ts_dt 1.0e-12 -ts_adapt_time_step_increase_delay 5 -snes_force_iteration -ts_max_time 1000.0 -ts_adapt_dt_max 2.0e-3 -ts_adapt_wnormtype INFINITY -ts_exact_final_time stepover -fieldsplit_0_pc_type sor -ts_max_snes_failures -1 -pc_fieldsplit_detect_coupling -ts_monitor -pc_type fieldsplit -fieldsplit_1_pc_type redundant -ts_max_steps 100</div><div class="" style="font-style: normal; font-variant-caps: normal; font-weight: normal; letter-spacing: normal; text-align: start; text-indent: 0px; text-transform: none; white-space: normal; word-spacing: 0px; text-decoration: none; font-family: Calibri, Arial, Helvetica, sans-serif; font-size: 12pt;"><br class=""></div><div class="" style="font-style: normal; font-variant-caps: normal; font-weight: normal; letter-spacing: normal; text-align: start; text-indent: 0px; text-transform: none; white-space: normal; word-spacing: 0px; text-decoration: none; font-family: Calibri, Arial, Helvetica, sans-serif; font-size: 12pt;">I started to compare the performance of the code without changing anything of the executable and simply adding "-snes_mf", I see a reduction of memory usage as expected and a benchmark that would usually take ~5min to run now takes ~50min. Reading the documentation I saw that there are a few option to play with the matrix free method like -snes_mf_err, -snes_mf_umin, or switching to -snes_mf_type wp. I used and modified the values of each of these options separately but never saw a sizable change in runtime, is it expected?</div><div class="" style="font-style: normal; font-variant-caps: normal; font-weight: normal; letter-spacing: normal; text-align: start; text-indent: 0px; text-transform: none; white-space: normal; word-spacing: 0px; text-decoration: none; font-family: Calibri, Arial, Helvetica, sans-serif; font-size: 12pt;"><br class=""></div><div class="" style="font-style: normal; font-variant-caps: normal; font-weight: normal; letter-spacing: normal; text-align: start; text-indent: 0px; text-transform: none; white-space: normal; word-spacing: 0px; text-decoration: none; font-family: Calibri, Arial, Helvetica, sans-serif; font-size: 12pt;">And are there other ways to make the matrix free method faster? I saw in the documentation that you can define your own per-conditioner for instance. Let me know if you need additional information about the PETSc setup in the application I use.<br class=""></div><div class="" style="font-style: normal; font-variant-caps: normal; font-weight: normal; letter-spacing: normal; text-align: start; text-indent: 0px; text-transform: none; white-space: normal; word-spacing: 0px; text-decoration: none; font-family: Calibri, Arial, Helvetica, sans-serif; font-size: 12pt;"><br class=""></div><div class="" style="font-style: normal; font-variant-caps: normal; font-weight: normal; letter-spacing: normal; text-align: start; text-indent: 0px; text-transform: none; white-space: normal; word-spacing: 0px; text-decoration: none; font-family: Calibri, Arial, Helvetica, sans-serif; font-size: 12pt;">Best,</div><div class="" style="font-style: normal; font-variant-caps: normal; font-weight: normal; letter-spacing: normal; text-align: start; text-indent: 0px; text-transform: none; white-space: normal; word-spacing: 0px; text-decoration: none; font-family: Calibri, Arial, Helvetica, sans-serif; font-size: 12pt;"><br class=""></div><div class="" style="font-style: normal; font-variant-caps: normal; font-weight: normal; letter-spacing: normal; text-align: start; text-indent: 0px; text-transform: none; white-space: normal; word-spacing: 0px; text-decoration: none; font-family: Calibri, Arial, Helvetica, sans-serif; font-size: 12pt;">Sophie</div></div></blockquote></div><br class=""></div></div></div><span id="x_x_cid:A09F8670-6BEA-4EAE-9B4F-81D5A43993BF" class=""><log_mf.txt></span><span id="x_x_cid:F0448F58-40FD-4152-81DA-D1C986F6A242" class=""><log_std.txt></span></div></blockquote></div></div></div></blockquote></div></div></div></div></blockquote></div><br class=""></div></body></html>