FemProcess
The GeMA Fem Process Plugin
Nonlinear transient solver options

Solving a nonlinear transient problem with 'transient nonlinear' scheme

The solver options to nonlinear transient solvers with 'transient nonlinear' scheme can be defined in GeMA by using code simillar to the following fragment:

local solverOptions = {
type = 'transient nonlinear',
timeMax = 100.0,
timeInitIncrement = 0.1,
timeMinIncrement = 0.001,
timeMaxIncrement = 20.0,
iterationsMax = 15,
eulerTheta = 1.0,
tolerance = {mechanic = 1.0e-04, hydraulic = 1.0e-04}
}

In addition, you can define some controls for time increment. In this case, the option timeControl is defined as true (default false). Then, you can modify the default time control parameters:

local solverOptions = {
timeControl = true, -- active TIME INCREMENTATION PARAMETERS
IncT = 3, -- minimum number of consecutive increments in which the time integration accuracy measure must be satisfied without any cutbacks to allow a time increment increase
IterationsL = 10, -- number of consecutive equilibrium iterations (without severe discontinuities) above which the size of the next increment will be reduced.
IterationsG = 4, -- maximum number of consecutive equilibrium iterations (without severe discontinuities) allowed in consecutive increments for the time increment to be increased
DLTM = 2, -- maximum time increment increase factor for diffusion-dominated processes
DLTG = 0.8, -- increase factor for the next time increment, as a ratio of the average integration accuracy measure over IncT increments to the corresponding tolerance
DLTB = 0.75, -- cutback factor for the next increment when too many equilibrium iterations (IterationsL) are used in the current increment
DLTF = 0.25, -- cutback factor used when the solution appears to be diverging
}

Solving a nonlinear transient problem with 'transient automatic time step' scheme

The solver options to nonlinear transient solvers with 'transient automatic time step' scheme can be defined in GeMA by using code simillar to the following fragment:

local solverOptions = {
type = 'transient automatic time step',
timeMax = 100.0,
timeInitIncrement = 0.1,
timeMinIncrement = 0.00001,
timeMaxIncrement = 20.0,
iterationsMax = 15,
DTOL = 1.0e-03,
ITOL = 1.0e-04,
}

fem.init(physicsList, numericSolver, solverOptions)
Description: Creates and initializes a transient nonlinear FEM solver intended to be used on a time loop as exemplified above. Returns a solver object that should be passed as parameter in calls to fem.Step().
Parameters: physicsList A string with the id of the physics responsible for providing the FEM process with the equations to be solved. Can also be a table with multiple physics ids. All physics objects must be tied to the same mesh.
numericSolver A string with the id of the numeric solver that will be used to solve the final linear system.
solverOptions The nonlinear transient solution escheme allows setting flags with options for the solver. These flags can be: requiered, optional and advanced.

A required Lua table with a set of flags defining options for the nonlinear transient solver accepts the following options:
- type: Defines the type of transient nonlinear solver available, two methods are implemented and a name or type is assigned: 'transient nonlinear and 'transient automatic time step', the first based on the Euler scheme and the second based on the automatic adaptive scheme with local error control proposed by Sloan & Abbo.
- timeMax: Sets the total time of analysis.
- timeInitIncrement: Sets the initial time increment.
- timeMinIncrement: An options that defines the minimum increment of time allowed for the analysis.
- timeMaxIncrement: An option that defines the maximum time increment allowed during the analysis.
- iterationMax: Sets the maximum number of iterations allowed in the Newton-Raphson scheme.
- tolerance: Options that allowed to define a tolerance, for the control of the error during the iterative process. Tolerance is defined by each physics involved in the analysis: mechanic, hydraulic, thermic and chemical. If the user does not set the required tolerance value, the default value of 1.0e-05 will be used.

An optional Lua table with a set of flags defining options for the solver. The transient linear solver accepts the following options:
- assemblerDofMode: A an optional string parameter that instructs the solver how fixed dofs should be treated. By default (assemblerDofMode empty or equal to "automatic"), the solver selects how to work. The option "add fixed dofs" tells the solver that global matrices should include lines for fixed dofs that will be removed prior to solving the system. The option "remove fixed dofs" tells the solver that the assembler should remove fixed dofs and don't include them on the assembled matrix at all. This is a superior option (and the default option for the linear solver) but at the moment can not be used when boundary condition application points change over the simulation time (it can be used though if bc values are changing, without changing the nodes / edges / faces to which they are applied).
- numThreads: An option that defines the desired number of worker threads that will be used. If this option is missing, or is equal to -1, the project configured maximum number of concurrent threads, maxThreads from the simulation options, will be used. A value greater than maxThreads will revert to that maximum. A value of zero disables threading, meaning that tasks will be executed by the main thread. A value of 1, although weird, can be used for testing: all tasks are executed by a single thread, different from the main one.
- numTasks: An option that defines the number of tasks that the assembler will use to process mesh elements. If this option is missing or is equal to 0, the default number of tasks, defNumTasks from the simulation options, will be used. A value greater than 0 specifies a fixed number of tasks while a value less than 0 specifies a multiple of the number of worker threads, so a value of 5 means exactly 5 tasks, while a value of -2 means that the number of tasks should be equal to twice the number of worker threads.
- cellPartitionStrategy: An option that defines the default strategy that will be used to partition mesh cells into tasks. If this option is missing, the default strategy, defCellPartitionStrategy from the simulation options, will be used. Accepted values are "default" and "sequential".
- printOptions: A table with a set of print options for the solver. This options can instruct the solver to print all of its internal matrices and are usefull when debugging new physics implementations or to fully understand the FEM process behavior.
Returns: The solver object.

Nonlinear transient solvers with 'transient nonlinear' and 'transient automatic time step' scheme can be orchestrated in GeMA by using code simillar to the following fragment

function ProcessScript()
-- output file
local file = io.prepareMeshFile('mesh', '$SIMULATIONDIR/out/$SIMULATIONNAME', 'nf', {'listNodalVariables'}, {'listGaussPointVariables'},{split = true, saveDisplacements = true})
-- Create the solver model
local solver = fem.init({'physicsId',}, 'numericSolverId', solverOptions)
local dt = solverOptions.timeInitIncrement
local FinalTime = solverOptions.timeMax
local Time = dt
local LastStep = false
while (Time <= FinalTime) do
-- Solve to one step
dtnew = fem.step(solver, dt)
-- add results to output file
io.addResultToMeshFile(file, Time)
-- Adjust time to guarantee that the last iteration will be on the
-- requested final simulation time
if (Time + dt > FinalTime and not LastStep) then
dt = FinalTime - Time
Time = FinalTime
if equal(dt, 0.0) then break end
LastStep = true
else
Time = Time + dt
end
end
io.closeMeshFile(file)
end

dtnew = fem.step(solver, dt)
Description: Executes a new time step of the transient nonlinear solver created by the previous call to fem.init(). Should be called inside a loop until reaching the maxTime sets of the simulation.
Parameters: solver The solver object returned by fem.initTransientSolver().
dt The time step that should be used for advancing the simulation time. The value is expected to be given in the current time simulation unit defined in a call to setCurrentTimeUnit(), or in seconds if this function was not called.
dtnew The new time step calculated, should be used for advancing the simulation time. setCurrentTimeUnit(), or in seconds if this function was not called.
Returns: Nothing.

Nonlinear transient solvers with 'transient nonlinear' scheme can be orchestrated in GeMA by using code simillar to the following fragment

function ProcessScript()
-- output file
local file = io.prepareMeshFile('mesh', '$SIMULATIONDIR/out/$SIMULATIONNAME', 'nf', {'listNodalVariables'}, {'listGaussPointVariables'},{split = true, saveDisplacements = true})
-- Create the solver model
local solver = fem.init({'physicsId',}, 'numericSolverId', solverOptions)
local dt = solverOptions.timeInitIncrement
local FinalTime = solverOptions.timeMax
local Time = dt
local LastStep = false
while (Time <= FinalTime) do
-- Solve to one step
local dtnew, err, conv, niter = fem.step(solver, dt, true)
if(conv) then
-- add results to output file
io.addResultToMeshFile(file, Time)
dt = dtnew
Time = Time + dt
else
setCurrentTime(Time)
dt = dtnew
end
-- Adjust time to guarantee that the last iteration will be on the
-- requested final simulation time
if (Time + dt > FinalTime and not LastStep) then
dt = FinalTime - Time
Time = FinalTime
if equal(dt, 0.0) then break end
LastStep = true
else
Time = Time + dt
end
end
io.closeMeshFile(file)
end

dtnew,err,conv,niter = fem.step(solver, dt, true)
Description: Executes a new time step of the transient nonlinear solver created by the previous call to fem.init(). Should be called inside a loop until reaching the maxTime sets of the simulation.
Parameters: solver The solver object returned by fem.initTransientSolver().
dt The time step that should be used for advancing the simulation time. The value is expected to be given in the current time simulation unit defined in a call to setCurrentTimeUnit(), or in seconds if this function was not called.
dtnew The new time step calculated, should be used for advancing the simulation time. The user can decide to use this new calculated time increment, or set it to another value.
err The error value calculated during the iterative process is returned in this variable.
conv A Boolean is returned if convergence was possible for the current time increment. Then the conv variable will be true or false.
niter An integer value is returned that indicates how many iterations were required for the time increment during the iterations of the Newton-Raphson scheme.
Returns: Nothing.