***linear_solver dd_feti
#
Description#
This keyword specifies that the FETI domain decomposition method is used to solve the linear system of equations involved in the global step of the Newton-Raphson algorithm. FETI stands for Finite Element Tearing and Interconnecting. Zset has several implementations of the FETI method, this keyword refers to the lastest implementation made by C. Bovet, A. Parret-Fréaud and P. Gosselet.
Syntax#
This solver accepts the following commands.
[**solver_config
predefined_config ]
[**precond
type ]
[**scaling
type ]
[**projector_schur
type ]
[**projector_scaling
type ]
[**local_solver
local_solver ]
[**iterative_solver
iterative_solver ]
[**acceleration
acceleration ]
**solver_config
activates some predefined configuration for the preconditionner (see
precond
andscaling
) and for the projector (seeprojector_schur
andprojector_scaling
). It can take these four valuessecure
,moderate
,fast
,aggressive
. The following table gives the sets. In most cases,secure
is the most robust choice. In term of iterations, it converges faster thanmoderate
, which converges faster thanfast
which in turn is faster thanaggressive
. However, in terms of computational time, the performance may be problem dependent as the cost per iteration is decreasing with respect to the previous ordeer. Default value issecure
.Configuration
precond
scaling
projector_schur
projector_scaling
secure
full
stiffness
full
stiffness
moderate
full
stiffness
none
topological
fast
lumped
stiffness
sup
stiffness
aggressive
lumped
stiffness
none
none
Please note that, if
aggressive
is chosen, the CG algorithm exploits the CG short recurrence for orthogonalization. It may lead to a loss of search directions orthogonality caused by numerical round-off error.**precond
specifies the operator used during preconditioning. type can take these three values
full
(i.e. dirichlet preconditioning \(\hat{\bf S} = \bf S = \bf K_{bb} - \bf K_{bi} \bf K_{ii}^{-1} \bf K_{ib}\)),lumped
(i.e. lumped preconditioning \(\hat{\bf S} = \bf K_{bb}\)),super_lumped
(i.e. super lumped preconditioning \(\hat{\bf S} = \operatorname{diag}\left(\bf K_{bb}\right)\)). Those three values are sorted in descending order regarding the both the efficiency and the computational cost of the preconditionning step. Thus, the choice is a balance between the total number of iteration needed to reach convergence and the unitary cost of an interation. Default value isfull
.**scaling
specifies the scaling used during preconditioning. type can take these three values
none
(i.e. no scaling),topological
(i.e. scaling based on multiplicity of interface nodes),stiffness
(i.e. scaling based on the diagonal stiffness matrix terms of interfaces nodes). Stiffness scaling is strongly advised whenever the computation involves material heterogeneities. Default value isstiffness
.**projector_schur
specifies the operator used to build the coarse space projector. type can take four values. The first one is
none
(i.e. standard orthogonal projector). The three others values arefull
,lumped
andsuper_lumped
. These three values produce an oblique projector. The meaning of these values is the same than in the**precond
keyword. Default value isfull
.**projector_scaling
specifies the scaling used to build the projector. It can take the same values than the
scaling
keyword. Default value isstiffness
.**local_solver
specifies the solver used to solve problems defined in each subdomain. The following direct solvers are availables:
mumps
anddissection
. The wrappersrigid_geometric
andadvanced_detection
, which handle kernel detection through a geometric approach, can also be used. As the proper kernel detection and treatment of the floating subdomains operator is a crucial aspect in the convergence of FETI methods, proper choice of the local solver is mandatory. The wrappersrigid_geometric
andadvanced_detection
are known to be more**iterative_solver
specifies the iterative solver used to solve the global problem defined at the interface. robust regarding kernel detection and therefore should be preferred in problem involving multi material or complex geometries. Options can be passed tolocal_solver
in a similar manner as in its standalone invocation, but please note that you have to add the keywordno_option
as argument tomumps
ordissection
if you want to use it without extra options (e.g.**local_solver mumps no_option
). The default value is set torigid_geometric
.**iterative_solver
specifies the iterative solver used to solve the global problem defined at the interface. The following values are available whithin
dd_feti
context:ddcg
. For symmetric positive definite problems, you may consider useddcg
, which is the implementation of a conjugate gradient algorithm. Available iterative solvers names and commands are detailed in the**iterative_solver
section on . Default value isddcg
.**acceleration
allows to reuse search space across Newton steps and/or iterations, it may be useful for ill conditioned nonlinear simulations. Please note that this feature is still experimental. Available acceleration names and commands are detailed in the
**acceleration
section
Example#
The example below shows the use of the FETI solver associated with the conjugate gradient iterative solver. Once again, there is default values for all keywords but you may be interested in modifying specific values.
***linear_solver dd_feti
**solver_config moderate
**iterative_solver ddcg
*precision 1.e-10
*output_every_iter 10
*max_iteration 1000
**local_solver mumps no_option
The example below shows the use the rigid_geometric
wrapper for null
space detection. It also shows detailed configuration.
***linear_solver dd_mpfeti
**precond full
**scaling topological
**projector_schur sumper_lumped
**projector_scaling topological
**iterative_solver ddcg
*precision 1.e-10
*output_every_iter 10
*max_iteration 1000
*double_orthogonalization
**local_solver rigid_geometric
*local_solver mumps no_option