TransientProblem

The main purpose of this class is to represent all relevant information for a transient analysis. This will include information defining the loading condition as well as various other pieces of information.

Note

This class should be created using the pyTACS.createTransientProblem method.

Options

Options can be set for TransientProblem at time of creation for the class in the pyTACS.createTransientProblem method or using the TransientProblem.setOption method. Current option values for a class instance can be printed out using the TransientProblem.printOption method. The following options their default values and descriptions are listed below:

+----------------------------------------+
|   TransientProblem default options:    |
+----------------------------------------+
'outputDir': './'
	 Output directory for F5 file writer.
'timeIntegrator': 'BDF'
	 Time integration scheme to use. Currently supports 'BDF' and 'DIRK'.
'integrationOrder': 2
	 Integration order for time marching scheme.
'L2Convergence': 1e-12
	 Absolute convergence tolerance for integrator based on l2 norm of residual.
'L2ConvergenceRel': 1e-12
	 Relative convergence tolerance for integrator based on l2 norm of residual.
'RBEStiffnessScaleFactor': 1000.0
	 Constraint matrix scaling factor used in RBE Lagrange multiplier stiffness matrix.
'RBEArtificialStiffness': 0.001
	 Artificial constant added to diagonals of RBE Lagrange multiplier stiffness matrix 
	 to stabilize preconditioner.
'jacAssemblyFreq': 1
	 How frequently to reassemble Jacobian during time integration process.
'writeSolution': True
	 Flag for suppressing all f5 file writing.
'numberSolutions': True
	 Flag for attaching solution counter index to f5 files.
'printTiming': False
	 Flag for printing out timing information for class procedures.
'printLevel': 0
	 Print level for integration solver.
	 Accepts:
		   0 : No printing.
		   1 : Print major iterations.
		 > 1 : Print major + minor iterations.

API Reference

class tacs.problems.TransientProblem(name, tInit, tFinal, numSteps, assembler, comm, outputViewer=None, meshLoader=None, isNonlinear=False, options=None)[source]

NOTE: This class should not be initialized directly by the user. Use pyTACS.createTransientProblem instead.

Parameters:
  • name (str) -- Name of this tacs problem

  • tInit (float) -- Starting time for transient problem integration

  • tFinal (float) -- Ending time for transient problem integration

  • numSteps (int) -- Number of time steps for transient problem integration

  • assembler (TACS.Assembler) -- Cython object responsible for creating and setting tacs objects used to solve problem

  • comm (mpi4py.MPI.Intracomm) -- The comm object on which to create the pyTACS object.

  • outputViewer (TACS.TACSToFH5) -- Cython object used to write out f5 files that can be converted and used for postprocessing.

  • meshLoader (pymeshloader.pyMeshLoader) -- pyMeshLoader object used to create the assembler.

  • options (dict) -- Dictionary holding problem-specific option parameters (case-insensitive).

setOption(name, value)[source]

Set a solver option value. The name is not case sensitive.

Parameters:
  • name (str) -- Name of option to modify

  • value (depends on option) -- New option value to set

getNumTimeSteps()[source]

Get the number of time steps used in time integration for this problem.

Returns:

numSteps -- Number of time steps.

Return type:

int

getTimeSteps()[source]

Get the discrete time step slices used in time integration.

Returns:

timeSteps -- Discrete time step slices used in time integration.

Return type:

numpy.ndarray[float]

getNumTimeStages()[source]

Get the number of time stages used for multi-stage time integration for this problem.

Returns:

numStages -- Number of time stages.

Return type:

int

getTimeStages(timeStep)[source]

Get the discrete time stage sub-intervals used in a multi-stage integration scheme.

Parameters:

timeStep (int) -- Time step index to get stage times for.

Returns:

timeStages -- Time step slices used to discretize this time step.

Return type:

numpy.ndarray[float]

addLoadToComponents(timeStep, compIDs, F, timeStage=None, averageLoad=False)[source]

This method is used to add a FIXED TOTAL LOAD on one or more components, defined by COMPIDs, at a specific time instance. The purpose of this routine is to add loads that remain fixed throughout an optimization. An example would be an engine load. This routine determines all the unique nodes in the FE model that are part of the requested components, then takes the total 'force' by F and divides by the number of nodes. This average load is then applied to the nodes.

Parameters:
  • timeStep (int) -- Time step index to apply load to.

  • compIDs (list[int] or int) -- The components with added loads. Use pyTACS.selectCompIDs method to determine this.

  • F (Numpy 1d or 2d array length (varsPerNodes) or (numCompIDs, varsPerNodes)) -- Vector(s) of 'force' to apply to each component. If only one force vector is provided, force will be copied uniformly across all components.

  • timeStage (int or None) -- Time stage index to apply load to. Default is None, which is applicable only for multi-step methods like BDF. For multi-stage methods like DIRK, this index must be specified.

  • averageLoad (bool) -- Flag to determine whether load should be split evenly across all components (True) or copied and applied individually to each component (False). Defaults to False.

Notes

The units of the entries of the 'force' vector F are not necessarily physical forces and their interpretation depends on the physics problem being solved and the dofs included in the model.

A couple of examples of force vector components for common problems are listed below:

In Heat Conduction with varsPerNode = 1

F = [Qdot] # heat rate

In Elasticity with varsPerNode = 3,

F = [fx, fy, fz] # forces

In Elasticity with varsPerNode = 6,

F = [fx, fy, fz, mx, my, mz] # forces + moments

In Thermoelasticity with varsPerNode = 4,

F = [fx, fy, fz, Qdot] # forces + heat rate

In Thermoelasticity with varsPerNode = 7,

F = [fx, fy, fz, mx, my, mz, Qdot] # forces + moments + heat rate

addLoadToNodes(timeStep, nodeIDs, F, timeStage=None, nastranOrdering=False)[source]

This method is used to add a fixed point load of F to the selected node IDs at a specified time instance.

Parameters:
  • timeStep (int) -- Time step index to apply load to.

  • nodeIDs (list[int]) -- The nodes IDs with added loads.

  • F (Numpy 1d or 2d array length (varsPerNodes) or (numNodeIDs, varsPerNodes)) -- Array of force vectors, one for each node. If only one force vector is provided, force will be copied uniformly across all nodes.

  • timeStage (int or None) -- Time stage index to apply load to. Default is None, which is applicable only for multi-step methods like BDF. For multi-stage methods like DIRK, this index must be specified.

  • nastranOrdering (bool) -- Flag signaling whether nodeIDs are in TACS (default) or NASTRAN (grid IDs in bdf file) ordering

Notes

The units of the entries of the 'force' vector F are not necessarily physical forces and their interpretation depends on the physics problem being solved and the dofs included in the model.

A couple of examples of force vector components for common problems are listed below:

In Heat Conduction with varsPerNode = 1

F = [Qdot] # heat rate

In Elasticity with varsPerNode = 3,

F = [fx, fy, fz] # forces

In Elasticity with varsPerNode = 6,

F = [fx, fy, fz, mx, my, mz] # forces + moments

In Thermoelasticity with varsPerNode = 4,

F = [fx, fy, fz, Qdot] # forces + heat rate

In Thermoelasticity with varsPerNode = 7,

F = [fx, fy, fz, mx, my, mz, Qdot] # forces + moments + heat rate

addLoadToRHS(timeStep, Fapplied, timeStage=None)[source]

This method is used to add a FIXED TOTAL LOAD directly to the right hand side vector given the equation below:

M*udotdot + K*u = f

Where:
  • K : Stiffness matrix for problem

  • u : State variables for problem

  • M : Mass matrix for problem

  • udotdot : Second time derivitive of state variables for problem

  • f : Right-hand side vector to add loads to

Parameters:
  • timeStep (int) -- Time step index to apply load to.

  • Fapplied (numpy.ndarray or tacs.TACS.Vec) -- Distributed array containing loads to applied to RHS of the problem.

  • timeStage (int or None) -- Time stage index to apply load to. Default is None, which is applicable only for multi-step methods like BDF. For multi-stage methods like DIRK, this index must be specified.

addTractionToComponents(timeStep, compIDs, tractions, timeStage=None, faceIndex=0)[source]

This method is used to add a FIXED TOTAL TRACTION on one or more components, defined by COMPIDs, at specified time instance. The purpose of this routine is to add loads that remain fixed throughout an optimization.

Parameters:
  • timeStep (int) -- Time step index to apply load to.

  • compIDs (list[int] or int) -- The components with added loads. Use pyTACS.selectCompIDs method to determine this.

  • tractions (numpy.ndarray length 1 or compIDs) -- Array of traction vectors for each component

  • timeStage (int or None) -- Time stage index to apply load to. Default is None, which is applicable only for multi-step methods like BDF. For multi-stage methods like DIRK, this index must be specified.

  • faceIndex (int) -- Indicates which face (side) of element to apply traction to. Note: not required for certain elements (i.e. shells)

addTractionToElements(timeStep, elemIDs, tractions, timeStage=None, faceIndex=0, nastranOrdering=False)[source]

This method is used to add a fixed traction to the selected element IDs at specified time instance. Tractions can be specified on an element by element basis (if tractions is a 2d array) or set to a uniform value (if tractions is a 1d array)

Parameters:
  • timeStep (int) -- Time step index to apply load to.

  • elemIDs (list[int]) -- The global element ID numbers for which to apply the traction.

  • tractions (numpy.ndarray 1d or 2d length varsPerNodes or (elemIDs, varsPerNodes)) -- Array of traction vectors for each element

  • timeStage (int or None) -- Time stage index to apply load to. Default is None, which is applicable only for multi-step methods like BDF. For multi-stage methods like DIRK, this index must be specified.

  • faceIndex (int) -- Indicates which face (side) of element to apply traction to. Note: not required for certain elements (i.e. shells)

  • nastranOrdering (bool) -- Flag signaling whether elemIDs are in TACS (default) or NASTRAN ordering

addPressureToComponents(timeStep, compIDs, pressures, timeStage=None, faceIndex=0)[source]

This method is used to add a FIXED TOTAL PRESSURE on one or more components, defined by COMPIDs, at specified time instance. The purpose of this routine is to add loads that remain fixed throughout an optimization. An example would be a fuel load.

Parameters:
  • timeStep (int) -- Time step index to apply load to.

  • compIDs (list[int] or int) -- The components with added loads. Use pyTACS.selectCompIDs method to determine this.

  • pressures (Numpy array length 1 or compIDs) -- Array of pressure values for each component

  • timeStage (int or None) -- Time stage index to apply load to. Default is None, which is applicable only for multi-step methods like BDF. For multi-stage methods like DIRK, this index must be specified.

  • faceIndex (int) -- Indicates which face (side) of element to apply pressure to. Note: not required for certain elements (i.e. shells)

addPressureToElements(timeStep, elemIDs, pressures, timeStage=None, faceIndex=0, nastranOrdering=False)[source]

This method is used to add a fixed presure to the selected element IDs at specified time instance. Pressures can be specified on an element by element basis (if pressures is an array) or set to a uniform value (if pressures is a scalar)

Parameters:
  • timeStep (int) -- Time step index to apply load to.

  • elemIDs (list[int]) -- The global element ID numbers for which to apply the pressure.

  • pressures (Numpy array length 1 or elemIDs) -- Array of pressure values for each element

  • timeStage (int or None) -- Time stage index to apply load to. Default is None, which is applicable only for multi-step methods like BDF. For multi-stage methods like DIRK, this index must be specified.

  • faceIndex (int) -- Indicates which face (side) of element to apply pressure to. Note: not required for certain elements (i.e. shells)

  • nastranOrdering (bool) -- Flag signaling whether elemIDs are in TACS (default) or NASTRAN ordering

addInertialLoad(timeStep, inertiaVector, timeStage=None)[source]

This method is used to add a fixed inertial load at a specified time step due to a uniform acceleration over the entire model. This is most commonly used to model gravity loads on a model.

Parameters:
  • timeStep (int) -- Time step index to apply load to.

  • inertiaVector (numpy.ndarray) -- Acceleration vector used to define inertial load.

  • timeStage (int or None) -- Time stage index to apply load to. Default is None, which is applicable only for multi-step methods like BDF. For multi-stage methods like DIRK, this index must be specified.

addCentrifugalLoad(timeStep, omegaVector, rotCenter, timeStage=None)[source]

This method is used to add a fixed centrifugal load at at specified time step due to a uniform rotational velocity over the entire model. This is most commonly used to model rotors, rolling aircraft, etc.

Parameters:
  • timeStep (int) -- Time step index to apply load to.

  • omegaVector (numpy.ndarray) -- Rotational velocity vector (rad/s) used to define centrifugal load.

  • rotCenter (numpy.ndarray) -- Location of center of rotation used to define centrifugal load.

  • timeStage (int or None) -- Time stage index to apply load to. Default is None, which is applicable only for multi-step methods like BDF. For multi-stage methods like DIRK, this index must be specified.

addLoadFromBDF(timeStep, loadID, timeStage=None, scale=1.0)[source]

This method is used to add a fixed load set defined in the BDF file to the problem at a specified time instance. Currently, only supports LOAD, FORCE, MOMENT, GRAV, RFORCE, PLOAD2, and PLOAD4.

Parameters:
  • timeStep (int) -- Time step index to apply load to.

  • loadID (int) -- Load identification number of load set in BDF file user wishes to add to problem.

  • timeStage (int or None) -- Time stage index to apply load to. Default is None, which is applicable only for multi-step methods like BDF. For multi-stage methods like DIRK, this index must be specified.

  • scale (float) -- Factor to scale the BDF loads by before adding to problem.

setInitConditions(vars=None, dvars=None, ddvars=None)[source]

Set the initial conditions associated with this problem

Parameters:
  • vars (float or numpy.ndarray or tacs.TACS.Vec) -- Initial conditions of the state variables

  • dvars (float or numpy.ndarray or tacs.TACS.Vec) -- Initial conditions of the first time-derivative of the state variables

  • ddvars (float or numpy.ndarray or tacs.TACS.Vec) -- Initial conditions of the second time-derivative of the state variables

solve()[source]

Solve the time integrated transient problem. The forces must already be set.

prepIterativeSolve()[source]

Prepare to solve the time integrated transient problem.

iterate(timeStep, timeStage=None, Fext=None)[source]

Iterate a single time instance in the time integrated transient problem. Useful for iterating an aeroelastic problem that is tightly coupled, where intermediate structural states need to be passed to an external fluid solver.

Requires the user to use an outer for-loop over the number of time steps/stages of the problem.

Parameters:
  • timeStep (int) -- Time step index to iterate

  • timeStage (int or None) -- Time stage index to iterate. If not None, solver must be multistage

  • Fext (tacs.TACS.Vec or numpy.ndarray or None) -- If Fext is not None, add this force vector to the loads applied at this time instance. Fext must be sized such that the flattened array is (numOwnedNodes*numVarsPerNode) in length. Useful for applying aerodynamic loads which change at every time instance in a tighly coupled aeroelastic solution.

Examples

>>> transientProblem.prepIterativeSolve()
>>> for step in range(transientProblem.numSteps + 1):
>>>     for stage in range(transientProblem.numStages):
>>>         # assume load array comes from an external CFD solver...
>>>         Fext = externalCFDSolver.solve(step, stage)
>>>         # iterate the transient problem in TACS
>>>         transientProblem.iterate(timeStep=step, timeStage=stage, Fext=Fext)
>>>         # get the necessary structural states for coupling
>>>         time, states, dstates, ddstates = transientProblem.getVariables(timeStep=step, timeStage=stage)
addFunction(funcName, funcHandle, compIDs=None, **kwargs)[source]

Generic method to add a function for TACS. It is intended to be reasonably generic since the user supplies the actual function handle to use. See the functions module for supported TACS eval functions.

Parameters:
  • funcName (str) -- The user-supplied name for the function. This will typically be a string that is meaningful to the user

  • funcHandle (tacs.TACS.Function) -- The function handle to use for creation. This must come from the functions module in tacs.

  • compIDs (list) -- List of compIDs to select.

  • **kwargs -- Any keyword arguments to be passed to the TACS function during setup.

evalFunctions(funcs, evalFuncs=None, ignoreMissing=False)[source]

This is the main routine for returning useful information from pytacs. The functions corresponding to the strings in evalFuncs are evaluated and updated into the provided dictionary.

Parameters:
  • funcs (dict) -- Dictionary into which the functions are saved.

  • evalFuncs (iterable object containing strings.) -- If not none, use these functions to evaluate.

  • ignoreMissing (bool) -- Flag to supress checking for a valid function. Please use this option with caution.

Examples

>>> funcs = {}
>>> transientProblem.solve()
>>> transientProblem.evalFunctions(funcs, ['mass'])
>>> funcs
>>> # Result will look like (if TransientProblem has name of 'c1'):
>>> # {'cl_mass':12354.10}
evalFunctionsSens(funcsSens, evalFuncs=None)[source]

This is the main routine for returning useful (sensitivity) information from problem. The derivatives of the functions corresponding to the strings in evalFuncs are evaluated and updated into the provided dictionary. The derivitives with respect to all design variables and node locations are computed.

Parameters:
  • funcsSens (dict) -- Dictionary into which the derivatives are saved.

  • evalFuncs (iterable object containing strings) -- The functions the user wants returned

Examples

>>> funcsSens = {}
>>> transientProblem.evalFunctionsSens(funcsSens, ['mass'])
>>> funcsSens
>>> # Result will look like (if TransientProblem has name of 'c1'):
>>> # {'c1_mass':{'struct':[1.234, ..., 7.89], 'Xpts':[3.14, ..., 1.59]}}
getVariables(timeStep, timeStage=None, states=None, dstates=None, ddstates=None)[source]

Return the current state values for the current problem

Parameters:
  • timeStep (int) -- Time step index to get state variables for.

  • timeStage (int or None) -- Time stage index to get state variables for.

  • states (tacs.TACS.Vec or numpy.ndarray or None) -- If states is not None, place the state variables into this array (optional).

  • dstates (tacs.TACS.Vec or numpy.ndarray or None) -- If dstates is not None, place the time derivative of the state variables into this array (optional).

  • ddstates (tacs.TACS.Vec or numpy.ndarray or None) -- If ddstates is not None, place the second time derivative of the state variables into this array (optional).

Returns:

  • time (float) -- The time at specified time instance

  • states (tacs.TACS.Vec or numpy.ndarray) -- The state variables.

  • dstates (tacs.TACS.Vec or numpy.ndarray or None) -- The time derivative of the state variables.

  • ddstates (tacs.TACS.Vec or numpy.ndarray or None) -- The second time derivative of the state variables.

zeroLoads()[source]

Zero all applied loads

writeSolution(outputDir=None, baseName=None, number=None, timeSteps=None)[source]

This is a generic shell function that writes the output file(s). The intent is that the user or calling program can call this function and pyTACS writes all the files that the user has defined. It is recommended that this function is used along with the associated logical flags in the options to determine the desired writing procedure

Parameters:
  • outputDir (str or None) -- Use the supplied output directory

  • baseName (str or None) -- Use this supplied string for the base filename. Typically only used from an external solver.

  • number (int or None) -- Use the user supplied number to index solution. Again, only typically used from an external solver

  • timeSteps (int or list[int] or None) -- Time step index or indices to get state variables for. If None, returns a solution for all time steps. Defaults to None.

dtype

alias of float64

getDesignVarRange()

get the lower/upper bounds for the design variables.

Returns:

  • xlb (numpy.ndarray) -- The design variable lower bound.

  • xub (numpy.ndarray) -- The design variable upper bound.

getDesignVars()

Get the current set of design variables for this problem.

Returns:

x -- The current design variable vector set in tacs.

Return type:

numpy.ndarray

getFunctionKeys()

Return a list of the current function key names

Returns:

funcNames -- List containing user-defined names for functions added so far.

Return type:

list[str]

getNodes()

Return the mesh coordinates of this problem.

Returns:

coords -- Structural coordinate in array of size (N * 3) where N is the number of structural nodes on this processor.

Return type:

numpy.ndarray

getNumCoordinates()

Return the number of mesh coordinates on this processor.

Returns:

ncoords -- Number of mesh coordinates on this processor.

Return type:

int

getNumDesignVars()

Return the number of design variables on this processor.

Returns:

ndvs -- Number of design variables on this processor.

Return type:

int

getNumOwnedNodes()

Get the number of nodes owned by this processor.

Returns:

nnodes -- Number of nodes on this processor.

Return type:

int

getNumVariables()

Return the number of degrees of freedom (states) that are on this processor

Returns:

nstate -- number of states.

Return type:

int

getOption(name)

Get a solver option value. The name is not case sensitive.

Parameters:

name (str) -- Name of option to get

getVarsPerNode()

Get the number of variables per node for the model.

Returns:

vpn -- Number of variables per node.

Return type:

int

property isNonlinear

The public interface for the isNonlinear attribute. Implemented as a property so that it is read-only.

classmethod printDefaultOptions()

Prints a nicely formatted dictionary of all the default solver options to the stdout

printModifiedOptions()

Prints a nicely formatted table of all the options that have been modified from their defaults

printOptions()

Prints a nicely formatted dictionary of all the current solver options to the stdout on the root processor

setDesignVars(x)

Update the design variables used by tacs.

Parameters:

x (numpy.ndarray or dict or tacs.TACS.Vec) -- The variables (typically from the optimizer) to set. It looks for variable in the self.varName attribute if in dict.

setNodes(Xpts)

Set the mesh coordinates of the structure.

Parameters:

Xpts (numpy.ndarray) -- Structural coordinate in array of size (N * 3) where N is the number of structural nodes on this processor.

setOptions(options)

Set multiple solver options at once. The names are not case sensitive.

Parameters:

options (dict) -- Dictionary of option names and values to set

setVarName(varName)

Set a name for the structural variables in pyOpt. Only needs to be changed if more than 1 pytacs object is used in an optimization

Parameters:

varName (str) -- Name of the structural variable used in addVarGroup().

writeSensFile(evalFuncs, tacsAim, proc: int = 0, root=0)

write an ESP/CAPS .sens file from the tacs aim Optional tacs_aim arg for TacsAim wrapper class object in root/tacs/caps2tacs/

Parameters:
  • evalFuncs (list[str]) -- names of TACS functions to be evaluated

  • tacsAim (tacs.caps2tacs.TacsAIM) -- class which handles the sensitivity file writing for ESP/CAPS shape derivatives

  • proc (int) -- which processor (in case of parallel tacsAIM instances) to write the sens file to

  • root (int) -- which processor we get the data from (usually proc 0)