Direct¶
In the direct Python interface for TACS, the user will be responsible for setting up and bookkeeping of most TACS objects. This includes objects like: the MeshLoader, Assembler, state vectors, FE matrices, etc. This approach allows for more visibility into the workings of the C++ source code, but can be daunting for new users. For a more simplified user interface where most of this setup and bookkeeping has been automated for the user, see pyTACS.
Workflow¶
The most common usage of TACS is to evaluate the values and gradients of desired structural functions with respect to specified design variables. Using the direct interface, this workflow proceeds as follows:
Load in a finite element model of the desired structure (in the form of a NASTRAN-style file) using an instance of the
MeshLoader
class.For each component of the loaded model, generate an element with the desired constitutive properties and design variables.
Create an instance of the
Assembler
class and apply boundary conditions.Solve the system and evaluate the functions and their gradients with respect to the design variables.
These function values and gradients can then be passed to an optimizer (such as ParOpt
)
in order to minimize the value of a particular function subject to some constraints.
Improved design variable values are iteratively computed by the optimizer and Step 4 is
repeated until the optimization criteria are satisfied.
Assembler¶
The Assembler
object can be created using the MeshLoader
or Creator
classes. It contains the methods that allow the user to solve the finite element problem and evaluate the desired structural functions and their gradients. Once an instance has been created, its typical usage is as follows:
Apply loads to the model with the
applyBCs()
.Create the solver using the
createVec()
andcreateFEMat()
functions.Use the
setDesignVars()
function to set design variable vales.Evaluate structural functions (e.g. Structural Mass, KSFailure) using the
evalFunctions()
call.The gradients of the functions with respect to the design variables can be evaluated using the adjoint method with the
addDVSens()
,addSVSens()
, andaddAdjointResProducts()
functions.
- class TACS.Assembler¶
- addAdjointResProducts(adjlist, dfdxlist, alpha=1.0, loadScale=1.0)¶
This function is collective on all TACSAssembler processes. This computes the product of the derivative of the residual w.r.t. the design variables with several adjoint vectors simultaneously. This saves computational time as the derivative of the element residuals can be reused for each adjoint vector. This function performs the same task as evalAdjointResProduct, but uses more memory than calling it for each adjoint vector.
adjoint: the array of adjoint vectors dvSens: the product of the derivative of the residuals and the adjoint num_dvs: the number of design variables loadScale: Scaling factor for the aux element contributions, by default 1
- addAdjointResXptSensProducts(adjlist, dfdXlist, alpha=1.0, loadScale=1.0)¶
This function is collective on all TACSAssembler processes. This computes the product of the derivative of the residual w.r.t. the node locations with several adjoint vectors simultaneously.
- addDVSens(funclist, dfdxlist, alpha=1.0)¶
Evaluate the derivative of a list of functions w.r.t. the design variables.
- addJacobianVecProduct(scale, alpha, beta, gamma, x, y, matOr=0, loadScale=1.0)¶
Compute the Jacobian-vector product
- addMatDVSensInnerProduct(scale, matType, psi, phi, dfdx)¶
Add the derivative of the inner product of the specified matrix with the input vectors to the design variable sensitivity vector A.
- addSVSens(funclist, dfdulist, alpha=1.0, beta=0.0, gamma=0.0)¶
Evaluate the derivative of the function w.r.t. the state variables.
function: the function pointer
vec: the derivative of the function w.r.t. the state variables
- addXptSens(funclist, dfdXlist, alpha=1.0)¶
Evaluate the derivative of a list of functions w.r.t. the node locations
- applyBCs(vec)¶
Apply boundary conditions to the vector
- applyMatBCs(mat)¶
Apply boundary conditions to the matrix
- assembleJacobian(alpha, beta, gamma, residual, A, matOr=0, loadScale=1.0)¶
Assemble the Jacobian matrix
This function assembles the global Jacobian matrix and residual. This Jacobian includes the contributions from all elements. The Dirichlet boundary conditions are applied to the matrix by zeroing the rows of the matrix associated with a boundary condition, and setting the diagonal to unity. The matrix assembly also performs any communication required so that the matrix can be used immediately after assembly.
alpha: coefficient on the variables beta: coefficient on the time-derivative terms gamma: coefficient on the second time derivative term residual: the residual of the governing equations A: the Jacobian matrix matOr: the matrix orientation NORMAL or TRANSPOSE loadScale: Scaling factor for the aux element contributions, by default 1
- assembleMatCombo(matType1, scale1, matType2, scale2, A, matOr=0, loadScale=1.0)¶
Assemble a combination of two matrices
- assembleMatType(matType, A, matOr=0, loadScale=1.0)¶
Assemble the Jacobian matrix
This function assembles the global Jacobian matrix and residual. This Jacobian includes the contributions from all elements. The Dirichlet boundary conditions are applied to the matrix by zeroing the rows of the matrix associated with a boundary condition, and setting the diagonal to unity. The matrix assembly also performs any communication required so that the matrix can be used immediately after assembly.
residual: the residual of the governing equations A: the Jacobian matrix alpha: coefficient on the variables beta: coefficient on the time-derivative terms gamma: coefficient on the second time derivative term matOr: the matrix orientation NORMAL or TRANSPOSE loadScale: Scaling factor for the aux element contributions, by default 1
- assembleRes(residual, loadScale=1.0)¶
Assemble the residual associated with the input load case.
This residual includes the contributions from element tractions set in the TACSSurfaceTraction class and any point loads. Note that the vector entries are zeroed first, and that the Dirichlet boundary conditions are applied after the assembly of the residual is complete.
rhs: the residual output loadScale: Scaling factor for the aux element contributions, by default 1
- computeReordering(order_type, mat_type)¶
Compute a reordering of the unknowns before initialize()
- copyVariables(vec=None, dvec=None, ddvec=None)¶
Set the values of the state variables
- static create(comm, varsPerNode, numOwnedNodes, numElements, numDependentNodes=0)¶
Static factory method for creating an instance of Assembler
- createDesignVec()¶
Create a distribute design variable vector
- createMat()¶
Create a distributed matrix
- createNodeVec()¶
Create a distributed node vector
- createSchurMat(order_type=4)¶
Create a parallel matrix specially suited for finite-element analysis.
On the first call, this computes a reordering with the scheme provided. On subsequent calls, the reordering scheme is reused s that all FEMats, created from the same TACSAssembler object have the same non-zero structure. This makes adding matrices together easier (which is required for eigenvalue computations.)
The first step is to determine the coupling nodes. (For a serial case there are no coupling nodes, so this is very simple!) Then, the nodes that are not coupled to other processes are determined. The coupling and non-coupling nodes are ordered separately. The coupling nodes must be ordered at the end of the block, while the local nodes must be ordered first. This type of constraint is not usually imposed in matrix ordering routines, so here we use a kludge. First, order all the nodes and determine the ordering of the coupling variables within the full set. Next, order the local nodes. Tis hopefully reduces the fill-ins required, although there is no firm proof to back that up.
The results from the reordering are placed in a set of objects. The matrix reordering is stored in feMatBIndices and feMatCIndices while two mapping objects are created that map the variables from the global vector to reordered matrix.
Mathematically this reordering can be written as follows,
A1 = (P A P^{T})
where P^{T} is a permutation of the columns (variables), while P is a permutation of the rows (equations).
- createVec()¶
Create a distributed vector.
Vector classes initialized by one TACS object, cannot be used by a second, unless they share are exactly the parallel layout.
- evalEnergies()¶
Evaluate the kinetic and potential energies
- evalFunctions(funclist)¶
Evaluate a list of TACS function
- getBcMap()¶
Get boundary conditions.
- getDesignVarRange(lb, ub)¶
Retrieve the design variable range.
This call is collective on all TACS processes. The ranges provided by individual objects may not be consistent (if someone provided incorrect data they could be.) Make a best guess; take the minimum upper bound and the maximum lower bound.
lowerBound: the lower bound on the design variables (output) upperBound: the upper bound on the design variables (output) numDVs: the number of design variables
- getDesignVars(x)¶
Collect all the design variable values assigned by this process
This code does not ensure consistency of the design variable values between processes. If the values of the design variables are inconsistent to begin with, the maximum design variable value is returned. Call setDesignVars to make them consistent.
Each process contains objects that maintain their own design variable values. Ensuring the consistency of the ordering is up to the user. Having multiply-defined design variable numbers corresponding to different design variables results in undefined behaviour.
- getElementData(num)¶
Return the element data associated with the element
- getElementNodes(num)¶
Get the node numbers associated with the given element
- getElements()¶
Get the elements
- getInitConditions(vec=None, dvec=None, ddvec=None)¶
Retrieve the initial conditions
- getMPIComm()¶
Get the MPI communicator
- getNodes(X)¶
Get the node locations
- getNumDependentNodes()¶
Return the number of dependent nodes
- getNumElements()¶
Return the number of elements
- getNumNodes()¶
Return the number of nodes in the TACSAssembler
- getNumOwnedNodes()¶
Return the number of owned nodes
- getOwnerRange()¶
Get the ranges of global node numbers owned by each processor
- getReordering()¶
Get the reordering
- getSimulationTime()¶
Retrieve the simulation time from TACS
- getVariables(vec=None, dvec=None, ddvec=None)¶
Set the values of the state variables
- getVarsPerNode()¶
Return the number of variables per node
- initialize()¶
Function to call after all the nodes and elements have been added into the created instance of TACS. This function need not be called when tacs is created using TACSCreator class.
- reorderVec(vec)¶
Reorder the vector based on the TACSAssembler reordering
- setAuxElements(elems=None)¶
Set the auxiliary elements
- setBCValuesFromVec(vec)¶
Set new Dirichlet BC values at nodes where BCs are imposed
This takes the new boundary condition values as the entries in the given vector where the Dirichlet boundary conditions are imposed.
- setBCs(vec)¶
Apply the Dirichlet boundary conditions to the state vector
- setDependentNodes(ptr, conn, weights)¶
Set the dependent node connectivity
- setDesignVars(x)¶
Set the design variables.
The design variable values provided must be the same on all processes for consistency. This call however, is not collective.
- setElementConnectivity(ptr, conn)¶
Set the connectivity
- setElements(elements)¶
Set the elements in to TACSAssembler
- setInitConditions(vec=None, dvec=None, ddvec=None)¶
Set the initial conditions as a vector
- setNodes(X)¶
Set the node locations
- setNumThreads(t)¶
Set the number of threads to use in computation
- setSimulationTime(time)¶
Set the simulation time within TACS
- setVariables(vec=None, dvec=None, ddvec=None)¶
Set the values of the state variables
- testElement(elemNum, print_level, dh=1e-06, rtol=1e-08, atol=0.1)¶
Test the implementation of the given element number.
This tests the stiffness matrix and various parts of the design-sensitivities: the derivative of the determinant of the Jacobian, the derivative of the strain w.r.t. the nodal coordinates, and the state variables and the derivative of the residual w.r.t. the design variables and nodal coordinates.
elemNum: the element number to test print_level: the print level to use
- testFunction(func, dh)¶
Test the implementation of the function.
This tests the state variable sensitivities and the design variable sensitivities of the function of interest. These sensitivities are computed based on a random perturbation of the input values. Note that a system of equations should be solved - or the variables should be set randomly before calling this function, otherwise this function may produce unrealistic function values.
Note that this function uses a central difference if the real code is compiled, and a complex step approximation if the complex version of the code is used.
func: the function to test num_dvs: the number of design variables to use dh: the step size to use
- zeroDDotVariables()¶
Zero the values of the 2nd time-derivatives of the state variables
- zeroDotVariables()¶
Zero the values of the time-derivatives of the state variables
- zeroVariables()¶
Zero the entries of the local variables
MeshLoader¶
MeshLoader
is an interface for reading in FEM data from
NASTRAN-style files (such as .bdf files).
The typical usage for a MeshLoader
class is as follows:
Create the object and call
scanBDFFile()
on the desired NASTRAN .bdf file.Retrieve the number of components using
getNumComponents()
Iterate through each component, create elements with desired constitutive properties and design variables. Set elements into the
MeshLoader
object and create theAssembler
object.
- class TACS.MeshLoader¶
- addAuxElement(aux, comp_num, elem)¶
Add the auxiliary element to the given component
- addFunctionDomain(func, comp_list)¶
Add the specified components to the domain of the function
- createTACS(varsPerNode, order_type=0, mat_type=2)¶
Create a distribtued version of TACS
- getBCs()¶
Return the boundary conditions associated with the file
- getComponentDescript(comp_num)¶
Return the component description
- getConnectivity()¶
Return the connectivity of the mesh
- getElementDescript(comp_num)¶
Retrieve the element description corresponding to the component number
- getNumComponents()¶
Return the number of components
- scanBDFFile(fname)¶
This scans a Nastran file - only scanning in information from the bulk data section
The only entries scanned are the entries beginning with elem_types and any GRID/GRID* entries
- setElement(comp_num, elem)¶
Set the element associated with a given component number
Creator¶
The Creator
object is similar to the MeshLoader
object, but
sets the nodes, elements, boundary conditions, etc. manually rather than loading them
from a NASTRAN-style file. This involves the use of the setNodes()
,
setElements()
, and setBoundaryConditions()
functions,
and finally the createTACS()
function which creates the Assembler
object.
- class TACS.Creator¶
- getElementPartition()¶
Retrieve the element partition
- setBoundaryConditions(nodes, bcptr=None, bcvars=None, bcvals=None)¶
Set the boundary conditions
- setElements(elements)¶
Set the elements
- setGlobalConnectivity(num_nodes, node_ptr, node_conn, id_nums)¶
Set the connectivity and the element id numbers
FrequencyAnalysis¶
The FrequencyAnalysis
object solves the natural frequency eigenproblem
and extracts the eigenvalues and eigenvectors (natural frequencies and mode shapes).
This could be used to, for example, minimize the mass of a beam by varying element
thicknesses with a lower bound constraint on its lowest natural frequency and an upper
bound constraint on the KSFailure function.
Integrator¶
The Integrator
class contains functions for solving the adjoint equations and governing equations forward in time. Classes for BDF, DIRK, ABM, and NBG integration inherit from this class.
- class TACS.Integrator¶
Class containing functions for solving the equations forward in time and adjoint.
- checkGradients(self, double dh)¶
Performs a FD/CSD verification of the gradients
- evalFunctions(self, funclist)¶
Evaluate a list of TACS function in time
- getAdjoint(self, int step_num, int func_num)¶
Get the adjoint vector at the given step
- getGradient(self, int func_num)¶
Get the time-dependent derivative of functionals
- getNumTimeSteps(self)¶
Get the number of time steps
- getStates(self, int time_step)¶
TACS state vectors are returned at the given time step
- getXptGradient(self, int func_num)¶
Get the time-dependent nodal derivatives of the functional
- initAdjoint(step_num)¶
initAdjoint(self, int step_num):
Initialize adjoint at the specified step
- integrate(self)¶
Integrates the governing equations forward in time
- integrateAdjoint(self)¶
Integrates the adjoint backwards in time
- iterate(self, int step_num, Vec forces=None)¶
Solve the nonlinear system at current time step
- iterateAdjoint(step_num, adjlist=None)¶
iterateAdjoint(self, int step_num, list adjlist=None):
Perform one iteration in reverse mode
- loadStates(step_num, prefix='')¶
Loads the states variables to disk. The string argument prefix can be used to put the binaries in a separate directory.
- persistStates(step_num, prefix='')¶
Writes the states variables to disk. The string argument prefix can be used to put the binaries in a separate directory.
- postAdjoint(self, int step_num)¶
Finish the calculations at the specified adjoint step
- setAbsTol(self, double atol)¶
Absolute tolerance of Newton solver
- setFH5(self, ToFH5 f5)¶
Configure the export of rigid bodies
- setFunctions(funcs, start_plane=-1, end_plane=-1)¶
Sets the functions for obtaining the derivatives.
- setInitNewtonDeltaFraction(self, double frac)¶
Parameter for globalization in Newton solver
- setJacAssemblyFreq(self, int freq)¶
How frequent to assemble the Jacobian for nonlinear solve
- setKrylovSubspaceMethod(self, KSM ksm)¶
Make TACS use this linear solver
- setMaxNewtonIters(max_newton_iters)¶
setMaxNewtonIters(self, int max_newton_iters):
Maximum iteration in Newton solver
- setOutputFrequency(self, int write_freq=0)¶
Configure how frequent to write f5 files
- setOutputPrefix(self, _prefix)¶
Output directory to use for f5 files
- setPrintLevel(self, int print_level, fname=None)¶
Level of print from TACSIntegrator 0: off 1: summary each step 2: summary each newton iteration
- setRelTol(self, double rtol)¶
Set the relative tolerance of Newton solver
- setTimeInterval(tinit, tfinal)¶
setTimeInterval(self, double tinit, double tfinal):
Set the time interval for the simulation
- setUseLapack(self, use_lapack)¶
Should TACSIntegrator use lapack for linear solve. This will be slow and need to be in serial mode only
- setUseSchurMat(self, int use_schur_mat, OrderingType order_type)¶
Use the TACSSchurMat for parallel execution