Solvers

Currently, EdelweissFE provides

  • a nonlinear implicit static solver (NIST),

  • a parallel nonlinear implicit static solver (NISTParallel),

  • a parallel nonlinear implicit static solver tuned for marmot elements (NISTParallelForMarmotElements),

  • and a parallel arc length solver (NISTPArcLength).

Choose the solver in the *solver definition:

*solver, name=mySolver, solver=NISTParallel
edelweissfe.config.solvers.getSolverByName(name)[source]

Get the class type of the requested solver.

Parameters

name (str) – The name of the solver to load.

Returns

The solver class type.

Return type

type

Currently available:

Option

Description

NIST

nonlinearimplicitstatic

NISTParallel

nonlinearimplicitstaticparallelmk2

NISTParallelForMarmotElements

nonlinearimplicitstaticparallel

NISTPArcLength

nonlinearimplicitstaticparallelarclength

NIST - Nonlinear Implicit Static

class edelweissfe.solvers.nonlinearimplicitstatic.NIST(jobInfo, journal, **kwargs)[source]

This is the Nonlinear Implicit STatic – solver.

Parameters
  • jobInfo – A dictionary containing the job information.

  • journal – The journal instance for logging.

Methods

applyDirichlet(timeStep, R, dirichlets)

Apply the dirichlet bcs on the residual vector Is called by solveStep() before solving the global equatuon system.

applyDirichletK(K, dirichlets)

Apply the dirichlet bcs on the global stiffness matrix Is called by solveStep() before solving the global system.

applyStepActionsAtIncrementStart(model, ...)

Called when all step actions should be applied at the start of a step.

applyStepActionsAtStepEnd(model, stepActions)

Called when all step actions should finish a step.

applyStepActionsAtStepStart(model, stepActions)

Called when all step actions should be appliet at the start a step.

assembleConstraints(constraints, U_np, dU, ...)

Loop over all elements, and evaluate them.

assembleLoads(nodeForces, distributedLoads, ...)

Assemble all loads into a right hand side vector.

assembleStiffnessCSR(K)

Construct a CSR matrix from VIJ format.

checkConvergence(R, ddU, F, ...)

Check the convergence, individually for each field, similar to Abaqus based on the current total flux residual and the field correction Is called by solveStep() to decide whether to continue iterating or stop.

checkDivergingSolution(...)

Check if the iterative solution scheme is diverging.

computeBodyForces(bodyForces, U_np, PExt, K, ...)

Loop over all body forces loads acting on elements, and evaluate them.

computeDistributedLoads(distributedLoads, ...)

Loop over all distributed loads acting on elements, and evaluate them.

computeElements(elements, U_np, dU, P, K, F, ...)

Loop over all elements, and evalute them.

computeSpatialAveragedFluxes(F)

Compute the spatial averaged flux for every field Is usually called by checkConvergence().

extrapolateLastIncrement(extrapolation, ...)

Depending on the current setting, extrapolate the solution of the last increment.

linearSolve(A, b)

Solve the linear equation system.

printResidualOutlierNodes(residualOutliers)

Print which nodes have the largest residuals.

solveIncrement(U_n, dU, P, K, stepActions, ...)

Standard Newton-Raphson scheme to solve for an increment.

solveStep(step, model, ...)

Public interface to solve for a step.

findDirichletIndices

applyDirichlet(timeStep, R, dirichlets)[source]

Apply the dirichlet bcs on the residual vector Is called by solveStep() before solving the global equatuon system.

Parameters
Returns

The modified residual vector.

Return type

DofVector

applyDirichletK(K, dirichlets)[source]

Apply the dirichlet bcs on the global stiffness matrix Is called by solveStep() before solving the global system. http://stackoverflux.com/questions/12129948/scipy-sparse-set-row-to-zeros

Parameters
Returns

The modified system matrix.

Return type

VIJSystemMatrix

applyStepActionsAtIncrementStart(model, timeStep, stepActions)[source]

Called when all step actions should be applied at the start of a step.

Parameters
applyStepActionsAtStepEnd(model, stepActions)[source]

Called when all step actions should finish a step.

Parameters
applyStepActionsAtStepStart(model, stepActions)[source]

Called when all step actions should be appliet at the start a step.

Parameters
assembleConstraints(constraints, U_np, dU, PExt, K, timeStep)[source]

Loop over all elements, and evaluate them. Is is called by solveStep() in each iteration.

Parameters
Returns

  • The modified external load vector.

  • The modified system matrix.

Return type

tuple[DofVector,VIJSystemMatrix,DofVector]

assembleLoads(nodeForces, distributedLoads, bodyForces, U_np, PExt, K, timeStep)[source]

Assemble all loads into a right hand side vector.

Parameters
Returns

  • The augmented external load vector.

  • The augmented system matrix.

Return type

tuple[DofVector,VIJSystemMatrix]

assembleStiffnessCSR(K)[source]

Construct a CSR matrix from VIJ format.

Parameters

K (VIJSystemMatrix) – The system matrix in VIJ format.

Returns

The system matrix in compressed sparse row format.

Return type

csr_matrix

checkConvergence(R, ddU, F, iterationCounter, residualHistory)[source]

Check the convergence, individually for each field, similar to Abaqus based on the current total flux residual and the field correction Is called by solveStep() to decide whether to continue iterating or stop.

Parameters
  • R (DofVector) – The current residual.

  • ddU (DofVector) – The current correction increment.

  • F (DofVector) – The accumulated fluxes.

  • iterationCounter (int) – The current iteration number.

  • residualHistory (dict) – The previous residuals.

Returns

  • True if converged.

  • The residual histories field wise.

Return type

tuple[bool,dict]

checkDivergingSolution(incrementResidualHistory, maxGrowingIter)[source]

Check if the iterative solution scheme is diverging.

Parameters
  • incrementResidualHistory (dict) – The dictionary containing the residual history of all fields.

  • maxGrowingIter (int) – The maximum allows number of growths of a residual during the iterative solution scheme.

Returns

True if solution is diverging.

Return type

bool

computeBodyForces(bodyForces, U_np, PExt, K, timeStep)[source]

Loop over all body forces loads acting on elements, and evaluate them. Assembles into the global external load vector and the system matrix.

Parameters
Returns

The augmented load vector and system matrix.

Return type

tuple[DofVector,VIJSystemMatrix]

computeDistributedLoads(distributedLoads, U_np, PExt, K, timeStep)[source]

Loop over all distributed loads acting on elements, and evaluate them. Assembles into the global external load vector and the system matrix.

Parameters
Returns

The augmented load vector and system matrix.

Return type

tuple[DofVector,VIJSystemMatrix]

computeElements(elements, U_np, dU, P, K, F, timeStep)[source]

Loop over all elements, and evalute them. Is is called by solveStep() in each iteration.

Parameters
  • elements (list) – The list of finite elements.

  • U_np (DofVector) – The current solution vector.

  • dU (DofVector) – The current solution increment vector.

  • P (DofVector) – The reaction vector.

  • K (VIJSystemMatrix) – The system matrix.

  • F (DofVector) – The vector of accumulated fluxes for convergence checks.

  • timeStep (TimeStep) – The time step.

Returns

  • The modified reaction vector.

  • The modified system matrix.

  • The modified accumulated flux vector.

Return type

tuple[DofVector,VIJSystemMatrix,DofVector]

computeSpatialAveragedFluxes(F)[source]

Compute the spatial averaged flux for every field Is usually called by checkConvergence().

Parameters

F (DofVector) – The accumulated flux vector.

Returns

A dictioary containg the spatial average fluxes for every field.

Return type

dict[str,float]

extrapolateLastIncrement(extrapolation, timeStep, dU, dirichlets, prevTimeStep, model)[source]

Depending on the current setting, extrapolate the solution of the last increment.

Parameters
  • extrapolation (str) – The type of extrapolation.

  • timeStep (TimeStep) – The current time step.

  • dU (DofVector) – The last solution increment.

  • dirichlets (list) – The list of active dirichlet boundary conditions.

  • lastIncrementSize – The size of the last increment.

  • prevTimeStep (TimeStep) –

Returns

  • The extrapolated solution increment.

  • True if an extrapolation was performed.

Return type

tuple[DofVector,bool]

linearSolve(A, b)[source]

Solve the linear equation system.

Parameters
  • A (csr_matrix) – The system matrix in compressed spare row format.

  • b (DofVector) – The right hand side.

Returns

The solution ‘x’.

Return type

ndarray

printResidualOutlierNodes(residualOutliers)[source]

Print which nodes have the largest residuals.

Parameters

residualOutliers (dict) – The dictionary containing the outlier nodes for every field.

solveIncrement(U_n, dU, P, K, stepActions, model, timeStep, prevTimeStep, extrapolation, maxIter, maxGrowingIter)[source]

Standard Newton-Raphson scheme to solve for an increment.

Parameters
  • Un – The old solution vector.

  • dU (DofVector) – The old solution increment.

  • P (DofVector) – The old reaction vector.

  • K (VIJSystemMatrix) – The system matrix to be used.

  • elements – The dictionary containing all elements.

  • stepActions (list) – The list of active step actions.

  • model (FEModel) – The model tree.

  • increment – The increment.

  • lastIncrementSize – The size of the previous increment.

  • extrapolation (str) – The type of extrapolation to be used.

  • maxIter (int) – The maximum number of iterations to be used.

  • maxGrowingIter (int) – The maximum number of growing residuals until the Newton-Raphson is terminated.

  • U_n (DofVector) –

  • timeStep (TimeStep) –

  • prevTimeStep (TimeStep) –

Returns

A tuple containing
  • the new solution vector

  • the solution increment

  • the new reaction vector

  • the number of required iterations

  • the history of residuals per field

Return type

tuple[DofVector,DofVector,DofVector,int,dict]

solveStep(step, model, fieldOutputController, outputmanagers)[source]

Public interface to solve for a step.

Parameters
  • stepNumber – The step number.

  • step – The dictionary containing the step definition.

  • stepActions – The dictionary containing all step actions.

  • model (FEModel) – The model tree.

  • fieldOutputController (FieldOutputController) – The field output controller.

  • outputmanagers (dict[str, edelweissfe.outputmanagers.base.outputmanagerbase.OutputManagerBase]) –

Return type

tuple[bool, edelweissfe.models.femodel.FEModel]

NISTParallelForMarmotElements - Nonlinear Implicit Static (parallel)

Parallel implementation of the NIST solver.

class edelweissfe.solvers.nonlinearimplicitstaticparallelmk2.NISTParallel(jobInfo, journal, **kwargs)

Methods

applyDirichlet(timeStep, R, dirichlets)

Apply the dirichlet bcs on the residual vector Is called by solveStep() before solving the global equatuon system.

applyDirichletK(K, dirichlets)

Apply the dirichlet bcs on the global stiffnes matrix Is called by solveStep() before solving the global sys.

applyStepActionsAtIncrementStart(model, ...)

Called when all step actions should be applied at the start of a step.

applyStepActionsAtStepEnd(model, stepActions)

Called when all step actions should finish a step.

applyStepActionsAtStepStart(model, stepActions)

Called when all step actions should be appliet at the start a step.

assembleConstraints(constraints, U_np, dU, ...)

Loop over all elements, and evaluate them.

assembleLoads(nodeForces, distributedLoads, ...)

Assemble all loads into a right hand side vector.

assembleStiffnessCSR(K)

Construct a CSR matrix from VIJ format.

checkConvergence(R, ddU, F, ...)

Check the convergence, individually for each field, similar to Abaqus based on the current total flux residual and the field correction Is called by solveStep() to decide whether to continue iterating or stop.

checkDivergingSolution(...)

Check if the iterative solution scheme is diverging.

computeBodyForces(bodyForces, U_np, PExt, K, ...)

Loop over all body forces loads acting on elements, and evaluate them.

computeDistributedLoads(distributedLoads, ...)

Loop over all distributed loads acting on elements, and evaluate them.

computeSpatialAveragedFluxes(F)

Compute the spatial averaged flux for every field Is usually called by checkConvergence().

extrapolateLastIncrement(extrapolation, ...)

Depending on the current setting, extrapolate the solution of the last increment.

linearSolve(A, b)

Solve the linear equation system.

printResidualOutlierNodes(residualOutliers)

Print which nodes have the largest residuals.

solveIncrement(U_n, dU, P, K, stepActions, ...)

Standard Newton-Raphson scheme to solve for an increment.

computeElements

findDirichletIndices

solveStep

applyDirichletK(K, dirichlets)

Apply the dirichlet bcs on the global stiffnes matrix Is called by solveStep() before solving the global sys. http://stackoverflux.com/questions/12129948/scipy-sparse-set-row-to-zeros

Cythonized version for speed!

Parameters
  • K – The system matrix.

  • dirichlets – The list of dirichlet boundary conditions.

Returns

The modified system matrix.

Return type

VIJSystemMatrix

computeElements(elements, Un1, dU, P, K, F, timeStep)

Loop over all elements, and evalute them. Is is called by solveStep() in each iteration.

Parameters
  • elements – The list of finite elements.

  • U_np – The current solution vector.

  • dU – The current solution increment vector.

  • P – The reaction vector.

  • K – The system matrix.

  • F – The vector of accumulated fluxes for convergence checks.

  • timeStep – The time step.

Returns

  • The modified reaction vector.

  • The modified system matrix.

  • The modified accumulated flux vector.

Return type

tuple[DofVector,VIJSystemMatrix,DofVector]

solveStep(step, model, fieldOutputController, outputmanagers)

Public interface to solve for a step.

Parameters
  • stepNumber – The step number.

  • step – The dictionary containing the step definition.

  • stepActions – The dictionary containing all step actions.

  • model – The model tree.

  • fieldOutputController – The field output controller.

NISTPArcLength - Nonlinear Implicit Static - Arc length

(Parallel) Arc Length Solver, based on the proposed approach in Jirásek/Bažant 2001. Replaces the NewtonRaphson scheme of the NISTParallel Solver.

class edelweissfe.solvers.nonlinearimplicitstaticparallelarclength.NISTPArcLength(jobInfo, journal, **kwargs)[source]

Methods

applyDirichlet(timeStep, R, dirichlets)

Apply the dirichlet bcs on the residual vector Is called by solveStep() before solving the global equatuon system.

applyDirichletK(K, dirichlets)

Apply the dirichlet bcs on the global stiffnes matrix Is called by solveStep() before solving the global sys.

applyStepActionsAtIncrementStart(model, ...)

Called when all step actions should be applied at the start of a step.

applyStepActionsAtStepEnd(model, stepActions)

Called when all step actions should finish a step.

applyStepActionsAtStepStart(model, stepActions)

Called when all step actions should be appliet at the start a step.

assembleConstraints(constraints, U_np, dU, ...)

Loop over all elements, and evaluate them.

assembleLoads(nodeForces, distributedLoads, ...)

Assemble all loads into a right hand side vector.

assembleStiffnessCSR(K)

Construct a CSR matrix from VIJ format.

checkConvergence(R, ddU, F, ...)

Check the convergence, individually for each field, similar to Abaqus based on the current total flux residual and the field correction Is called by solveStep() to decide whether to continue iterating or stop.

checkDivergingSolution(...)

Check if the iterative solution scheme is diverging.

computeBodyForces(bodyForces, U_np, PExt, K, ...)

Loop over all body forces loads acting on elements, and evaluate them.

computeDistributedLoads(distributedLoads, ...)

Loop over all distributed loads acting on elements, and evaluate them.

computeSpatialAveragedFluxes(F)

Compute the spatial averaged flux for every field Is usually called by checkConvergence().

extrapolateLastIncrement(extrapolation, ...)

Depending on the current setting, extrapolate the solution of the last increment.

linearSolve(A, b)

Solve the linear equation system.

printResidualOutlierNodes(residualOutliers)

Print which nodes have the largest residuals.

solveIncrement(U_n, dU, P, K, stepActions, ...)

Arc length method scheme to solve for an increment.

solveStep(step, model, ...)

Public interface to solve for a step.

computeElements

findDirichletIndices

applyStepActionsAtStepEnd(model, stepActions)[source]

Called when all step actions should finish a step.

For the arc length controlled simulation, external loads at the end of a step depend on the arc length parameter. Hence, for subsequent steps, the new reference loads depend on this parameter. We communicate this to the step actions before the parent solver tells them that the step is finished.

Parameters
extrapolateLastIncrement(extrapolation, timeStep, dU, dirichlets, prevTimeStep, model, dLambda=None)[source]

Depending on the current setting, extrapolate the solution of the last increment. Also accounts for the arc length loading parameter.

Parameters
  • extrapolation (str) – The type of extrapolation.

  • timeStep (TimeStep) – The current time increment.

  • dU (DofVector) – The last solution increment.

  • dirichlets (list) – The list of active dirichlet boundary conditions.

  • prevTimeStep (TimeStep) – The previous TimeStep

  • model (FEModel) –

  • dLambda (float) –

Returns

  • The extrapolated solution increment.

  • True if an extrapolation was performed.

Return type

tuple[DofVector,bool]

solveIncrement(U_n, dU, P, K, stepActions, model, timeStep, prevTimeStep, extrapolation, maxIter, maxGrowingIter)[source]

Arc length method scheme to solve for an increment.

Parameters
  • U_n (DofVector) – The old solution vector.

  • dU (DofVector) – The old solution increment.

  • P (DofVector) – The old reaction vector.

  • K (VIJSystemMatrix) – The system matrix to be used.

  • elements – The dictionary containing all elements.

  • activeStepActions – The list of active step actions.

  • constraints – The dictionary containing all elements.

  • timeStep (TimeStep) – The current time step.

  • prevTimeStep (TimeStep) – The previous time step.

  • extrapolation (str) – The type of extrapolation to be used.

  • maxIter (int) – The maximum number of iterations to be used.

  • maxGrowingIter (int) – The maximum number of growing residuals until the Newton-Raphson is terminated.

  • stepActions (list) –

Returns

A tuple containing
  • the new solution vector

  • the solution increment

  • the new reaction vector

  • the number of required iterations

  • the history of residuals per field

Return type

tuple[DofVector,DofVector,DofVector,int,dict]

solveStep(step, model, fieldOutputController, outputmanagers)[source]

Public interface to solve for a step.

Parameters
  • stepNumber – The step number.

  • step (dict) – The dictionary containing the step definition.

  • stepActions – The dictionary containing all step actions.

  • model (FEModel) – The model tree.

  • fieldOutputController (FieldOutputController) – The field output controller.

  • outputmanagers (dict[str, edelweissfe.outputmanagers.base.outputmanagerbase.OutputManagerBase]) –