Predefined solutions¶
-
class
solution.predefined_solutions.
BucklingAnalysisEigenPart
(prb, staticAnalysisPart, name=None, printFlag=0, soeType: str = None, solverType: str = None, shift: float = 0.0)¶ Bases:
solution.predefined_solutions.SolutionProcedure
Eigenvalue part of a linear buckling analysis.
-
setup
()¶ Defines the solution procedure in the finite element problem object.
-
-
class
solution.predefined_solutions.
BucklingAnalysisStaticPart
(prb, name=None, constraintHandlerType='transformation', maxNumIter=100, convergenceTestTol=1e-09, printFlag=0, numSteps=1, numberingMethod='rcm', convTestType=None, soeType: str = None, solverType: str = None, solnAlgorithmType='krylov_newton_soln_algo')¶ Bases:
solution.predefined_solutions.SolutionProcedure
Static analysis part of a linear buckling analysis.
-
setup
()¶ Defines the solution procedure in the finite element problem object.
-
-
class
solution.predefined_solutions.
FrequencyAnalysis
(prb, name=None, printFlag=0, systemPrefix='sym_band', numberingMethod='rcm', shift: float = None)¶ Bases:
solution.predefined_solutions.SolutionProcedure
Return a natural frequency computation procedure.
-
setup
()¶ Defines the solution procedure in the finite element problem object.
-
-
class
solution.predefined_solutions.
IllConditioningAnalysis
(prb, name=None, printFlag=0, systemPrefix='band_arpack', shift: float = 0.0)¶ Bases:
solution.predefined_solutions.IllConditioningAnalysisBase
Procedure to obtain the modes associated with very small eigenvalues of the stiffness matrix.
-
class
solution.predefined_solutions.
IllConditioningAnalysisBase
(prb, name=None, printFlag=0, systemPrefix='sym_band_eigen', shift: float = None, numberingMethod='rcm')¶ Bases:
solution.predefined_solutions.SolutionProcedure
Base class for ill-conditioning solution procedures.
-
setup
()¶ Defines the solution procedure in the finite element problem object.
-
-
class
solution.predefined_solutions.
LineSearchBase
(prb, name, constraintHandlerType, maxNumIter, convergenceTestTol, printFlag, numSteps, numberingMethod, convTestType, lineSearchMethod, soeType, solverType)¶ Bases:
solution.predefined_solutions.SolutionProcedure
Base class for line search solution aggregations.
-
class
solution.predefined_solutions.
LinearBucklingAnalysis
(prb, numModes, constraintHandlerType='transformation', numberingMethod='rcm', convTestType='norm_disp_incr_conv_test', convergenceTestTol=1e-08, maxNumIter=1000, soeType='band_gen_lin_soe', solverType='band_gen_lin_lapack_solver', solnAlgorithmType='krylov_newton_soln_algo', eigenSOEType='band_arpackpp_soe', eigenSolverType='band_arpackpp_solver')¶ Bases:
object
Linear buckling analysis.
Variables: - staticPart – solution for the static problem.
- eigenPart – solution for the eigen problem.
- numModes – number of buckling modes to compute.
-
analysisSetup
()¶ Define the analysis object.
Parameters: analysisType – type of the analysis to perform.
-
setup
()¶ Defines the solution procedure in the finite element problem object.
-
soluControlSetup
()¶ Defines the solution control object.
-
solve
()¶ Compute the linear buckling modes.
-
class
solution.predefined_solutions.
PenaltyModifiedNewton
(prb, name=None, maxNumIter=150, convergenceTestTol=1e-09, printFlag=0, numSteps=1, numberingMethod='rcm', convTestType='relative_total_norm_disp_incr_conv_test')¶ Bases:
solution.predefined_solutions.PenaltyModifiedNewtonBase
Static solution procedure with a modified Newton algorithm and a penalty constraint handler.
-
setup
()¶ Defines the solution procedure in the finite element problem object.
-
-
class
solution.predefined_solutions.
PenaltyModifiedNewtonBase
(prb, name, maxNumIter, convergenceTestTol, printFlag, numSteps, numberingMethod, convTestType='relative_total_norm_disp_incr_conv_test', soeType='sparse_gen_col_lin_soe', solverType='super_lu_solver')¶ Bases:
solution.predefined_solutions.SolutionProcedure
Base class for penalty modified Newton solution aggregation.
-
setup
()¶ Defines the solution procedure in the finite element problem object.
-
-
class
solution.predefined_solutions.
PenaltyModifiedNewtonMUMPS
(prb, name=None, maxNumIter=150, convergenceTestTol=1e-09, printFlag=0, numSteps=1, numberingMethod='rcm', convTestType='relative_total_norm_disp_incr_conv_test')¶ Bases:
solution.predefined_solutions.PenaltyModifiedNewtonBase
Static solution procedure with a modified Newton algorithm, a penalty constraint handler and a MUMPS (parallel sparse direct solver) solver.
-
setup
()¶ Defines the solution procedure in the finite element problem object.
-
-
class
solution.predefined_solutions.
PenaltyModifiedNewtonUMF
(prb, name=None, maxNumIter=150, convergenceTestTol=1e-09, printFlag=0, numSteps=1, numberingMethod='rcm', convTestType='relative_total_norm_disp_incr_conv_test')¶ Bases:
solution.predefined_solutions.PenaltyModifiedNewtonBase
Static solution procedure with a modified Newton algorithm, a penalty constraint handler and a UMF (Unsimmetric multi-frontal method) solver.
-
setup
()¶ Defines the solution procedure in the finite element problem object.
-
-
class
solution.predefined_solutions.
PenaltyNewmarkNewtonRapshon
(prb, name=None, maxNumIter=10, convergenceTestTol=1e-09, printFlag=0, numSteps=1, numberingMethod='rcm', convTestType='norm_disp_incr_conv_test')¶ Bases:
solution.predefined_solutions.SolutionProcedure
Newmark solution procedure with a Newton Raphson algorithm and a penalty constraint handler.
-
setup
()¶ Defines the solution procedure in the finite element problem object.
-
-
class
solution.predefined_solutions.
PenaltyNewtonLineSearch
(prb, name=None, maxNumIter=150, convergenceTestTol=1e-09, printFlag=0, numSteps=1, numberingMethod='rcm', convTestType='relative_total_norm_disp_incr_conv_test', lineSearchMethod='regula_falsi_line_search')¶ Bases:
solution.predefined_solutions.PenaltyNewtonLineSearchBase
Static solution procedure with a Newton line search algorithm and a penalty constraint handler.
-
setup
()¶ Defines the solution procedure in the finite element problem object.
-
-
class
solution.predefined_solutions.
PenaltyNewtonLineSearchBase
(prb, name, maxNumIter, convergenceTestTol, printFlag, numSteps, numberingMethod, convTestType, lineSearchMethod, soeType, solverType)¶ Bases:
solution.predefined_solutions.LineSearchBase
Base class for penalty Newton line search solution aggregation.
-
setup
()¶ Defines the solution procedure in the finite element problem object.
-
-
class
solution.predefined_solutions.
PenaltyNewtonLineSearchMUMPS
(prb, name=None, maxNumIter=150, convergenceTestTol=1e-09, printFlag=0, numSteps=1, numberingMethod='rcm', convTestType='relative_total_norm_disp_incr_conv_test', lineSearchMethod='regula_falsi_line_search')¶ Bases:
solution.predefined_solutions.PenaltyNewtonLineSearchBase
Static solution procedure with a Newton line search algorithm, a penalty constraint handler and a MUMPS (parallel sparse direct solver) solver.
-
setup
()¶ Defines the solution procedure in the finite element problem object.
-
-
class
solution.predefined_solutions.
PenaltyNewtonLineSearchUMF
(prb, name=None, maxNumIter=150, convergenceTestTol=1e-09, printFlag=0, numSteps=1, numberingMethod='rcm', convTestType='relative_total_norm_disp_incr_conv_test', lineSearchMethod='regula_falsi_line_search')¶ Bases:
solution.predefined_solutions.PenaltyNewtonLineSearchBase
Static solution procedure with a Newton line search algorithm, a penalty constraint handler and a UMF (Unsimmetric multi-frontal method) solver.
-
setup
()¶ Defines the solution procedure in the finite element problem object.
-
-
class
solution.predefined_solutions.
PenaltyNewtonRaphson
(prb, name=None, maxNumIter=10, convergenceTestTol=1e-09, printFlag=0, numSteps=1, numberingMethod='rcm', convTestType='norm_unbalance_conv_test')¶ Bases:
solution.predefined_solutions.PenaltyNewtonRaphsonBase
Return a static solution procedure with a Newton Raphson algorithm and a penalty constraint handler.
-
setup
()¶ Defines the solution procedure in the finite element problem object.
-
-
class
solution.predefined_solutions.
PenaltyNewtonRaphsonBase
(prb, name, maxNumIter, convergenceTestTol, printFlag, numSteps, numberingMethod, convTestType, soeType, solverType)¶ Bases:
solution.predefined_solutions.SolutionProcedure
Base class for penalty modified Newton solution aggregation.
-
setup
()¶ Defines the solution procedure in the finite element problem object.
-
-
class
solution.predefined_solutions.
PenaltyNewtonRaphsonMUMPS
(prb, name=None, maxNumIter=150, convergenceTestTol=1e-09, printFlag=0, numSteps=1, numberingMethod='rcm', convTestType='norm_unbalance_conv_test')¶ Bases:
solution.predefined_solutions.PenaltyNewtonRaphsonBase
Static solution procedure with a Newton algorithm, a penalty constraint handler and a MUMPS (parallel sparse direct solver) solver.
-
setup
()¶ Defines the solution procedure in the finite element problem object.
-
-
class
solution.predefined_solutions.
PenaltyNewtonRaphsonUMF
(prb, name=None, maxNumIter=150, convergenceTestTol=1e-09, printFlag=0, numSteps=1, numberingMethod='rcm', convTestType='norm_unbalance_conv_test')¶ Bases:
solution.predefined_solutions.PenaltyNewtonRaphsonBase
Static solution procedure with a Newton algorithm, a penalty constraint handler and a UMF (Unsimmetric multi-frontal method) solver.
-
setup
()¶ Defines the solution procedure in the finite element problem object.
-
-
class
solution.predefined_solutions.
PenaltyStaticLinearBase
(prb, name=None, printFlag=0, numSteps=1, numberingMethod='rcm', soeType='band_spd_lin_soe', solverType='band_spd_lin_lapack_solver')¶ Bases:
solution.predefined_solutions.SolutionProcedure
Base class for linear static solution algorithm with a penalty constraint handler.
-
setup
()¶ Defines the solution procedure in the finite element problem object.
-
-
class
solution.predefined_solutions.
PlainKrylovNewton
(prb, name=None, maxNumIter=150, convergenceTestTol=1e-09, printFlag=0, numSteps=1, numberingMethod='simple', convTestType='energy_incr_conv_test', maxDim=6)¶ Bases:
solution.predefined_solutions.SolutionProcedure
- KrylovNewton algorithm object which uses a Krylov subspace
- accelerator to accelerate the convergence of the modified newton method.
See appendix C. Analysis model script of the document “Finite Element Modeling of Gusset Plate Failure Using Opensees” Andrew J. Walker. Oregon State University
-
setup
()¶ Defines the solution procedure in the finite element problem object.
-
class
solution.predefined_solutions.
PlainLinearNewmark
(prb, timeStep, name=None, maxNumIter=10, convergenceTestTol=1e-09, printFlag=0, numSteps=1, numberingMethod='simple')¶ Bases:
solution.predefined_solutions.SolutionProcedure
Return a linear Newmark solution algorithm with a plain constraint handler.
-
setup
()¶ Defines the solution procedure in the finite element problem object.
-
solve
(calculateNodalReactions=False, includeInertia=False, reactionCheckTolerance=1e-12)¶ Compute the solution (run the analysis).
Parameters: - calculateNodalReactions – if true calculate reactions at nodes.
- includeInertia – if true calculate reactions including inertia effects.
- reactionCheckTolerance – tolerance when checking reaction values.
-
-
class
solution.predefined_solutions.
PlainNewtonRaphson
(prb, name=None, maxNumIter=10, convergenceTestTol=1e-09, printFlag=0, numSteps=1, numberingMethod='rcm', convTestType='norm_unbalance_conv_test')¶ Bases:
solution.predefined_solutions.SolutionProcedure
Newton-Raphson solution algorithm with a plain constraint handler.
-
setup
()¶ Defines the solution procedure in the finite element problem object.
-
-
class
solution.predefined_solutions.
PlainNewtonRaphsonBandGen
(prb, name=None, maxNumIter=10, convergenceTestTol=1e-09, printFlag=0, numSteps=1, numberingMethod='rcm', convTestType='norm_unbalance_conv_test')¶ Bases:
solution.predefined_solutions.SolutionProcedure
Newton-Raphson solution algorithm with a s plain constraint handler and a band general SOE solver.
-
setup
()¶ Defines the solution procedure in the finite element problem object.
-
-
class
solution.predefined_solutions.
PlainStaticModifiedNewton
(prb, name=None, maxNumIter=10, convergenceTestTol=1e-09, printFlag=0, numSteps=1, numberingMethod='rcm', convTestType='relative_total_norm_disp_incr_conv_test')¶ Bases:
solution.predefined_solutions.SolutionProcedure
Static solution procedure with a modified Newton solution algorithm with a plain constraint handler.
-
setup
()¶ Defines the solution procedure in the finite element problem object.
-
-
class
solution.predefined_solutions.
SimpleLagrangeStaticLinear
(prb, name=None, maxNumIter=10, convergenceTestTol=1e-09, printFlag=0, numSteps=1, numberingMethod='rcm')¶ Bases:
solution.predefined_solutions.SolutionProcedure
Linear static solution algorithm with a Lagrange constraint handler.
-
setup
()¶ Defines the solution procedure in the finite element problem object.
-
-
class
solution.predefined_solutions.
SimpleStaticLinear
(prb, name=None, printFlag=0, numSteps=1, numberingMethod='rcm', soeType='band_spd_lin_soe', solverType='band_spd_lin_lapack_solver')¶ Bases:
solution.predefined_solutions.PenaltyStaticLinearBase
Return a linear static solution algorithm with a penalty constraint handler.
-
class
solution.predefined_solutions.
SimpleStaticLinearMUMPS
(prb, name=None, printFlag=0, numSteps=1, numberingMethod='rcm', soeType='mumps_soe', solverType='mumps_solver')¶ Bases:
solution.predefined_solutions.PenaltyStaticLinearBase
Return a linear static solution algorithm with a penalty constraint handler.
-
class
solution.predefined_solutions.
SimpleStaticLinearUMF
(prb, name=None, printFlag=0, numSteps=1, numberingMethod='rcm', soeType='umfpack_gen_lin_soe', solverType='umfpack_gen_lin_solver')¶ Bases:
solution.predefined_solutions.PenaltyStaticLinearBase
Return a linear static solution algorithm with a penalty constraint handler.
-
class
solution.predefined_solutions.
SimpleTransformationStaticLinear
(prb, name=None, maxNumIter=10, convergenceTestTol=1e-09, printFlag=0, numSteps=1)¶ Bases:
solution.predefined_solutions.SolutionProcedure
Linear static solution algorithm with a transformation constraint handler.
-
setup
()¶ Defines the solution procedure in the finite element problem object.
-
-
class
solution.predefined_solutions.
SolutionProcedure
(name=None, constraintHandlerType='plain', maxNumIter=10, convergenceTestTol=1e-09, printFlag=0, numSteps=1, numberingMethod='rcm', convTestType=None, soeType: str = None, solverType: str = None, shift: float = None)¶ Bases:
object
Variables: - solutionStrategy – object that aggregates all the component objects that define the type of analysis (constraint handler, DOF_Numberer, integrator, solution algorithm, system of equations, convergence test)
- cHandler – constraint handler. Determines how the constraint equations are enforced in the analysis, how it handles the boundary conditions/imposed displacements
- numberer – DOF numberer. Determines the mapping between equation numbers and degrees of freedom (how DOF are numbered)
- integ –
integrator, used to:
- determine the predictive step for time t+dt
- specify the tangent matrix and residual vector at any iteration
- determine the corrective step based on the displacement increment
- solAlgo –
solution algorithm, which determines the sequence of steps taken to solve the non-linear equation:
- linear
- Newton Raphson (the tangent is updated at each iteration)
- Newton with Line Search
- Modified Newton Raphson (the tangent is not updated at
each iteration) - Kyrlov-Newton - BFGS, for symetric systems - Broyden
- soe – sparse system of equations: band general, band SPD, profile SPD, sparse general, umfPack, sparse SPD
- solver – solver for the system of equations.
- analysis – determines what type of analysis is to be performed
- convergenceTestTol – convergence tolerance (defaults to 1e-9)
- printFlag – if not zero print convergence results on each step.
- numSteps – number of steps to use in the analysis (useful only when loads are variable in time).
- numberingMethod – numbering method (plain or reverse Cuthill-McKee or alterntive minimum degree).
- convTestType – convergence test type for non linear analysis (norm unbalance,…).
- soeType – type of the system of equations object.
- solverType – type of the solver.
- maxNumIter – maximum number of iterations (defauts to 10)
- numSteps – number of steps to use in the analysis (useful only when loads are variable in time).
- solu –
- solCtrl –
- modelWrapper – model representation for the analysis.
- shift – shift-and-invert mode (used with ARPACK).
-
analysisSetup
(analysisType='static_analysis')¶ Define the analysis object.
Parameters: analysisType – type of the analysis to perform.
-
clear
()¶ Wipe out the solution procedure.
-
constraintHandlerSetup
()¶ Define the constraint handler and return a reference to it.
-
getConvergenceTest
()¶ Return the convergence test.
-
getModelWrapperName
()¶ Return the name for the model wrapper.
-
getSolutionStrategyName
()¶ Return the name for the model wrapper.
-
integratorSetup
(integratorType='load_control_integrator')¶ Define the type of integrator to use in the analysis.
Parameters: integratorType – type of integrator to use.
-
modelWrapperSetup
()¶ Defines the model wrapper.
-
resetLoadCase
()¶ Remove previous load from the domain.
-
setPenaltyFactors
(alphaSP=1000000000000000.0, alphaMP=1000000000000000.0)¶ Define the penalty factors to use with the penalty constraint handler.
Parameters: - alphaSP – penalty factor on single points constraints (defaults to 1e15).
- alphaMP – penalty factor on multi-poing constraints (defaults to 1e15).
-
setup
()¶ Defines the solution procedure in the finite element problem object.
-
soluControlSetup
()¶ Defines the solution control object.
-
solutionAlgorithmSetup
(solAlgType='linear_soln_algo', integratorType='load_control_integrator')¶ Define the solution strategy.
Parameters: - solAlgType – type of the solution algorithm (linear, Newton, modified Newton, …)
- integratorType – type of integrator to use.
-
solve
(calculateNodalReactions=False, includeInertia=False, reactionCheckTolerance=1e-12)¶ Compute the solution (run the analysis).
Parameters: - calculateNodalReactions – if true calculate reactions at nodes.
- includeInertia – if true calculate reactions including inertia effects.
- reactionCheckTolerance – tolerance when checking reaction values.
-
solveComb
(combName, calculateNodalReactions=False, includeInertia=False, reactionCheckTolerance=1e-12)¶ Obtains the solution for the combination argument.
Parameters: - combName – name of the combination to obtain the response for.
- calculateNodalReactions – if true calculate reactions at nodes.
- includeInertia – if true calculate reactions including inertia effects.
- reactionCheckTolerance – tolerance when checking reaction values.
-
sysOfEqnSetup
()¶ Defines the solver to use for the resulting system of equations.
-
class
solution.predefined_solutions.
TransformationNewtonLineSearch
(prb, name=None, maxNumIter=150, convergenceTestTol=1e-09, printFlag=0, numSteps=1, numberingMethod='rcm', convTestType='relative_total_norm_disp_incr_conv_test', lineSearchMethod='regula_falsi_line_search')¶ Bases:
solution.predefined_solutions.TransformationNewtonLineSearchBase
Static solution procedure with a Newton line search algorithm and a transformation constraint handler.
-
setup
()¶ Defines the solution procedure in the finite element problem object.
-
-
class
solution.predefined_solutions.
TransformationNewtonLineSearchBase
(prb, name, maxNumIter, convergenceTestTol, printFlag, numSteps, numberingMethod, convTestType, lineSearchMethod, soeType, solverType)¶ Bases:
solution.predefined_solutions.LineSearchBase
Base class for transformation Newton line search solution aggregation.
-
setup
()¶ Defines the solution procedure in the finite element problem object.
-
-
class
solution.predefined_solutions.
TransformationNewtonLineSearchMUMPS
(prb, name=None, maxNumIter=150, convergenceTestTol=1e-09, printFlag=0, numSteps=1, numberingMethod='rcm', convTestType='relative_total_norm_disp_incr_conv_test', lineSearchMethod='regula_falsi_line_search')¶ Bases:
solution.predefined_solutions.TransformationNewtonLineSearchBase
Static solution procedure with a Newton line search algorithm, a transformation constraint handler and a MUMPS (parallel sparse direct solver) solver.
-
setup
()¶ Defines the solution procedure in the finite element problem object.
-
-
class
solution.predefined_solutions.
TransformationNewtonLineSearchUMF
(prb, name=None, maxNumIter=150, convergenceTestTol=1e-09, printFlag=0, numSteps=1, numberingMethod='rcm', convTestType='relative_total_norm_disp_incr_conv_test', lineSearchMethod='regula_falsi_line_search')¶ Bases:
solution.predefined_solutions.TransformationNewtonLineSearchBase
Static solution procedure with a Newton line search algorithm, a transformation constraint handler and a UMF (Unsimmetric multi-frontal method) solver.
-
setup
()¶ Defines the solution procedure in the finite element problem object.
-
-
class
solution.predefined_solutions.
ZeroEnergyModes
(prb, name=None, printFlag=0, systemPrefix='sym_band_eigen', shift: float = None)¶ Bases:
solution.predefined_solutions.IllConditioningAnalysisBase
Procedure to obtain the zero energy modes of the finite element model.
-
solution.predefined_solutions.
frequency_analysis
(prb, systemPrefix='sym_band', shift: float = None)¶ Return a solution procedure that computes the natural frequencies of the model.
-
solution.predefined_solutions.
ill_conditioning_analysis
(prb)¶ Return a solution procedure that computes the modes that correspond to ill-conditioned degrees of freedom.
-
solution.predefined_solutions.
penalty_modified_newton
(prb, mxNumIter=10, convergenceTestTol=0.0001, printFlag=0, convTestType='relative_total_norm_disp_incr_conv_test')¶ Return a simple static modified Newton solution procedure.
Parameters: - maxNumIter – maximum number of iterations (defauts to 10)
- convergenceTestTol – convergence tolerance (defaults to 1e-9)
- printFlag – print message on each iteration
-
solution.predefined_solutions.
penalty_newton_line_search_mumps
(prb, mxNumIter=150, convergenceTestTol=1e-09, printFlag=0, convTestType='relative_total_norm_disp_incr_conv_test', lineSearchMethod='regula_falsi_line_search')¶ Return a simple static modified Newton solution procedure.
Parameters: - maxNumIter – maximum number of iterations (defauts to 10)
- convergenceTestTol – convergence tolerance (defaults to 1e-9)
- printFlag – print message on each iteration
-
solution.predefined_solutions.
penalty_newton_raphson
(prb, mxNumIter=10, convergenceTestTol=0.0001, printFlag=0)¶ Return a penalty Newton-Raphson solution procedure.
Parameters: - maxNumIter – maximum number of iterations (defauts to 10)
- convergenceTestTol – convergence tolerance (defaults to 1e-9)
- printFlag – print message on each iteration
-
solution.predefined_solutions.
plain_krylov_newton
(prb, mxNumIter=300, convergenceTestTol=1e-09, printFlag=0, maxDim=6)¶ Return a plain Krylov Newton solution procedure.
Parameters: - maxNumIter – maximum number of iterations (defauts to 300)
- convergenceTestTol – convergence tolerance (defaults to 1e-9)
- printFlag – print message on each iteration
- maxDim – max number of iterations until the tangent is reformed and the acceleration restarts (default = 6).
-
solution.predefined_solutions.
plain_newton_raphson
(prb, mxNumIter=10)¶
-
solution.predefined_solutions.
plain_newton_raphson_band_gen
(prb, mxNumIter=10)¶
-
solution.predefined_solutions.
plain_static_modified_newton
(prb, mxNumIter=10, convergenceTestTol=0.01, printFlag=0)¶ Return a simple static modified Newton solution procedure.
Parameters: - maxNumIter – maximum number of iterations (defauts to 10)
- convergenceTestTol – convergence tolerance (defaults to 1e-9)
- printFlag – print message on each iteration
-
solution.predefined_solutions.
simple_static_linear
(prb)¶ Return a simple static linear solution procedure.
-
solution.predefined_solutions.
solveCombEstat2ndOrderLin
(combName, solutionProcedure)¶
-
solution.predefined_solutions.
solveStaticLinearComb
(combName, solutionProcedure)¶
-
solution.predefined_solutions.
solveStaticNoLinCase
(combName)¶
-
solution.predefined_solutions.
transformation_newton_line_search_mumps
(prb, mxNumIter=150, convergenceTestTol=1e-09, printFlag=0, convTestType='relative_total_norm_disp_incr_conv_test')¶ Return a simple static modified Newton solution procedure.
Parameters: - maxNumIter – maximum number of iterations (defauts to 10)
- convergenceTestTol – convergence tolerance (defaults to 1e-9)
- printFlag – print message on each iteration
-
solution.predefined_solutions.
zero_energy_modes
(prb)¶ Return a solution procedure that computes the zero energy modes of the model.