<div dir="ltr">And don't use -pc_gamg_parallel_coarse_grid_solver<div>You can use that in production but for debugging use -mg_coarse_pc_type svd</div><div>Also, use -options_left and remove anything that is not used.</div><div>(I am puzzled, I see -pc_type gamg not -pc_type fieldsplit)</div><div><br></div><div>Mark</div><div><br></div></div><br><div class="gmail_quote gmail_quote_container"><div dir="ltr" class="gmail_attr">On Mon, Jun 16, 2025 at 6:40 AM Matthew Knepley <<a href="mailto:knepley@gmail.com">knepley@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"><div dir="ltr">On Sun, Jun 15, 2025 at 9:46 PM hexioafeng <<a href="mailto:hexiaofeng@buaa.edu.cn" target="_blank">hexiaofeng@buaa.edu.cn</a>> wrote:</div><div class="gmail_quote"><blockquote class="gmail_quote" style="margin:0px 0px 0px 0.8ex;border-left:1px solid rgb(204,204,204);padding-left:1ex"><div><div><font face="Verdana" size="4">Hello,</font></div><div><font face="Verdana" size="4"><br></font></div><div><font face="Verdana" size="4">Here are the options and outputs:</font></div><div><font face="Verdana" size="4"><br></font></div><div><div><font face="Verdana" size="4">options:</font></div><div><font face="Verdana" size="4"><br></font></div><div><font face="Verdana" size="4">-ksp_type cg -pc_type gamg -pc_gamg_parallel_coarse_grid_solver  -pc_fieldsplit_detect_saddle_point  -pc_fieldsplit_type schur -pc_fieldsplit_schur_precondition selfp -fieldsplit_1_mat_schur_complement_ainv_type lump -pc_fieldsplit_schur_fact_type full -fieldsplit_0_ksp_type preonly -fieldsplit_0_pc_type gamg -fieldsplit_0_mg_coarse_pc_type_type svd -fieldsplit_1_ksp_type preonly -fieldsplit_1_pc_type bjacobi -fieldsplit_1_sub_pc_type sor -ksp_view  -ksp_monitor_true_residual  -ksp_converged_reason  -fieldsplit_0_mg_levels_ksp_monitor_true_residual  -fieldsplit_0_mg_levels_ksp_converged_reason  -fieldsplit_1_ksp_monitor_true_residual  -fieldsplit_1_ksp_converged_reason </font></div></div></div></blockquote><div><br></div><div>This option was wrong:</div><div><br></div><div>  -fieldsplit_0_mg_coarse_pc_type_type svd</div><div><br></div><div>from the output, we can see that it should have been</div><div><br></div><div>  -fieldsplit_0_mg_coarse_sub_pc_type_type svd</div><div><br></div><div>  THanks,</div><div><br></div><div>     Matt</div><div> </div><blockquote class="gmail_quote" style="margin:0px 0px 0px 0.8ex;border-left:1px solid rgb(204,204,204);padding-left:1ex"><div><div><div><font face="Verdana" size="4">output:</font></div><div><font face="Verdana" size="4"><br></font></div><div><font face="Verdana" size="4">0 KSP unpreconditioned resid norm 2.777777777778e+01 true resid norm 2.777777777778e+01 ||r(i)||/||b|| 1.000000000000e+00</font></div><div><font face="Verdana" size="4">Linear solve did not converge due to DIVERGED_PC_FAILED iterations 0</font></div><div><font face="Verdana" size="4">               PC failed due to SUBPC_ERROR</font></div><div><font face="Verdana" size="4">KSP Object: 1 MPI processes</font></div><div><font face="Verdana" size="4">  type: cg</font></div><div><font face="Verdana" size="4">  maximum iterations=200, initial guess is zero</font></div><div><font face="Verdana" size="4">  tolerances:  relative=1e-06, absolute=1e-12, divergence=1e+30</font></div><div><font face="Verdana" size="4">  left preconditioning</font></div><div><font face="Verdana" size="4">  using UNPRECONDITIONED norm type for convergence test</font></div><div><font face="Verdana" size="4">PC Object: 1 MPI processes</font></div><div><font face="Verdana" size="4">  type: gamg</font></div><div><font face="Verdana" size="4">    type is MULTIPLICATIVE, levels=2 cycles=v</font></div><div><font face="Verdana" size="4">      Cycles per PCApply=1</font></div><div><font face="Verdana" size="4">      Using externally compute Galerkin coarse grid matrices</font></div><div><font face="Verdana" size="4">      GAMG specific options</font></div><div><font face="Verdana" size="4">        Threshold for dropping small values in graph on each level =</font></div><div><font face="Verdana" size="4">        Threshold scaling factor for each level not specified = 1.</font></div><div><font face="Verdana" size="4">        AGG specific options</font></div><div><font face="Verdana" size="4">          Symmetric graph false</font></div><div><font face="Verdana" size="4">          Number of levels to square graph 1</font></div><div><font face="Verdana" size="4">          Number smoothing steps 1</font></div><div><font face="Verdana" size="4">        Complexity:    grid = 1.00176</font></div><div><font face="Verdana" size="4">  Coarse grid solver -- level -------------------------------</font></div><div><font face="Verdana" size="4">    KSP Object: (mg_coarse_) 1 MPI processes</font></div><div><font face="Verdana" size="4">      type: preonly</font></div><div><font face="Verdana" size="4">      maximum iterations=10000, initial guess is zero</font></div><div><font face="Verdana" size="4">      tolerances:  relative=1e-05, absolute=1e-50, divergence=10000.</font></div><div><font face="Verdana" size="4">      left preconditioning</font></div><div><font face="Verdana" size="4">      using NONE norm type for convergence test</font></div><div><font face="Verdana" size="4">    PC Object: (mg_coarse_) 1 MPI processes</font></div><div><font face="Verdana" size="4">      type: bjacobi</font></div><div><font face="Verdana" size="4">        number of blocks = 1</font></div><div><font face="Verdana" size="4">        Local solver is the same for all blocks, as in the following KSP and PC objects on rank 0:</font></div><div><font face="Verdana" size="4">      KSP Object: (mg_coarse_sub_) 1 MPI processes</font></div><div><font face="Verdana" size="4">        type: preonly</font></div><div><font face="Verdana" size="4">        maximum iterations=1, initial guess is zero</font></div><div><font face="Verdana" size="4">        tolerances:  relative=1e-05, absolute=1e-50, divergence=10000.</font></div><div><font face="Verdana" size="4">        left preconditioning</font></div><div><font face="Verdana" size="4">        using NONE norm type for convergence test</font></div><div><font face="Verdana" size="4">      PC Object: (mg_coarse_sub_) 1 MPI processes</font></div><div><font face="Verdana" size="4">        type: lu</font></div><div><font face="Verdana" size="4">          out-of-place factorization</font></div><div><font face="Verdana" size="4">          tolerance for zero pivot 2.22045e-14</font></div><div><font face="Verdana" size="4">          using diagonal shift on blocks to prevent zero pivot [INBLOCKS]</font></div><div><font face="Verdana" size="4">          matrix ordering: nd</font></div><div><font face="Verdana" size="4">          factor fill ratio given 5., needed 1.</font></div><div><font face="Verdana" size="4">            Factored matrix follows:</font></div><div><font face="Verdana" size="4">              Mat Object: 1 MPI processes</font></div><div><font face="Verdana" size="4">                type: seqaij</font></div><div><font face="Verdana" size="4">                rows=7, cols=7</font></div><div><font face="Verdana" size="4">                package used to perform factorization: petsc</font></div><div><font face="Verdana" size="4">                total: nonzeros=45, allocated nonzeros=45</font></div><div><font face="Verdana" size="4">                  using I-node routines: found 3 nodes, limit used is 5</font></div><div><font face="Verdana" size="4">        linear system matrix = precond matrix:</font></div><div><font face="Verdana" size="4">        Mat Object: 1 MPI processes</font></div><div><font face="Verdana" size="4">          type: seqaij</font></div><div><font face="Verdana" size="4">          rows=7, cols=7</font></div><div><font face="Verdana" size="4">          total: nonzeros=45, allocated nonzeros=45</font></div><div><font face="Verdana" size="4">          total number of mallocs used during MatSetValues calls=0</font></div><div><font face="Verdana" size="4">            using I-node routines: found 3 nodes, limit used is 5</font></div><div><font face="Verdana" size="4">      linear system matrix = precond matrix:</font></div><div><font face="Verdana" size="4">      Mat Object: 1 MPI processes</font></div><div><font face="Verdana" size="4">        type: mpiaij</font></div><div><font face="Verdana" size="4">        rows=7, cols=7</font></div><div><font face="Verdana" size="4">        total: nonzeros=45, allocated nonzeros=45</font></div><div><font face="Verdana" size="4">        total number of mallocs used during MatSetValues calls=0</font></div><div><font face="Verdana" size="4">          using nonscalable MatPtAP() implementation</font></div><div><font face="Verdana" size="4">          using I-node (on process 0) routines: found 3 nodes, limit used is 5</font></div><div><font face="Verdana" size="4">  Down solver (pre-smoother) on level 1 -------------------------------</font></div><div><font face="Verdana" size="4">    KSP Object: (mg_levels_1_) 1 MPI processes</font></div><div><font face="Verdana" size="4">      type: chebyshev</font></div><div><font face="Verdana" size="4">        eigenvalue estimates used:  min = 0., max = 0.</font></div><div><font face="Verdana" size="4">        eigenvalues estimate via gmres min 0., max 0.</font></div><div><font face="Verdana" size="4">        eigenvalues estimated using gmres with translations  [0. 0.1; 0. 1.1]</font></div><div><font face="Verdana" size="4">        KSP Object: (mg_levels_1_esteig_) 1 MPI processes</font></div><div><font face="Verdana" size="4">          type: gmres</font></div><div><font face="Verdana" size="4">            restart=30, using Classical (unmodified) Gram-Schmidt Orthogonalization with no iterative refinement</font></div><div><font face="Verdana" size="4">            happy breakdown tolerance 1e-30</font></div><div><font face="Verdana" size="4">          maximum iterations=10, initial guess is zero</font></div><div><font face="Verdana" size="4">          tolerances:  relative=1e-12, absolute=1e-50, divergence=10000.</font></div><div><font face="Verdana" size="4">          left preconditioning</font></div><div><font face="Verdana" size="4">          using PRECONDITIONED norm type for convergence test</font></div><div><font face="Verdana" size="4">        PC Object: (mg_levels_1_) 1 MPI processes</font></div><div><font face="Verdana" size="4">          type: sor</font></div><div><font face="Verdana" size="4">            type = local_symmetric, iterations = 1, local iterations = 1, omega = 1.</font></div><div><font face="Verdana" size="4">          linear system matrix = precond matrix:</font></div><div><font face="Verdana" size="4">          Mat Object: 1 MPI processes</font></div><div><font face="Verdana" size="4">            type: mpiaij</font></div><div><font face="Verdana" size="4">            rows=624, cols=624</font></div><div><font face="Verdana" size="4">            total: nonzeros=25536, allocated nonzeros=25536</font></div><div><font face="Verdana" size="4">            total number of mallocs used during MatSetValues calls=0</font></div><div><font face="Verdana" size="4">              using I-node (on process 0) routines: found 336 nodes, limit used is 5</font></div><div><font face="Verdana" size="4">        estimating eigenvalues using noisy right hand side</font></div><div><font face="Verdana" size="4">      maximum iterations=2, nonzero initial guess</font></div><div><font face="Verdana" size="4">      tolerances:  relative=1e-05, absolute=1e-50, divergence=10000.</font></div><div><font face="Verdana" size="4">      left preconditioning</font></div><div><font face="Verdana" size="4">      using NONE norm type for convergence test</font></div><div><font face="Verdana" size="4">    PC Object: (mg_levels_1_) 1 MPI processes</font></div><div><font face="Verdana" size="4">      type: sor</font></div><div><font face="Verdana" size="4">        type = local_symmetric, iterations = 1, local iterations = 1, omega = 1.      linear system matrix = precond matrix:</font></div><div><font face="Verdana" size="4">      Mat Object: 1 MPI processes</font></div><div><font face="Verdana" size="4">        type: mpiaij</font></div><div><font face="Verdana" size="4">        rows=624, cols=624</font></div><div><font face="Verdana" size="4">        total: nonzeros=25536, allocated nonzeros=25536</font></div><div><font face="Verdana" size="4">        total number of mallocs used during MatSetValues calls=0</font></div><div><font face="Verdana" size="4">          using I-node (on process 0) routines: found 336 nodes, limit used is 5  Up solver (post-smoother) same as down solver (pre-smoother)</font></div><div><font face="Verdana" size="4">  linear system matrix = precond matrix:</font></div><div><font face="Verdana" size="4">  Mat Object: 1 MPI processes</font></div><div><font face="Verdana" size="4">    type: mpiaij</font></div><div><font face="Verdana" size="4">    rows=624, cols=624</font></div><div><font face="Verdana" size="4">    total: nonzeros=25536, allocated nonzeros=25536</font></div><div><font face="Verdana" size="4">    total number of mallocs used during MatSetValues calls=0</font></div><div><font face="Verdana" size="4">      using I-node (on process 0) routines: found 336 nodes, limit used is 5</font></div></div><div><font face="Verdana" size="4"><br></font></div><div><font face="Verdana" size="4"><br></font></div><div><font face="Verdana" size="4">Best regards,</font></div><div><font face="Verdana" size="4"><br></font></div><div><font face="Verdana" size="4">Xiaofeng</font></div><br id="m_2839857572054762463m_-3003253161800872436lineBreakAtBeginningOfMessage"><div><br><blockquote type="cite"><div>On Jun 14, 2025, at 07:28, Barry Smith <<a href="mailto:bsmith@petsc.dev" target="_blank">bsmith@petsc.dev</a>> wrote:</div><br><div><div><div><br></div>  Matt,<div><br></div><div>   Perhaps we should add options -ksp_monitor_debug and -snes_monitor_debug that turn on all possible monitoring for the (possibly) nested solvers and all of their converged reasons also? Note this is not completely trivial because each preconditioner will have to supply its list based on the current solver options for it.</div><div><br></div><div>   Then we won't need to constantly list a big string of problem specific monitor options to ask the user to use.</div><div><br></div><div>  Barry</div><div><br></div><div><br></div><div><br id="m_2839857572054762463m_-3003253161800872436lineBreakAtBeginningOfMessage"><div><br><blockquote type="cite"><div>On Jun 13, 2025, at 9:09 AM, Matthew Knepley <<a href="mailto:knepley@gmail.com" target="_blank">knepley@gmail.com</a>> wrote:</div><br><div><div dir="ltr"><div dir="ltr">On Thu, Jun 12, 2025 at 10:55 PM hexioafeng <<a href="mailto:hexiaofeng@buaa.edu.cn" target="_blank">hexiaofeng@buaa.edu.cn</a>> wrote:</div><div class="gmail_quote"><blockquote class="gmail_quote" style="margin:0px 0px 0px 0.8ex;border-left:1px solid rgb(204,204,204);padding-left:1ex"><div><div><font face="Verdana" size="4">Dear authors,</font></div><div><font face="Verdana" size="4"><br></font></div><div><font face="Verdana" size="4">I tried <b>-pc_type game -pc_gamg_parallel_coarse_grid_solver</b> and <b>-pc_type field split -pc_fieldsplit_detect_saddle_point -fieldsplit_0_ksp_type pronely -fieldsplit_0_pc_type game -fieldsplit_0_mg_coarse_pc_type sad -fieldsplit_1_ksp_type pronely -fieldsplit_1_pc_type Jacobi _fieldsplit_1_sub_pc_type for</b> , both options got the KSP_DIVERGE_PC_FAILED error.</font></div></div></blockquote><div><br></div><div>With any question about convergence, we need to see the output of</div><div><br></div><div>  -ksp_view -ksp_monitor_true_residual -ksp_converged_reason -fieldsplit_0_mg_levels_ksp_monitor_true_residual -fieldsplit_0_mg_levels_ksp_converged_reason -fieldsplit_1_ksp_monitor_true_residual -fieldsplit_1_ksp_converged_reason</div><div><br></div><div>and all the error output.</div><div><br></div><div>  Thanks,</div><div><br></div><div>     Matt</div><div> </div><blockquote class="gmail_quote" style="margin:0px 0px 0px 0.8ex;border-left:1px solid rgb(204,204,204);padding-left:1ex"><div><div><font face="Verdana" size="4">Thanks,</font></div><div><font face="Verdana" size="4"><br></font></div><div><font face="Verdana" size="4">Xiaofeng</font></div><br id="m_2839857572054762463m_-3003253161800872436m_-8146634516906487715lineBreakAtBeginningOfMessage"><div><br><blockquote type="cite"><div>On Jun 12, 2025, at 20:50, Mark Adams <<a href="mailto:mfadams@lbl.gov" target="_blank">mfadams@lbl.gov</a>> wrote:</div><br><div><br><br style="font-family:Helvetica;font-size:12px;font-style:normal;font-variant-caps:normal;font-weight:400;letter-spacing:normal;text-align:start;text-indent:0px;text-transform:none;white-space:normal;word-spacing:0px;text-decoration:none"><div class="gmail_quote" style="font-family:Helvetica;font-size:12px;font-style:normal;font-variant-caps:normal;font-weight:400;letter-spacing:normal;text-align:start;text-indent:0px;text-transform:none;white-space:normal;word-spacing:0px;text-decoration:none"><div dir="ltr" class="gmail_attr">On Thu, Jun 12, 2025 at 8:44 AM Matthew Knepley <<a href="mailto:knepley@gmail.com" target="_blank">knepley@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"><div dir="ltr">On Thu, Jun 12, 2025 at 4:58 AM Mark Adams <<a href="mailto:mfadams@lbl.gov" target="_blank">mfadams@lbl.gov</a>> wrote:</div><div class="gmail_quote"><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>Adding this to the PETSc mailing list,</div><br><div class="gmail_quote"><div dir="ltr" class="gmail_attr">On Thu, Jun 12, 2025 at 3:43 AM hexioafeng <<a href="mailto:hexiaofeng@buaa.edu.cn" target="_blank">hexiaofeng@buaa.edu.cn</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><div><font size="4" face="Verdana"><br></font></div><div><font size="4" face="Verdana">Dear Professor,</font></div><div><font size="4" face="Verdana"><br></font></div><div><font size="4" face="Verdana">I hope this message finds you well.</font></div><div><font size="4" face="Verdana"><br></font></div><div><font size="4" face="Verdana">I am an employee at a CAE company and a heavy user of the PETSc library. I would like to thank you for your contributions to PETSc and express my deep appreciation for your work.</font></div><div><font size="4" face="Verdana"><br></font></div><div><font size="4" face="Verdana">Recently, I encountered some difficulties when using PETSc to solve structural mechanics problems with Lagrange multiplier constraints. After searching extensively online and reviewing several papers, I found your previous paper titled "<b>Algebraic multigrid methods for constrained linear systems with applications to contact problems in solid mechanics</b>" seems to be the most relevant and helpful. </font></div><div><font size="4" face="Verdana"><br></font></div><div><font size="4" face="Verdana">The stiffness matrix I'm working with,<span> </span><b>K</b>, is a block saddle-point matrix of the form (A00 A01; A10 0), where<span> </span><b>A00 is singular</b>—just as described in your paper, and different from many other articles . I have a few questions regarding your work and would greatly appreciate your insights:</font></div><div><font size="4" face="Verdana"><br></font></div><div><font size="4" face="Verdana">1. Is the<span> </span><b>AMG/KKT</b><span> </span>method presented in your paper available in PETSc? I tried using<span> </span><b>CG+GAMG</b><span> </span>directly but received a<span> </span><b>KSP_DIVERGED_PC_FAILED</b><span> </span>error. I also attempted to use<span> </span><b>CG+PCFIELDSPLIT</b><span> </span>with the following options:  </font></div></div></blockquote><div><br></div><div>No</div><div> </div><blockquote class="gmail_quote" style="margin:0px 0px 0px 0.8ex;border-left:1px solid rgb(204,204,204);padding-left:1ex"><div><div><font size="4" face="Verdana">   </font></div><div><font size="4" face="Verdana">   <span> </span>-pc_type fieldsplit -pc_fieldsplit_detect_saddle_point -pc_fieldsplit_type schur -pc_fieldsplit_schur_precondition selfp -pc_fieldsplit_schur_fact_type full -fieldsplit_0_ksp_type preonly -fieldsplit_0_pc_type gamg -fieldsplit_1_ksp_type preonly -fieldsplit_1_pc_type bjacobi  </font></div><div><font size="4" face="Verdana">   </font></div><div><font size="4" face="Verdana">   Unfortunately, this also resulted in a<span> </span><b>KSP_DIVERGED_PC_FAILED</b><span> </span>error. Do you have any suggestions?</font></div><div><font size="4" face="Verdana"><br></font></div><div><font size="4" face="Verdana">2. In your paper, you compare the method with<span> </span><b>Uzawa</b>-type approaches. To my understanding, Uzawa methods typically require A00 to be invertible. How did you handle the singularity of A00 to construct an M-matrix that is invertible?</font></div><div><font size="4" face="Verdana"><br></font></div></div></blockquote><div><br></div>You add a regularization term like A01 * A10 (like springs). See the paper or any reference to augmented lagrange or Uzawa</div><div class="gmail_quote"><br><div><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><div><font size="4" face="Verdana"></font></div><div><font size="4" face="Verdana">3. Can i implement the AMG/KKT method in your paper using existing<span> </span><b>AMG APIs</b>? Implementing a production-level AMG solver from scratch would be quite challenging for me, so I’m hoping to utilize existing AMG interfaces within PETSc or other packages.</font></div><div><font size="4" face="Verdana"><br></font></div></div></blockquote><div><br></div><div>You can do Uzawa and make the regularization matrix with matrix-matrix products. Just use AMG for the A00 block.</div><div><br></div><div> </div><blockquote class="gmail_quote" style="margin:0px 0px 0px 0.8ex;border-left:1px solid rgb(204,204,204);padding-left:1ex"><div><div><font size="4" face="Verdana"></font></div><div><font size="4" face="Verdana">4. For saddle-point systems where A00 is singular, can you recommend any more robust or efficient solutions? Alternatively, are you aware of any open-source software packages that can handle such cases out-of-the-box?</font></div><div><font size="4" face="Verdana"><br></font></div></div></blockquote><div><br></div>No, and I don't think PETSc can do this out-of-the-box, but others may be able to give you a better idea of what PETSc can do.</div><div class="gmail_quote">I think PETSc can do Uzawa or other similar algorithms but it will not do the regularization automatically (it is a bit more complicated than just A01 * A10)</div></div></blockquote><div><br></div><div>One other trick you can use is to have</div><div><br></div><div> <span> </span>-fieldsplit_0_mg_coarse_pc_type svd</div><div><br></div><div>This will use SVD on the coarse grid of GAMG, which can handle the null space in A00 as long as the prolongation does not put it back in. I have used this for the Laplacian with Neumann conditions and for freely floating elastic problems.</div><div><br></div></div></div></blockquote><div><br></div><div>Good point.</div><div>You can also use -pc_gamg_parallel_coarse_grid_solver to get GAMG to use a on level iterative solver for the coarse grid.</div><div> </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 class="gmail_quote"><div></div><div> <span> </span>Thanks,</div><div><br></div><div>     Matt</div><div> </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 class="gmail_quote">Thanks,</div><div class="gmail_quote">Mark<br><blockquote class="gmail_quote" style="margin:0px 0px 0px 0.8ex;border-left:1px solid rgb(204,204,204);padding-left:1ex"><div><div><font size="4" face="Verdana"></font></div><div><font size="4" face="Verdana"><br></font></div><div><font size="4" face="Verdana">Thank you very much for taking the time to read my email. Looking forward to hearing from you.</font></div><div><font size="4" face="Verdana"><br></font></div><div><font size="4" face="Verdana"><br></font></div><div><font size="4" face="Verdana"><br></font></div><div><font size="4" face="Verdana">Sincerely,  </font></div><div><font size="4" face="Verdana"><br></font></div><div><font size="4" face="Verdana">Xiaofeng He</font></div><div><font size="4" face="Verdana">-----------------------------------------------------</font></div><div><font size="4" face="Verdana"><br></font></div><div><font size="4" face="Verdana">Research Engineer</font></div><div><font size="4" face="Verdana"><br></font></div><div><font size="4" face="Verdana">Internet Based Engineering, Beijing, China</font></div><div><div><br></div></div></div></blockquote></div></div></blockquote></div><div><br clear="all"></div><div><br></div><span class="gmail_signature_prefix">--<span> </span></span><br><div dir="ltr" class="gmail_signature"><div dir="ltr"><div><div dir="ltr"><div><div dir="ltr"><div>What most experimenters take for granted before they begin their experiments is infinitely more interesting than any results to which their experiments lead.<br>-- Norbert Wiener</div><div><br></div><div><a href="https://urldefense.us/v3/__http://www.cse.buffalo.edu/*knepley/__;fg!!G_uCfscf7eWS!f-YJSzthRa7atIa1xs1GPHW53hGIqSenvp1eO2kDsSyf4jv1_Vp0kL9Lg8pyyPeG8al4Im8XlLqGRHw1FxYh$" target="_blank">https://www.cse.buffalo.edu/~knepley/</a></div></div></div></div></div></div></div></div></blockquote></div></div></blockquote></div><br></div></blockquote></div><div><br clear="all"></div><div><br></div><span class="gmail_signature_prefix">-- </span><br><div dir="ltr" class="gmail_signature"><div dir="ltr"><div><div dir="ltr"><div><div dir="ltr"><div>What most experimenters take for granted before they begin their experiments is infinitely more interesting than any results to which their experiments lead.<br>-- Norbert Wiener</div><div><br></div><div><a href="https://urldefense.us/v3/__http://www.cse.buffalo.edu/*knepley/__;fg!!G_uCfscf7eWS!f-YJSzthRa7atIa1xs1GPHW53hGIqSenvp1eO2kDsSyf4jv1_Vp0kL9Lg8pyyPeG8al4Im8XlLqGRHw1FxYh$" target="_blank">https://www.cse.buffalo.edu/~knepley/</a><br></div></div></div></div></div></div></div></div>
</div></blockquote></div><br></div></div></div></blockquote></div><br></div></blockquote></div><div><br clear="all"></div><div><br></div><span class="gmail_signature_prefix">-- </span><br><div dir="ltr" class="gmail_signature"><div dir="ltr"><div><div dir="ltr"><div><div dir="ltr"><div>What most experimenters take for granted before they begin their experiments is infinitely more interesting than any results to which their experiments lead.<br>-- Norbert Wiener</div><div><br></div><div><a href="https://urldefense.us/v3/__http://www.cse.buffalo.edu/*knepley/__;fg!!G_uCfscf7eWS!dYETsi-moODALE1tmLrk5pxFKF9l552nNiC0cBgsCQ9ebugJWHtsNYa0QBS5Gmws9J_VC_Iec3Nx0c1FgNl1$" target="_blank">https://www.cse.buffalo.edu/~knepley/</a><br></div></div></div></div></div></div></div></div>
</blockquote></div>