[CIG-SHORT] Linear convergence for nonplanar faults

Brad Aagaard baagaard at usgs.gov
Tue Nov 8 16:17:59 PST 2016


Things I would check:

1. Without a fault, check convergence with using ML for the preconditioner.

2. Add the fault and with zero prescribed slip (FaultCohesiveKin), 
examine the convergence with the solver_fault_fieldsplit.cfg parameters.

3. Use fault friction (FaultCohesiveDyn) but under conditions when you 
don't have slip. You should see similar solver behavior to 2. I would 
set the ksp_atol parameter to 2 orders of magnitude smaller than 
zero_tolerance and snes_atol to 2 order of magnitude greater than 
zero_tolerance.

4. If 1-3 look good, then run under conditions when you expect slip. The 
linear convergence should be similar, but you will need multiple 
nonlinear iteractions to converge.

Regards,
Brad



On 11/08/2016 02:47 PM, Josimar Alves da Silva wrote:
> Dear Charles and Brad,
>
> Thank you for following up on my question. I implemented the changes you
> suggested and it did not improve the convergence. I am attaching the log
> file for you reference.
>
> All the other parameters of the previous files were kept unchanged,
> except the length_scale = 500*m.
>
> Let me know if you have any other suggestion.
> Thank you
> Josimar
>
> On Tue, Nov 8, 2016 at 4:46 PM, Brad Aagaard <baagaard at usgs.gov
> <mailto:baagaard at usgs.gov>> wrote:
>
>     Josimar,
>
>     Your solver settings in pylithapp.cfg are obsolete. To use
>     fieldsplit and the custom preconditioner, please use the solver
>     settings in share/settings/solver_fault_fieldsplit.cfg (shown
>     below). Note that instead of 0/1 for the field, we now use names
>     (displacement/lagrange_multiplier).
>
>     These settings will only work well if there is a fault. Without a
>     fault, you can just use ML as the preconditioner. See
>     examples/3d/hex8/step02.cfg (just use pc_type = ml).
>
>     Regards,
>     Brad
>
>
>     [pylithapp.problem.formulation]
>     split_fields = True
>     matrix_type = aij
>     use_custom_constraint_pc = True
>
>     [pylithapp.petsc]
>     fs_pc_type = fieldsplit
>     fs_pc_use_amat =
>     fs_pc_fieldsplit_type = multiplicative
>     fs_fieldsplit_displacement_pc_type = ml
>     fs_fieldsplit_lagrange_multiplier_pc_type = jacobi
>     fs_fieldsplit_displacement_ksp_type = preonly
>     fs_fieldsplit_lagrange_multiplier_ksp_type = preonly
>
>
>
>
>     On 11/08/2016 11:00 AM, Josimar Alves da Silva wrote:
>
>         Dear Brad, Matt and Charles,
>
>         I am running Pylith for a 2D friction model with gravity for a
>         nonplanar
>         fault. I am finding that the linear solver does not converge
>         when the
>         fault is included on the model. When the fault is not included
>         on the
>         model the linear solver takes roughly 500 to 1000 iterations to
>         converge. Please refer to the attached "log_no_fault.cfg" file.
>
>         The three links below show previous discussions on linear
>         convergence:
>
>         http://lists.geodynamics.org/pipermail/cig-short/2012-September/001123.html
>         <http://lists.geodynamics.org/pipermail/cig-short/2012-September/001123.html>
>
>
>         http://lists.geodynamics.org/pipermail/cig-short/2008-December/000453.html
>         <http://lists.geodynamics.org/pipermail/cig-short/2008-December/000453.html>
>
>
>         http://lists.geodynamics.org/pipermail/cig-short/2013-March/001320.html
>         <http://lists.geodynamics.org/pipermail/cig-short/2013-March/001320.html>
>
>
>
>         As pointed out by Charles previously and also discussed on the links
>         above, overall, it seems that there are 3 main causes for
>         convergence
>         issues: 1) Poor quality mesh; 2) incorrect solver settings; 3) Poor
>         setup/initial conditions.
>
>
>         For my set up, I am unable to make the linear solver to converge.
>
>
>         This is what I have done so far, based on recommendations from
>         previous
>         discussions:
>
>
>         1) Make sure that the mesh has condition number and aspect
>         ratios <1.5.
>         Please see attached presentation for histograms.
>
>         2) used the solver parameters shown
>         on share/settings/solver_fault_fieldsplit.cfg
>
>
>         I would like to know if there is anything else that I can do in
>         order to
>         make the linear solver to converge. Any other advices to improve my
>         simulation would be very much appreciated.
>
>
>         You can find attached the pylith .cfg files that I am using and
>         also the
>         log files for the case without the fault (log_no_fault.out) and
>         with the
>         fault (log_with_fault.out).
>
>
>         I run my simulation using:
>
>
>         1) pylith no_fault.cfg > log_nor_fault.out
>
>         2) pylith no_fault.cfg   fault.cfg > log_with_fault.out
>
>
>
>         Thank you in advance for you help,
>
>         Josimar
>
>
>
>
>         _______________________________________________
>         CIG-SHORT mailing list
>         CIG-SHORT at geodynamics.org <mailto:CIG-SHORT at geodynamics.org>
>         http://lists.geodynamics.org/cgi-bin/mailman/listinfo/cig-short
>         <http://lists.geodynamics.org/cgi-bin/mailman/listinfo/cig-short>
>
>
>     _______________________________________________
>     CIG-SHORT mailing list
>     CIG-SHORT at geodynamics.org <mailto:CIG-SHORT at geodynamics.org>
>     http://lists.geodynamics.org/cgi-bin/mailman/listinfo/cig-short
>     <http://lists.geodynamics.org/cgi-bin/mailman/listinfo/cig-short>
>
>
>
>
> _______________________________________________
> CIG-SHORT mailing list
> CIG-SHORT at geodynamics.org
> http://lists.geodynamics.org/cgi-bin/mailman/listinfo/cig-short
>



More information about the CIG-SHORT mailing list