Thanks Jens, we'll try this next week.


On 02/01/2014 04:06 AM, Jens Lohne Eftang wrote:
> Hi,
>
> I used these arguments for 3D elasticity and they seemed to work well
> also in combination with DirichletBoundary (why does that break symmetry?)
>
> -ksp_type cg
> -pc_type fieldsplit
> -pc_fieldsplit_block_size 3
> -fieldsplit_pc_type hypre
> -fieldsplit_pc_hypre_type boomeramg
> -fieldsplit_pc_hypre_boomeramg_strong_threshold 0.7
> -pc_fieldsplit_0 0,1,2
> -ksp_atol 1e-11
> -ksp_rtol 1e-11
> -pc_fieldsplit_type symmetric_multiplicative
>
> It is key to use the fieldsplit preconditioner in order to be able to
> exploit smoothness within each field component. If it's 2D elasticity I
> believe you'd have to modify to -pc_fieldsplit_block_size 2 and
> -pc_fieldsplit_0 0,1.
>
> I was able to solve a problem with tens of millions dofs in 90 CG
> iterations using these parameters.
>
> Best,
> Jens
>
> On 01/31/2014 11:11 PM, David Knezevic wrote:
>> On 01/31/2014 12:36 PM, John Peterson wrote:
>>>> Hello,
>>>>
>>>> I'm trying to solve a linear elasticity problem based on a tetraheadral
>>>> mesh with over a million degrees of freedom.
>>>>
>>>> I've been using the LU solver in MUMPS but I'm now forced to switch to an
>>>> iterative solver because of the memory requirements of that solver.
>>>>
>>>> CG and BoomerAMG in hypre seem to be a good combination, but I've not been
>>>> able to achieve convergence. The residual never goes below 1 even after
>>>> several thousand iterations and shows high oscillations.
>>>>
>>>> The parameters I'm using are
>>>> -ksp_type cg
>>>> -ksp_monitor
>>>> -ksp_converged_reason
>>>> -pc_type hypre
>>>> -pc_hypre_type boomeramg
>>>> -pc_hypre_boomeramg_strong_threshold 0.7
>>>> -log_summary
>>>>
>>> First thing I would think of is that the problem isn't actually symmetric
>>> even though logically it should be?
>>>
>>> Does it work any better with -ksp_type bcgs?
>>>
>>>
>>>
>>>> -eps_monitor
>>>> -st_type sinvert
>>>>
>>> Is this an eigenvalue problem too?
>> Dana is working with me on this. It's not an eigenproblem, so
>> "-eps_monitor -st_type sinvert" isn't meant to be there (but it those
>> options would've just been ignored I believe).
>>
>> Regarding symmetry: That's a good thought. The PDE is linear elasticity,
>> hence it's symmetric. But we enforce boundary conditions using
>> DirchletBoundary, so that'll ruin symmetry on the Dirichlet rows/cols.
>> Not sure if that's the problem here or not, but Dana will try "-ksp_type
>> bcgs" to see if that helps convergence.
>>
>> Thanks,
>> David
>>
>>
>>
>>
>> ------------------------------------------------------------------------------
>> WatchGuard Dimension instantly turns raw network data into actionable
>> security intelligence. It gives you real-time visual feedback on key
>> security issues and trends.  Skip the complicated setup - simply import
>> a virtual appliance and go from zero to informed in seconds.
>> http://pubads.g.doubleclick.net/gampad/clk?id=123612991&iu=/4140/ostg.clktrk
>> _______________________________________________
>> Libmesh-users mailing list
>> [email protected]
>> https://lists.sourceforge.net/lists/listinfo/libmesh-users
>
> ------------------------------------------------------------------------------
> WatchGuard Dimension instantly turns raw network data into actionable
> security intelligence. It gives you real-time visual feedback on key
> security issues and trends.  Skip the complicated setup - simply import
> a virtual appliance and go from zero to informed in seconds.
> http://pubads.g.doubleclick.net/gampad/clk?id=123612991&iu=/4140/ostg.clktrk
> _______________________________________________
> Libmesh-users mailing list
> [email protected]
> https://lists.sourceforge.net/lists/listinfo/libmesh-users


------------------------------------------------------------------------------
WatchGuard Dimension instantly turns raw network data into actionable 
security intelligence. It gives you real-time visual feedback on key
security issues and trends.  Skip the complicated setup - simply import
a virtual appliance and go from zero to informed in seconds.
http://pubads.g.doubleclick.net/gampad/clk?id=123612991&iu=/4140/ostg.clktrk
_______________________________________________
Libmesh-users mailing list
[email protected]
https://lists.sourceforge.net/lists/listinfo/libmesh-users

Reply via email to