StaticProblem

The main purpose of this class is to represent all relevant information for a static 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.createStaticProblem method.

Options

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

+----------------------------------------+
|     StaticProblem default options:     |
+----------------------------------------+
'outputDir': './'
	 Output directory for F5 file writer.
'linearSolver': 'GMRES'
	 Krylov subspace method to use for linear solver. Currently only supports 'GMRES'
'nonlinearSolver': 'Continuation'
	 Convergence method to use for nonlinear solver. Currently only supports 'Continuation'
'orderingType': 3
	 Ordering type to use for matrix partitioning.
	 Acceptable values are:
		 tacs.TACS.NATURAL_ORDER = 0
		 tacs.TACS.RCM_ORDER = 1
		 tacs.TACS.ND_ORDER = 3
		 tacs.TACS.TACS_AMD_ORDER = 4
		 tacs.TACS.MULTICOLOR_ORDER = 5
'PCFillLevel': 1000
	 Preconditioner fill level.
'PCFillRatio': 20.0
	 Preconditioner fill ratio.
'subSpaceSize': 10
	 Subspace size for Krylov solver.
'nRestarts': 15
	 Max number of restarts for Krylov solver.
'flexible': True
	 Flag for whether the preconditioner is flexible.
'L2Convergence': 1e-12
	 Absolute convergence tolerance for linear solver based on l2 norm of residual.
'L2ConvergenceRel': 1e-12
	 Relative convergence tolerance for linear solver 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.
'useMonitor': False
	 Flag for whether to attach a debug monitor to the linear solver.
'monitorFrequency': 10
	 Print frequency for sub iterations of linear solver.
'writeSolution': True
	 Flag for suppressing all f5 file writing.
'writeNLIterSolutions': False
	 Flag to save a solution file at every nonlinear solver iterations.
'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 nonlinear solver.
	 Accepts:
		   0 : No printing.
		   1 : Print nonlinear iterations.

Nonlinear solvers

When you create an assembler using a nonlinear element type or constitutive model, any static problems you create will automatically setup a nonlinear solver. A continuation solver is used to control an adaptive load incrementation process, and a Newton solver is used to solve the nonlinear system of equations at each load increment. The options for these solvers should be set directly to problem.nonlinearSolver and problem.nonlinearSolver.innerSolver. See ContinuationSolver and NewtonSolver for more information.

API Reference

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

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

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

  • assembler (tacs.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.TACS.TACSToFH5) -- Cython object used to write out f5 files that can be converted and used for postprocessing.

  • meshLoader (tacs.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

property loadScale

This is a scaling factor applied to all forcing terms

Forcing terms includes both the user supplied force vector and the forcing terms coming from aux elements in the TACS assembler (e.g inertial, centrifugal forces)

Returns:

The current load scale

Return type:

float or complex

setLoadScale(value)[source]

Set the scaling applied to external loads

Parameters:

value (float or complex) -- Value to set the load scale to

getLoadScale()[source]

Get the current load scale

This function is required by the continuation solver

Returns:

Current load scale

Return type:

float

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.

setDesignVars(x)[source]

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.

setNodes(Xpts)[source]

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.

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

This method is used to add a FIXED TOTAL LOAD on one or more components, defined by COMPIDs. 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:
  • compIDs (list[int] or int) -- The components with added loads. Use pyTACS selectCompIDs method to determine this.

  • F (numpy.ndarray 1d or 2d 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.

  • 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(nodeIDs, F, nastranOrdering=False)[source]

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

Parameters:
  • 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.

  • 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(Fapplied)[source]

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

K*u = f

Where:
  • K : Stiffness matrix for problem

  • u : State variables for problem

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

Parameters:

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

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

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

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

  • tractions (Numpy array length 1 or compIDs) -- Array of traction vectors for each component

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

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

This method is used to add a fixed traction to the selected element IDs. 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:
  • elemIDs (list[int]) -- The global element ID numbers for which to apply the traction.

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

  • 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(compIDs, pressures, faceIndex=0)[source]

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

Parameters:
  • 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

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

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

This method is used to add a fixed presure to the selected element IDs. 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:
  • 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

  • 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(inertiaVector)[source]

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

Parameters:

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

addCentrifugalLoad(omegaVector, rotCenter, firstOrder=False)[source]

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

Parameters:
  • 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.

  • firstOrder (bool, optional) -- Whether to use first order approximation for centrifugal load, which computes the force in the displaced position. By default False

addLoadFromBDF(loadID, scale=1.0)[source]

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

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

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

solve(Fext=None)[source]

Solution of the static problem for current load set. The forces must already be set.

Parameters:

Fext (numpy.ndarray or tacs.TACS.Vec, optional) -- Distributed array containing additional loads (ex. aerodynamic forces for aerostructural coupling) to applied to RHS of the static problem, by default None

updateJacobian(res=None)[source]

Update the Jacobian (a.k.a stiffness) matrix

The Jacobian will only actually be updated if the _jacobianUpdateRequired flag is set to True.

Parameters:

res (tacs.TACS.Vec, optional) -- If provided, the residual is also computed and stored in this vector

updatePreconditioner()[source]

Update the Jacobian (a.k.a stiffness) matrix preconditioner

By default, the static problem uses a full LU factorization of the Jacobian matrix as a preconditioner, which means that this step is typically the most expensive operation in the solution process.

The preconditioner will only actually be updated if the _preconditionerUpdateRequired flag is set to True. This occurs whenever the Jacobian is updated.

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 = {}
>>> staticProblem.solve()
>>> staticProblem.evalFunctions(funcs, ['mass'])
>>> funcs
>>> # Result will look like (if StaticProblem 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 = {}
>>> staticProblem.evalFunctionsSens(funcsSens, ['mass'])
>>> funcsSens
>>> # Result will look like (if StaticProblem has name of 'c1'):
>>> # {'c1_mass':{'struct':[1.234, ..., 7.89], 'Xpts':[3.14, ..., 1.59]}}
addSVSens(evalFuncs, svSensList)[source]

Add the state variable partial sensitivity to the ADjoint RHS for given evalFuncs

Parameters:
  • evalFuncs (list[str]) -- The functions the user wants returned

  • svSensList (list[tacs.TACS.Vec] or list[numpy.ndarray]) -- List of sensitivity vectors to add partial sensitivity to

addDVSens(evalFuncs, dvSensList, scale=1.0)[source]

Add partial sensitivity contribution due to design vars for evalFuncs

Parameters:
  • evalFuncs (list[str]) -- The functions the user wants returned

  • dvSensList (list[tacs.TACS.Vec] or list[numpy.ndarray]) -- List of sensitivity vectors to add partial sensitivity to

  • scale (float) -- Scalar to multiply partial sensitivity by. Defaults to 1.0

addAdjointResProducts(adjointlist, dvSensList, scale=-1.0)[source]

Add the adjoint product contribution to the design variable sensitivity arrays

Parameters:
  • adjointlist (list[tacs.TACS.Vec] or list[numpy.ndarray]) -- List of adjoint vectors for residual sensitivity product

  • dvSensList (list[tacs.TACS.Vec] or list[numpy.ndarray]) -- List of sensitivity vectors to add product to

  • scale (float) -- Scalar to multiply product by. Defaults to -1.0

addXptSens(evalFuncs, xptSensList, scale=1.0)[source]

Add partial sensitivity contribution due to nodal coordinates for evalFuncs

Parameters:
  • evalFuncs (list[str]) -- The functions the user wants returned

  • xptSensList (list[tacs.TACS.Vec] or list[numpy.ndarray]) -- List of sensitivity vectors to add partial sensitivity to

  • scale (float) -- Scalar to multiply partial sensitivity by. Defaults to 1.0

addAdjointResXptSensProducts(adjointlist, xptSensList, scale=-1.0)[source]

Add the adjoint product contribution to the nodal coordinates sensitivity arrays

Parameters:
  • adjointlist (list[tacs.TACS.Vec] or list[numpy.ndarray]) -- List of adjoint vectors for residual sensitivity product

  • xptSensList (list[tacs.TACS.Vec] or list[numpy.ndarray]) -- List of sensitivity vectors to add product to

  • scale (float) -- Scalar to multiply product by. Defaults to -1.0

getResidual(res, Fext=None)[source]

This routine is used to evaluate directly the structural residual. Only typically used with aerostructural analysis.

Parameters:
  • res (tacs.TACS.Vec or numpy.ndarray) -- If res is not None, place the residuals into this array.

  • Fext (tacs.TACS.Vec or numpy.ndarray, optional) -- Distributed array containing additional loads (ex. aerodynamic forces for aerostructural coupling) to applied to RHS of the static problem.

getForces(externalForceVec, internalForceVec, Fext=None)[source]

Compute the internal and external forces acting on the structure

The computations here are based on the residual equation:

r(u, loadScale) = -(Fint(u) + loadScale * Fext(u))

Thus, the internal forces are given by:

Fint(u) = -r(u, 0)

And the external forces are given by:

Fext(u) = -r(u, 1) - Fi(u)

Parameters:
  • externalForceVec (TACS BVec or numpy array) -- Vector/array to store external forces in

  • internalForceVec (TACS BVec or numpy array) -- Vector/array to store internal forces in

  • Fext (TACS BVec or numpy array, optional) -- Distributed array containing additional loads (ex. aerodynamic forces for aerostructural coupling) to applied to RHS of the static problem.

getJacobian()[source]

Get the problem's Jacobian in sciPy sparse matrix format

Returns:

K -- A tuple of 2 scipy.sparse.bsr_matrices (A, B) if Jacobian is a TACSParallelMat, or 4 scipy.sparse.bsr_matrices (A, B, C, D) if Jacobian is a TACSSchurMat

Return type:

(scipy.sparse.bsr_matrix, scipy.sparse.bsr_matrix) or (scipy.sparse.bsr_matrix, scipy.sparse.bsr_matrix, scipy.sparse.bsr_matrix, scipy.sparse.bsr_matrix)

addTransposeJacVecProduct(phi, prod, scale=1.0)[source]

Adds product of transpose Jacobian and input vector into output vector as shown below: prod += scale * J^T . phi

Parameters:
  • phi (tacs.TACS.Vec or numpy.ndarray) -- Input vector to product with the transpose Jacobian.

  • prod (tacs.TACS.Vec or numpy.ndarray) -- Output vector to add Jacobian product to.

  • scale (float) -- Scalar used to scale Jacobian product by.

zeroVariables()[source]

Zero all the tacs solution b-vecs

zeroLoads()[source]

Zero all applied loads

solveAdjoint(rhs, phi)[source]

Solve the structural adjoint.

Parameters:
  • rhs (tacs.TACS.Vec or numpy.ndarray) -- right hand side vector for adjoint solve

  • phi (tacs.TACS.Vec or numpy.ndarray) -- BVec or numpy array into which the adjoint is saved

getVariables(states=None)[source]

Return the current state values for the problem

Parameters:

states (tacs.TACS.Vec or numpy.ndarray) -- Vector to place current state variables into (optional)

Returns:

states -- current state vector

Return type:

numpy.ndarray

setVariables(states)[source]

Set the structural states for current load case.

Parameters:

states (numpy.ndarray) -- Values to set. Must be the size of getNumVariables()

getOutputFileName(outputDir=None, baseName=None, number=None)[source]

Figure out a base path/name for output files

Parameters:
  • outputDir (str, optional) -- Directory to write file to, by default uses the 'outputDir' option

  • baseName (str, optional) -- Case name, by default uses self.name

  • number (int, optional) -- A number to append to the filename, by default uses the call-count of the problem

Returns:

Full path to output file (excluding any extension)

Return type:

str

writeSolution(outputDir=None, baseName=None, number=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

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

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)

writeSolutionHistory(outputDir=None, baseName=None, number=None)[source]

Write the nonlinear solver history to a file

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

writeLoadToBDF(bdfFile, loadCaseID)[source]

Write loads from problem to NASTRAN BDF file. NOTE: To get correct loads, solve method should be called before this method.

Parameters:
  • bdfFile (str or pyNastran.bdf.bdf.BDF or None) -- Name of file to write BDF file to. Only required on root proc, can be None otherwise.

  • loadCaseID (int) -- NASTARAN loadcase ID to assign loads to in BDF.