MACH

The MACH interface provides a specialized wrapper for coupling TACS with MDOLab's MACH framework codes, particularly for aerostructural analysis and multidisciplinary optimization. It is built around the StructProblem class, which adapts a TACS StaticProblem to conform to the baseclasses.StructProblem interface expected by MACH framework codes.

Overview

The MACH interface is designed for scenarios where TACS is integrated with external computational libraries, such as:

  • ADflow — aerodynamic solver for aerostructural coupling

  • pyGeo — geometric parametrization and design variables

  • pyoptsparse — gradient-based optimization

The typical workflow is:

  1. Create a pyTACS assembler and a StaticProblem.

  2. Set up any external forces on the StaticProblem.

  3. Wrap the StaticProblem with StructProblem.

  4. Optionally attach a DVGeometry object for geometric design variables.

  5. At each optimization step: call solve(), then call evalFunctions() / evalFunctionsSens().

Key Features

  • External Force Integration: The Fext property holds the aerodynamic coupling force vector. It is automatically included in the residual during solve().

  • Aitken Acceleration: solve() optionally applies Aitken relaxation to stabilize fixed-point iteration in aerostructural coupling loops.

  • Adjoint Capabilities: The class manages all adjoint state vectors needed for coupled sensitivity analysis, including phi, pLoad, dSdu, and adjRHS.

  • Geometric Design Variables: When a DVGeometry object is provided (via setDVGeo()), node coordinates are automatically updated before each analysis call and geometric sensitivities are folded in by evalFunctionsSens().

  • pyoptsparse Integration: addVariablesPyOpt() and addConstraintsPyOpt() register structural design variables and constraints directly with a pyoptsparse optimization problem.

  • Force File I/O: External coupling forces can be saved to and loaded from BDF-format files for post-processing or restarting aerostructural analyses.

Managed State Vectors

During coupled computation, both in forward and reverse mode, the following state vectors are managed by StructProblem:

Property

Description

Fext

External aerodynamic coupling force. Set by the aero solver before each solve() call.

phi

Structural adjoint vector \(\phi\). Updated by solveAdjoint().

pLoad

Aerodynamic adjoint contribution to the structural adjoint RHS: \([dA/du]^T \psi\). Supplied by the aero solver.

dSdu

Product \([dS/du]^T \phi\) from a coupled structural block. Supplied by the coupling layer.

adjRHS

Assembled adjoint RHS: \(dI/du - [dA/du]^T\psi - [dS/du]^T\phi\). Built by assembleAdjointRHS().

dIdu

Partial derivative of the objective with respect to state variables. Set by setAdjointRHS() and BCs are applied automatically.

matVecRHS / matVecSolve

Scratch vectors used for matrix–vector products in the coupled direct/adjoint system.

StructProblem Class

class tacs.mach.struct_problem.StructProblem(staticProblem, FEAAssembler, DVGeo=None, loadFile=None)[source]

MACH StructProblem wrapper for pyTACS StaticProblem class

Initialize the StructProblem.

Parameters:
  • staticProblem (tacs.problems.StaticProblem) -- StaticProblem object used for modeling and solving static cases.

  • FEAAssembler (tacs.pytacs.pyTACS) -- pyTACS assembler object.

  • DVGeo (pygeo.DVGeometry, optional) -- Object responsible for manipulating the geometry design variables. If None, no geometric design variables will be used.

  • loadFile (str, optional) -- Filename of the (static) external load file. Should be generated from pyAerostructure.

property name

Get the name of the structural problem.

Returns:

Name of the structural problem.

Return type:

str

property callCounter: int

Get the call counter for the underlying static problem.

Returns:

Number of times the static problem's solve method has been called.

Return type:

int

property Fext

Get the external aerodynamic force vector.

Returns:

External aerodynamic force vector as a numpy array.

Return type:

numpy.ndarray

property phi

Get the adjoint vector.

Returns:

Adjoint vector as a numpy array.

Return type:

numpy.ndarray

property pLoad

The contribution of the aerodynamic adjoint to the structural adjoint RHS. This term is given by dAdu^T*psi and is subtracted from the structural adjoint RHS.

Returns:

Aerodynamic adjoint contribution to the structural adjoint RHS.

Return type:

numpy.ndarray

property dSdu

Get the product of the derivative of structural residual with respect to state variables and the structural adjoint vector.

Returns:

Product of the derivative of structural residual with respect to state variables and the structural adjoint vector.

Return type:

numpy.ndarray

property adjRHS

Get the adjoint right-hand side vector.

Returns:

Adjoint RHS vector as a numpy array.

Return type:

numpy.ndarray

property dIdu

Get the derivative of the function with respect to state variables.

Returns:

Derivative vector as a numpy array with boundary conditions applied.

Return type:

numpy.ndarray

property matVecRHS

Get the matrix-vector product right-hand side vector.

Returns:

Matrix-vector RHS vector as a numpy array.

Return type:

numpy.ndarray

property matVecSolve

Get the matrix-vector solve vector.

Returns:

Matrix-vector solve vector as a numpy array.

Return type:

numpy.ndarray

getVarName()[source]

Get name for the design variables in pyOpt. Only needed if more than 1 pytacs object is used in an optimization

Returns:

varName -- Name of the design variables used in setDesignVars() dict.

Return type:

str

setDVGeo(DVGeo, pointSetKwargs=None)[source]

Set the DVGeometry object that will manipulate 'geometry' in this object. Note that TACS doe not strictly need a DVGeometry object, but if optimization with geometric changes is desired, then it is required.

Parameters:
  • dvGeo (A DVGeometry object.) -- Object responsible for manipulating the constraints that this object is responsible for.

  • pointSetKwargs (dict) -- Keyword arguments to be passed to the DVGeo addPointSet call. Useful for DVGeometryMulti, specifying FFD projection tolerances, etc

Examples

>>> StructProblem.setDVGeo(DVGeo)
getDVGeo()[source]

Get the DVGeometry object.

Returns:

DVGeometry -- DVGeometry object.

Return type:

pygeo.parameterization.BaseDVGeometry or None

setDesignVars(x)[source]

Set the variables in the x-dict for this object.

Parameters:

x (dict) -- Dictionary of variables which may or may not contain the design variable names this object needs

addConstraint(constr)[source]

Add pyTACS constraint object to StructProblem.

Parameters:

constr (tacs.constraints.TACSConstraint) -- pyTACS Constraint object

addVariablesPyOpt(optProb)[source]

Add the current set of variables to the optProb object.

Parameters:

optProb (pyOpt_optimization class) -- Optimization problem definition to which variables are added

addConstraintsPyOpt(optProb, nonLinear=True, linear=True)[source]

Add any linear constraints that were generated during setup to the specified pyOpt problem.

Parameters:
  • optProb (Optimization instance) -- Optimization problem object to add constraints to

  • nonLinear (bool) -- Flag to include non-linear constraints.

  • linear (bool) -- Flag to include linear constraints.

getNodes()[source]

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

solve(damp=1.0, useAitkenAcceleration=False, dampLB=0.2, loadScale=1.0)[source]

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

Parameters:
  • damp (float) -- Value to use to damp the solution update.

  • useAitkenAcceleration (bool) -- Flag to use aitkenAcceleration. Only applicable for aerostructural problems

  • dampLB (float) -- Lower bound for the damping parameter used in Aitken acceleration.

  • loadScale (float) -- Value to scale external loads by. Only useful for load step approach on nonlinear problems.

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

Evaluate the current functions

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

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

evalConstraints(fcon, evalCons=None, ignoreMissing=False, nonLinear=True, linear=False)[source]

Evaluate values for constraints. The constraints corresponding to the strings in evalCons are evaluated and updated into the provided dictionary.

Parameters:
  • fcon (dict) -- Dictionary into which the constraints are saved.

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

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

  • nonLinear (bool) -- Flag to include non-linear constraints.

  • linear (bool) -- Flag to include linear constraints.

evalFunctionsSens(funcsSens, evalFuncs=None)[source]

Evaluate the sensitivity of the desired functions

Parameters:
  • funcsSens (dict) -- Dictionary into which the function sensitivities are saved

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

evalConstraintsSens(fconSens, evalCons=None, nonLinear=True, linear=False)[source]

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

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

  • evalCons (iterable object containing strings) -- The constraints the user wants returned

  • nonLinear (bool) -- Flag to include non-linear constraints.

  • linear (bool) -- Flag to include linear constraints.

getVarsPerNode()[source]

The number of degrees of freedom used at each output location.

Returns:

Number of degrees of freedom of each node in the domain.

Return type:

int

getNumNodes()[source]

Get the number of nodes on this processor.

Returns:

Number of nodes on this processor.

Return type:

int

getOrigDesignVars()[source]

Get an array holding the original values for all design variables added to TACS at time of initialization.

Returns:

Array containing all design variable values across all processors.

Return type:

numpy.ndarray

getDesignVarRange()[source]

Get arrays containing the lower and upper bounds for the design variables, in the form needed by pyoptsparse.

Returns:

Lower and upper bounds for the design variables.

Return type:

tuple of numpy.ndarray

getDesignVarScales()[source]

Get an array containing the scaling factors for the design variables, in the form needed by pyoptsparse. method.

Returns:

Scaling values for design variables.

Return type:

numpy.ndarray

getNumDesignVars()[source]

Get total number of structural design variables across all processors.

Returns:

Total number of structural design variables.

Return type:

int

getStates(states=None)[source]

Return the current state values for the problem.

Parameters:

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

Returns:

Current state vector.

Return type:

numpy.ndarray

setStates(states=None)[source]

Set the current state values for the problem.

Parameters:

states (tacs.TACS.Vec or numpy.ndarray) -- Vector to replace current state variables with.

getExternalForce(Fext=None)[source]

Return the current external coupling force vector for the problem.

Parameters:

Fext (tacs.TACS.Vec or numpy.ndarray, optional) -- Vector to place current external coupling force vector into.

Returns:

Current external coupling force vector.

Return type:

numpy.ndarray

setExternalForce(Fext)[source]

Set the external coupling force vector for the problem.

Parameters:

Fext (tacs.TACS.Vec or numpy.ndarray) -- External coupling force vector to set.

getStateNorm()[source]

Get the norm of the current state vector.

Returns:

Norm of the current state vector.

Return type:

float

getResidual()[source]

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

Returns:

Computed residuals for problem.

Return type:

numpy.ndarray

getResNorms()[source]

Return the initial, starting and final residual norms. Note that the same norms are used for both solution and adjoint computations.

Returns:

Initial, starting, and final residual norms.

Return type:

tuple of float

setResNorms(initNorm=None, startNorm=None, finalNorm=None)[source]

Set one of these norms if not None. This is typically only used with aerostructural analysis.

Parameters:
  • initNorm (float, optional) -- Initial residual norm to set.

  • startNorm (float, optional) -- Starting residual norm to set.

  • finalNorm (float, optional) -- Final residual norm to set.

zeroVectors()[source]

Zero all the TACS b-vectors.

property evalFuncs

Get strings for added evaluation functions.

Returns:

List of function names that can be evaluated.

Return type:

list of str

getAdjoint(func)[source]

Return the adjoint values for func if they exist. Otherwise just return zeros.

Parameters:

func (str) -- Name of the function.

Returns:

Adjoint values for the specified function.

Return type:

numpy.ndarray

setAdjoint(adjoint, func=None)[source]

Set externally supplied adjoint values.

Parameters:
  • adjoint (float, numpy array) -- An array of size getStateSize() to be copied to the structural adjoint variables

  • func (str) -- Name of function to set adjoint for

setAdjointRHS(func)[source]

Set the adjoint right-hand side for given function.

Parameters:

func (str) -- Name of the function to set adjoint RHS for.

getdSduVec(inVec)[source]

Evaluate the direct residual.

Parameters:

inVec (numpy.ndarray) -- Input vector.

Returns:

Direct residual vector.

Return type:

numpy.ndarray

getdSduTVec(inVec)[source]

Evaluate the adjoint residual.

Parameters:

inVec (numpy.ndarray) -- Input vector.

Returns:

Adjoint residual vector.

Return type:

numpy.ndarray

getdIdXpt(func)[source]

Get the (partial) derivative of the structural function with respect to all structural nodes for the current structProblem.

Parameters:

func (str) -- Name of the function to compute derivatives for.

Returns:

Derivative vector with respect to structural nodes.

Return type:

tacs.TACS.Vec

getdIdXdv(func)[source]

Get the (partial) derivative of the structural function with respect to all structural design variables.

Parameters:

func (str) -- Name of the function to compute derivatives for.

Returns:

Derivative vector with respect to structural design variables.

Return type:

tacs.TACS.Vec

getdIduTransProd(vecT, evalFuncs=None)[source]

Perform the transpose matvec of the getdIduProd function and return a TACS vector. The result is stored locally.

We assume that evalFuncs is valid.

Parameters:
  • vecT (dict) -- Dictionary containing transpose vector values.

  • evalFuncs (list of str, optional) -- List of function names to evaluate.

Returns:

Transpose matrix-vector product result.

Return type:

numpy.ndarray

getdIdXdvTransProd(vecT, evalFuncs=None)[source]

Perform the transpose matvec of the getdIdXdvProd function and return a vector.

We assume that structProblem is set and that handles is valid.

Note that handles contains the list of function handles, whereas the names are contained in evalFuncs.

Parameters:
  • vecT (dict) -- Dictionary containing transpose vector values.

  • evalFuncs (list of str, optional) -- List of function names to evaluate.

Returns:

Dictionary containing transpose matrix-vector product results.

Return type:

dict

globalNKPreCon(inVec)[source]

This function is ONLY used as a preconditioner to the global Aero-Structural system.

Parameters:

inVec (numpy.ndarray) -- Input vector for preconditioning.

Returns:

Preconditioned output vector.

Return type:

numpy.ndarray

globalAdjointPreCon(inVec)[source]

This function is ONLY used as a preconditioner to the global Aero-Structural adjoint system. For the linear symmetric case this is actually identical to the NKPreCon function above.

Parameters:

inVec (numpy.ndarray) -- Input vector for preconditioning.

Returns:

Preconditioned output vector.

Return type:

numpy.ndarray

globalDirectPreCon(inVec)[source]

This function is ONLY used as a preconditioner to the global Aero-Structural direct system. For the linear symmetric case this is actually identical to the NKPreCon function above.

Parameters:

inVec (numpy.ndarray) -- Input vector for preconditioning.

Returns:

Preconditioned output vector.

Return type:

numpy.ndarray

assembleAdjointRHS()[source]

Compute the final adjoint RHS: Adjoint RHS should be: dIdu - dAdu^T*psi - dSdu^T*phi

solveAdjoint(damp=1.0)[source]

Solve the structural adjoint.

Parameters:

damp (float) -- A damping variable for adjoint update. Typically only used in multidisciplinary analysis

getdRdXptPhi(funcs)[source]

Get the result of \([dR/dX_{nodes}]^T \phi\) for each of the functions in the function list.

This is the total sensitivity calculation.

Parameters:

funcs (list of str) -- List of function names.

Returns:

List of sensitivity products for each function.

Return type:

list of tacs.TACS.Vec

getdRdXdvPhi(funcs)[source]

Get the result of \([dR/dX_{dv}]^T \phi\) for the total sensitivity calculation.

Parameters:

funcs (list of str) -- List of function names.

Returns:

List of sensitivity products for each function.

Return type:

list of tacs.TACS.Vec

getdRdXdvTransProd()[source]

Perform the transpose matvec of getdRdXdvProd. The result is returned as a dict of numpy arrays.

Returns:

Dictionary containing transpose matrix-vector product results.

Return type:

dict

convertDesignVecToDict(dvVec)[source]

Convert a design vector to a dictionary format.

Parameters:

dvVec (tacs.TACS.Vec or numpy.ndarray) -- Design vector to convert.

Returns:

Dictionary containing the design vector with variable name as key.

Return type:

dict

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

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

This function writes the external coupling loads stored in Fext to a file. This is typically used to save aero loads from an aerostructural run.

Parameters:
  • outputDir (str, optional) -- Output directory for the force file.

  • baseName (str, optional) -- Base name for the force file. If None, uses problem name + "_external_forces".

  • number (int, optional) -- Number to append to the filename.

readExternalForceFile(fileName)[source]

Reads in a force file and sets the external forces in the structural problem. This is typically used to read in loads saved from an aerostructural run.

Parameters:

fileName (str) -- Filename for force file. Should have .dat or .bdf extension.