BucklingProblem¶
The main purpose of this class is to represent all relevant information for a linearized buckling analysis.
Note
This class should be created using the
pyTACS.createBucklingProblem
method.
Options¶
Options can be set for BucklingProblem
at time of creation for the class in the
pyTACS.createBucklingProblem
method or using the
BucklingProblem.setOption
method. Current option values for a class
instance can be printed out using the BucklingProblem.printOption
method.
The following options, their default values and descriptions are listed below:
+----------------------------------------+
| BucklingProblem default options: |
+----------------------------------------+
'outputDir': './'
Output directory for F5 file writer.
'L2Convergence': 1e-12
Absolute convergence tolerance for Eigenvalue solver based on l2 norm of residual.
'L2ConvergenceRel': 1e-12
Relative convergence tolerance for Eigenvalue 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.
'subSpaceSize': 10
Subspace size for Krylov solver used by Eigenvalue solver.
'nRestarts': 15
Max number of resets for Krylov solver used by Eigenvalue solver.
'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 Eigenvalue solver.
Accepts:
0 : No printing.
1 : Print major iterations.
> 1 : Print major + minor iterations.
API Reference¶
- class tacs.problems.BucklingProblem(name, sigma, numEigs, assembler, comm, outputViewer=None, meshLoader=None, isNonlinear=False, options=None)[source]¶
NOTE: This class should not be initialized directly by the user. Use pyTACS.createBucklingProblem instead.
- Parameters:
name (str) -- Name of this tacs problem
sigma (float) -- Guess for the lowest eigenvalue. This corresponds to the lowest buckling load factor.
numEigs (int) -- Number of eigenvalues to solve for
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
- setValName(valName)[source]¶
Set a name for the eigenvalues. Only needs to be changed if more than 1 pytacs object is used in an optimization
- Parameters:
valName (str) -- Name of the eigenvalue output used in evalFunctions().
- getNumEigs()[source]¶
Get the number of eigenvalues requested from solver for this problem.
- Returns:
numEigs -- Number of eigenvalues.
- Return type:
- evalFunctions(funcs, evalFuncs=None, ignoreMissing=False)[source]¶
Evaluate eigenvalues for problem. The functions corresponding to the integers in evalFuncs are evaluated and updated into the provided dictionary.
- Parameters:
Examples
>>> funcs = {} >>> bucklingProblem.solve() >>> bucklingProblem.evalFunctions(funcs, 'eigsm.0') >>> funcs >>> # Result will look like (if bucklingProblem has name of 'c1'): >>> # {'c1_eigsm.0':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 = {} >>> bucklingProblem.evalFunctionsSens(funcsSens, 'eigsm.0') >>> funcsSens >>> # Result will look like (if bucklingProblem has name of 'c1'): >>> # {'c1_eigsm.0':{'struct':[1.234, ..., 7.89], 'Xpts':[3.14, ..., 1.59]}}
- 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:
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.
- solve(Fext=None, u0=None)[source]¶
Solve the eigenvalue problem. The forces must already be set.
- Parameters:
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.
u0 (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. Alternate to Fext.
- getVariables(index, states=None)[source]¶
Return the current state values for one mode of the current problem
- Parameters:
index (int) -- Mode index to return solution for.
states (tacs.TACS.Vec or numpy.ndarray or None) -- Place eigenvector for mode into this array (optional).
- Returns:
eigVal (float) -- Eigenvalue for mode corresponds to buckling load factor
states (numpy.ndarray) -- Eigenvector for mode
- getModalError(index)[source]¶
Return the error associated with a particular mode
- Parameters:
index (int) -- Mode index to return solution for
- addXptSens(indices, xptSensList, scale=1.0)[source]¶
Add partial sensitivity contribution due to nodal coordinates for eigenvalue
- addDVSens(indices, dvSensList, scale=1.0)[source]¶
Add partial sensitivity contribution due to design variables for eigenvalue
- evalSVSens(indices, svSensList)[source]¶
Add partial sensitivity contribution due to state variables for eigenvalue
- Parameters:
indices (list[int]) -- Mode indices to return sensitivity for.
svSensList (list[tacs.TACS.Vec] or list[numpy.ndarray]) -- List of sensitivity vectors to add partial sensitivity to
- writeSolution(outputDir=None, baseName=None, number=None, indices=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
indices (int or list[int] or None) -- Mode index or indices to get state variables for. If None, returns a solution for all modes. 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:
- getFunctionKeys()¶
Return a list of the current function key names
- 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:
- getNumCoordinates()¶
Return the number of mesh coordinates on this processor.
- Returns:
ncoords -- Number of mesh coordinates on this processor.
- Return type:
- getNumDesignVars()¶
Return the number of design variables on this processor.
- Returns:
ndvs -- Number of design variables on this processor.
- Return type:
- getNumOwnedNodes()¶
Get the number of nodes owned by this processor.
- Returns:
nnodes -- Number of nodes on this processor.
- Return type:
- getNumVariables()¶
Return the number of degrees of freedom (states) that are on this processor
- Returns:
nstate -- number of states.
- Return type:
- 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:
- 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)