constitutive module

The tacs.constitutive contains classes responsible for defining constitutive behaviors for elements.

Material classes

Most constitutive classes require a material properties class to setup. The following classes are available in TACS:

class tacs.constitutive.MaterialProperties

Bases: object

This class stores the mechanical and thermal material properties for isotropic and anisotropic materials.

The goal of this class is to store a set of material properties that can be queried by constitutive classes for beams, shells, plane stress and solid elements. The minimum set of properties consists of isotropic mechanical properties, with zero thermal expansion.

There are several different ways of specifying material properties. The following describes several of the possible ways and the appropriate situations:

'rho' + 'specific_heat' + 'kappa':

Specifies the density, specific heat, and thermal conductivity for an isotropic material. This is appropriate for 2D or 3D heat conduction analysis.

'rho' + 'specific_heat' + 'kappa1' + 'kappa2':

Specifies the density, specific heat, and thermal conductivity for an orthotropic material. This is appropriate for 2D heat conduction analysis.

'rho' + 'specific_heat' + 'kappa1' + 'kappa2' + 'kappa3':

Specifies the density, specific heat, and thermal conductivity for an orthotropic material. This is appropriate for 3D heat conduction analysis.

'rho' + 'E' + 'nu' + 'ys':

Specifies the density, Youngs' modulus, Poisson's ratio, and yield strength for an isotropic material. This is appropriate for 2D or 3D elastic structural analysis.

'rho' + 'E1' + 'E2' + 'nu12' + 'G12' + 'T1' + 'T2' + 'C1' + 'C2' + 'S12':

Specifies the density, Youngs' moduli, Poisson's ratios, and strength values for an orthotropic material. This is appropriate for 2D elastic structural analysis.

'rho' + 'E1' + 'E2' + 'E3' + 'nu12' + 'nu13' + 'nu23' + 'G12' + 'G13' + 'G23' + 'T1' + 'T2' + 'T3' + 'C1' + 'C2' + 'C3' + 'S12' + 'S13' + 'S23':

Specifies the density, Youngs' moduli, Poisson's ratios, and strength values for an orthotropic material. This is appropriate for 3D elastic structural analysis.

'rho' + 'specific_heat' + 'kappa' + 'alpha' + 'E' + 'nu' + 'ys':

Specifies the density, specific heat, thermal conductivity, thermal expansion, Youngs' modulus, Poisson's ratio, and yield strength for an isotropic material. This is appropriate for 2D or 3D thermoelastic structural analysis.

'rho' + 'specific_heat' + 'kappa1' + 'kappa2' + 'alpha1' + 'alpha2' + 'E1' + 'E2' + 'nu12' + 'G12' + 'T1' + 'T2' + 'C1' + 'C2' + 'S12':

Specifies the density, specific heat, thermal conductivity, thermal expansion, Youngs' moduli, Poisson's ratios, and strength values for an orthotropic material. This is appropriate for 2D thermoelastic structural analysis.

'rho' + 'specific_heat' + 'kappa1' + 'kappa2' + 'kappa3' + 'alpha1' + 'alpha2' + 'alpha3' + 'E1' + 'E2' + 'E3' + 'nu12' + 'nu13' + 'nu23' + 'G12' + 'G13' + 'G23' + 'T1' + 'T2' + 'T3' + 'C1' + 'C2' + 'C3' + 'S12' + 'S13' + 'S23':

Specifies the density, specific heat, thermal conductivity, thermal expansion, Youngs' moduli, Poisson's ratios, and strength values for an orthotropic material. This is appropriate for 3D thermoelastic structural analysis.

All parameters are optional and specified using a keyword arg format.

Parameters:
  • rho (float or complex, optional) -- The material density (keyword argument).

  • specific_heat (float or complex, optional) -- The material specific heat (keyword argument).

  • kappa (float or complex, optional) -- The material isotropic thermal conductivity (keyword argument).

  • kappa1 (float or complex, optional) -- The material orthotropic thermal conductivity in the '1' direction (keyword argument).

  • kappa2 (float or complex, optional) -- The material orthotropic thermal conductivity in the '2' direction (keyword argument).

  • kappa3 (float or complex, optional) -- The material orthotropic thermal conductivity in the '3' direction (keyword argument).

  • alpha (float or complex, optional) -- The material isotropic thermal expansion coeeficient (keyword argument).

  • alpha1 (float or complex) -- The material orthotropic thermal expansion coefficient in the '1' direction (keyword argument).

  • alpha2 (float or complex, optional) -- The material orthotropic thermal expansion coefficient in the '2' direction (keyword argument).

  • alpha3 (float or complex, optional) -- The material orthotropic thermal expansion coefficient in the '3' direction (keyword argument).

  • E (float or complex, optional) -- The material isotropic Youngs' modulus (keyword argument).

  • E1 (float or complex, optional) -- The material orthotropic Youngs' modulus in the '1' direction (keyword argument).

  • E2 (float or complex, optional) -- The material orthotropic Youngs' modulus in the '2' direction (keyword argument).

  • E3 (float or complex, optional) -- The material orthotropic Youngs' modulus in the '3' direction (keyword argument).

  • nu (float or complex, optional) -- The material isotropic Poisson's ratio (keyword argument).

  • nu12 (float or complex, optional) -- The material orthotropic Poisson's ratio in the '1-2' plane (keyword argument).

  • nu13 (float or complex, optional) -- The material orthotropic Poisson's ratio in the '1-3' plane (keyword argument).

  • nu23 (float or complex, optional) -- The material orthotropic Poisson's ratio in the '2-3' plane (keyword argument).

  • G (float or complex, optional) -- The material isotropic shear modulus (keyword argument).

  • G12 (float or complex, optional) -- The material orthotropic shear modulus in the '1-2' plane (keyword argument).

  • G13 (float or complex, optional) -- The material orthotropic shear modulus in the '1-3' plane (keyword argument).

  • G23 (float or complex, optional) -- The material orthotropic shear modulus in the '2-3' plane (keyword argument).

  • ys (float or complex, optional) -- The material isotropic yield strength (keyword argument).

  • T1 (float or complex, optional) -- The material orthotropic tension strength in the '1' direction (keyword argument).

  • T2 (float or complex, optional) -- The material orthotropic tension strength in the '2' direction (keyword argument).

  • T3 (float or complex, optional) -- The material orthotropic tension strength in the '3' direction (keyword argument).

  • C1 (float or complex, optional) -- The material orthotropic compression strength in the '1' direction (keyword argument).

  • C2 (float or complex, optional) -- The material orthotropic compression strength in the '2' direction (keyword argument).

  • C3 (float or complex, optional) -- The material orthotropic compression strength in the '3' direction (keyword argument).

  • S12 (float or complex, optional) -- The material orthotropic shear strength in the '1-2' plane (keyword argument).

  • S13 (float or complex, optional) -- The material orthotropic shear strength in the '1-3' plane (keyword argument).

  • S23 (float or complex, optional) -- The material orthotropic shear strength in the '2-3' plane (keyword argument).

generateBDFCard()

Generate pyNASTRAN card class based on current design variable values.

Returns:

pyNastran card holding material information

Return type:

card (pyNastran.bdf.cards.materials.MAT1 or pyNastran.bdf.cards.materials.MAT8)

getMaterialProperties()

Return a dictionary of the material properties

Returns:

Dictionary holding material property information

Return type:

mat (dict)

getNastranID()

Get material ID assigned in NASTRAN card for this object.

Returns:

ID number associated with this object's NASTRAN card

Return type:

id (int)

setDensity(rho)

Set the density property values

Parameters:

rho (float or complex) -- The material density

setNastranID(id)

Set material ID to be used in NASTRAN card for this object. Should be set before generateBDFCard is called.

Parameters:

id (int) -- ID number to associate with this object's NASTRAN card

setSpecificHeat(specific_heat)

Set the density property values

Parameters:

rho (float or complex) -- The material specific heat

class tacs.constitutive.OrthotropicPly

Bases: object

The following class holds the material stiffness and strength properties for an orthotropic ply. This class is used by several constitutive classes within TACS.

The interaction coefficient for the Tsai-Wu failure criterion is set to zero by default. If a value of C, the failure stress under combined in-plane loading, is supplied, the interaction coefficient is determined. Be careful - the value can easily fall outside acceptable bounds - these are tested during initialization.

Parameters:
  • ply_thickness (float or complex) -- The ply thickness.

  • props (MaterialProperties) -- The ply material property.

  • max_strain_criterion (bool) -- Flag to determine if max strain strength criterion is to be used. Defaults to False (i.e. use max strength).

getMaterialProperties()

Get the MaterialProperties class associated with this object

Returns:

TACS material property class associated with object.

Return type:

prop (tacs.constitutive.MaterialProperties)

Constitutive classes

Most Element classes require a Constitutive class in their setup procedure. These objects typically defines mass, stiffness, failure and buckling behaviors for the elements that reference them. The following Constitutive classes are available in TACS:

class tacs.constitutive.BasicBeamConstitutive

Bases: BeamConstitutive

Timoshenko theory based constitutive object for an general beam.

Note

TACS uses a negative sign convention in the product of inertia definition, for example:

Iyz = -int[y * z * dA]

The moments of area are always positive, as usual:

Iyy = int[(z^2 * dA]

Parameters:
  • props (MaterialProperties) -- The material property.

  • A (float or complex) -- Beam cross-sectional area (keyword argument). Defaults to 0.0.

  • J (float or complex) -- Beam polar moment of area about x-axis (keyword argument). Defaults to 0.0.

  • Iy (float or complex) -- Beam area moment of area about y-axis (keyword argument). Defaults to 0.0.

  • Iz (float or complex) -- Beam area moment of area about z-axis (keyword argument). Defaults to 0.0.

  • Iyz (float or complex) -- Beam product of area in yz-plane (keyword argument). Defaults to 0.0.

  • ky (float or complex) -- Shear correction factor in y-direction (keyword argument). Defaults to 5/6.

  • kz (float or complex) -- Shear correction factor in z-direction (keyword argument). Defaults to 5/6.

generateBDFCard()

Generate pyNASTRAN card class based on current design variable values.

Returns:

pyNastran card holding property information

Return type:

card (pyNastran.bdf.cards.properties.bars.PBAR)

class tacs.constitutive.BladeStiffenedShellConstitutive

Bases: StiffenedShellConstitutive

This constitutive class models a shell stiffened with T-shaped stiffeners. The stiffeners are not explicitly modelled. Instead, their stiffness is "smeared" across the shell.

Parameters:
  • panelPly (tacs.constitutive.OrthotropicPly) -- Ply model to use for the panel

  • stiffenerPly (tacs.constitutive.OrthotropicPly) -- Ply model to use for the stiffener

  • panelLength (float or complex) -- Panel length DV value

  • stiffenerPitch (float or complex) -- Stiffener pitch DV value

  • panelThick (float or complex) -- Panel thickness DV value

  • panelPlyAngles (numpy.ndarray[float or complex]) -- Array of ply angles in the panel

  • panelPlyFracs (numpy.ndarray[float or complex]) -- Array of ply fractions in the panel

  • stiffenerHeight (float or complex) -- Stiffener height DV value

  • stiffenerThick (float or complex) -- Stiffener thickness DV value

  • stiffenerPlyAngles (numpy.ndarray[float or complex]) -- Array of ply angles for the stiffener

  • stiffenerPlyFracs (numpy.ndarray[float or complex]) -- Array of ply fractions for the stiffener

  • kcorr (float or complex, optional) -- Shear correction factor, defaults to 5.0/6.0

  • flangeFraction (float, optional) -- Ratio of the stiffener base width to the stiffener height. Defaults to 1.0

  • panelLengthNum (int, optional) -- Panel lenth DV number, passing a negative value tells TACS not to treat this as a DV. Defaults to -1

  • stiffenerPitchNum (int, optional) -- Stiffener pitch DV number, passing a negative value tells TACS not to treat this as a DV. Defaults to -1

  • panelThickNum (int, optional) -- Panel thickness DV number, passing a negative value tells TACS not to treat this as a DV. Defaults to -1

  • panelPlyFracNums (numpy.ndarray[np.intc], optional) -- Array of ply fraction DV numbers in the panel, passing negative values tells TACS not to treat that ply fraction as a DV. Defaults to -1's

  • stiffenerHeightNum (int, optional) -- Stiffener height DV number, passing a negative value tells TACS not to treat this as a DV. Defaults to -1

  • stiffenerThickNum (int, optional) -- Stiffener thickness DV number, passing a negative value tells TACS not to treat this as a DV. Defaults to -1

  • stiffenerPlyFracNums (numpy.ndarray[numpy.intc], optional) -- Array of ply fraction DV numbers for the stiffener, passing negative values tells TACS not to treat that ply fraction as a DV. Defaults to -1's

Raises:
  • ValueError -- Raises error if panelPlyAngles, panelPlyFracs, or panelPlyFracNums do not have same number of entries

  • ValueError -- Raises error if stiffenerPlyAngles, stiffenerPlyFracs, or stiffenerPlyFracNums do not have same number of entries

class tacs.constitutive.BucklingGP

Bases: GaussianProcess

Gaussian Process ML model to predict N_11,cr^* non-dimensional buckling loads of global axial modes. Local axial mode predictions can also be made with gamma = 0 and xi, rho0 computed for the local panel. The ML model uses non-dimensional inputs and outputs as follows:

log(N_11,cr^*) = my_axial_GP.predict_mean_test_data(log(1+xi), log(rho0), log(1+gamma), log(1+10^3 * zeta)) log(N_12,cr^*) = my_shear_GP.predict_mean_test_data(log(1+xi), log(rho0), log(1+gamma), log(1+10^3 * zeta))

The axialGP model accomodoates making multiple test point predictions in parallel as do many ML models.

Parameters:
  • Xtrain (np.ndarray) -- a rank-1 tensor (4*N_train,) of the training dataset non-dim inputs namely [log(1+xi), log(rho0), log(1+gamma), log(1+10^3 * zeta)] for each training point or [log_xi1, log_rho01, log_gamma1, log_zeta1, log_xi2, ...]

  • alpha (np.ndarray) -- a rank-1 tensor (N_train,) containing the training weights for the ML model (computed by python training in ml_buckling repo)

  • theta (np.ndarray) -- a rank-1 tensor (13,) containing each of the model hyperparameters in the kernel function

classmethod from_csv(csv_file, theta_csv)

Construct an BucklingGP from a csv_files in the ml_buckling repo (or your own dataset csv files) which contain the training weights of the ML model and the optimal hyperparameters theta.

This is the method commonly used to construct the BucklingGP. Namely using the ml_buckling

The common construction of this class from the ml_buckling repo, https://github.com/smdogroup/ml_buckling, and is the following:

axial_gp = BucklingGP.from_csv(

csv_file=mlb.axialGP_csv, theta_csv=mlb.axial_theta_csv

) shear_gp = BucklingGP.from_csv(

csv_file=mlb.shearGP_csv, theta_csv=mlb.shear_theta_csv

)

Parameters:
  • csv_file (filepath or str) -- path to a csv_file that holds the Xtrain training data with 5 columns for each entry of Xtrain namely [log(1+xi), log(rho0), log(1+gamma), log(1+10^3 * zeta), alphaTrain].

  • theta_csv (filepath or str) -- path to a csv file that holds the theta hyperparameters in one column with 13 entries.

  • Returns -- BucklingGP object

class tacs.constitutive.CompositeShellConstitutive

Bases: ShellConstitutive

This constitutive class defines the stiffness properties for a composite laminate first-order shear deformation theory (FSDT) shell type element.

Parameters:
  • ply_list (list[OrthotropicPly]) -- List of ply properties in layup.

  • ply_thicknesses (numpy.ndarray[float or complex]) -- Array of ply thicknesses in layup.

  • ply_angles (numpy.ndarray[float or complex]) -- Array of ply angles (in radians) in layup.

  • kcorr (float or complex, optional) -- FSDT shear correction factor. Defaults to 5.0/6.0.

  • tOffset (float or complex, optional) -- Offset distance of reference plane (where nodes are located) relative to thickness mid-plane. Measured in fraction of shell thickness. A value of 0.5 places the reference plane at the top of the plate, a value of 0.0 at the plate mid-plane, and a value of -0.5 at the bottom of the plate. Defaults to 0.0.

generateBDFCard()

Generate pyNASTRAN card class based on current design variable values.

Returns:

pyNastran card holding property information

Return type:

card (pyNastran.bdf.cards.properties.shell.PCOMP)

class tacs.constitutive.DOFSpringConstitutive

Bases: GeneralSpringConstitutive

This is the base class for the traditional spring constitutive objects with no dof coupling. Assumes 6 dofs.

Parameters:

k (array-like[float or complex]) -- Stiffness values for all 6 dofs.

generateBDFCard()

Generate pyNASTRAN card class based on current design variable values.

Returns:

pyNastran card holding property information

Return type:

card (pyNastran.bdf.cards.properties.bush.PBUSH)

class tacs.constitutive.GPBladeStiffenedShellConstitutive

Bases: StiffenedShellConstitutive

" This constitutive class models a shell stiffened with T-shaped stiffeners. The stiffeners are not explicitly modelled. Instead, their stiffness is "smeared" across the shell.

This class is a subclass of the BladeStiffenedShellConstitutive base class and includes higher-fidelity buckling constraints using Gaussian Processes for Machine Learning. The trained ML models are located on the repo, https://github.com/smdogroup/ml_buckling.

For the methods below, consider a panel with dimensions a the panel length, b the panel width,

h the panel thickness, and stiffeners along the length or 1-direction.

The Dij for axial and shear modes of the panel use the panel laminate design, with the centroid shifted towards the stiffener

only for the D11 stiffness for the global modes (the local modes use D11 at the center of the skin).

The stiffener has a lateral spacing s_p the stiffener pitch, stiffener height h_s, stiffener thickness t_s. For more information on this constitutive class, see our paper https://arc.aiaa.org/doi/abs/10.2514/6.2024-3981,

"Machine Learning to Improve Buckling Predictions for Efficient Structural Optimization of Aircraft Wings" by Sean Engelstad, Brian Burke, Graeme Kennedy.

Parameters:
  • panelPly (tacs.constitutive.OrthotropicPly) -- Ply model to use for the panel

  • stiffenerPly (tacs.constitutive.OrthotropicPly) -- Ply model to use for the stiffener

  • panelLength (float or complex) -- Panel length DV value

  • stiffenerPitch (float or complex) -- Stiffener pitch DV value

  • panelThick (float or complex) -- Panel thickness DV value

  • panelPlyAngles (numpy.ndarray[float or complex]) -- Array of ply angles in the panel

  • panelPlyFracs (numpy.ndarray[float or complex]) -- Array of ply fractions in the panel

  • stiffenerHeight (float or complex) -- Stiffener height DV value

  • stiffenerThick (float or complex) -- Stiffener thickness DV value

  • stiffenerPlyAngles (numpy.ndarray[float or complex]) -- Array of ply angles for the stiffener

  • stiffenerPlyFracs (numpy.ndarray[float or complex]) -- Array of ply fractions for the stiffener

  • panelWidth (float or complex) -- Panel width DV value

  • kcorr (float or complex, optional) -- Shear correction factor, defaults to 5.0/6.0

  • flangeFraction (float, optional) -- Ratio of the stiffener base width to the stiffener height. Defaults to 1.0

  • panelLengthNum (int, optional) -- Panel lenth DV number, passing a negative value tells TACS not to treat this as a DV. Defaults to -1

  • stiffenerPitchNum (int, optional) -- Stiffener pitch DV number, passing a negative value tells TACS not to treat this as a DV. Defaults to -1

  • panelThickNum (int, optional) -- Panel thickness DV number, passing a negative value tells TACS not to treat this as a DV. Defaults to -1

  • panelPlyFracNums (numpy.ndarray[np.intc], optional) -- Array of ply fraction DV numbers in the panel, passing negative values tells TACS not to treat that ply fraction as a DV. Defaults to -1's

  • stiffenerHeightNum (int, optional) -- Stiffener height DV number, passing a negative value tells TACS not to treat this as a DV. Defaults to -1

  • stiffenerThickNum (int, optional) -- Stiffener thickness DV number, passing a negative value tells TACS not to treat this as a DV. Defaults to -1

  • stiffenerPlyFracNums (numpy.ndarray[numpy.intc], optional) -- Array of ply fraction DV numbers for the stiffener, passing negative values tells TACS not to treat that ply fraction as a DV. Defaults to -1's

  • panelWidthNum (int, optional) -- Panel width DV number, passing a negative value tells TACS not to treat this as a DV. Defaults to -1

  • CPTstiffenerCrippling (bool) -- whether to use CPT (classical plate theory) for stiffener crippling buckling predictions in the closed-form solution case (default False)

  • panelGPs (TACSPanelGPs) -- container for the three GP models for each panel / TACS component

Raises:
  • ValueError -- Raises error if panelPlyAngles, panelPlyFracs, or panelPlyFracNums do not have same number of entries

  • ValueError -- Raises error if stiffenerPlyAngles, stiffenerPlyFracs, or stiffenerPlyFracNums do not have same number of entries

nondimCriticalGlobalAxialLoad(rho_0, xi, gamma, zeta=0.0)

predict the non-dimensional buckling load N11,cr^* for the global axial mode of a panel. Uses the closed-form solution if PanelGPs.axialGP is None or the axialGP ML model if PanelGPs.axialGP is not None Here D11 is for the centroid of the panel and stiffener, and the other Dij are for the panel at the panel center.

Parameters:
  • rho_0 (float) -- the affine aspect ratio of the panel, rho_0 = a/b * 4thRoot(D22 / D11)

  • xi (float) -- the laminate isotropy, xi = (D12 + 2 * D66) / sqrt(D11 * D22)

  • gamma (float) -- the stiffener stiffness ratio, gamma = E_1s * A_s / (s_p * D11)

  • zeta (float) -- the transverse shear parameters, zeta = A11 / A66 * (h/b)^2

  • Returns --

    N11,cr^* (float)the non-dimensional output buckling load where the dimensional buckling load N11,cr is given by

    N11,cr = (N11,cr^*) * pi^2 * sqrt(D11 * D22) / b^2 / (1 + delta) with delta = E1s * As / (E1p * sp * h)

nondimCriticalGlobalShearLoad(rho_0, xi, gamma, zeta=0.0)

predict the non-dimensional buckling load N12,cr^* for the global shear mode of a panel. Uses the closed-form solution if PanelGPs.shearGP is None or the shearGP ML model if PanelGPs.shearGP is not None. Here D11 is for the centroid of the panel and stiffener, and the other Dij are for the panel at the panel center.

Parameters:
  • rho_0 (float) -- the affine aspect ratio of the panel, rho_0 = a/b * 4thRoot(D22 / D11)

  • xi (float) -- the laminate isotropy, xi = (D12 + 2 * D66) / sqrt(D11 * D22)

  • gamma (float) -- the stiffener stiffness ratio, gamma = E_1s * A_s / (s_p * D11)

  • zeta (float) -- the transverse shear parameters, zeta = A11 / A66 * (h/b)^2

  • Returns --

    N12,cr^* (float)the non-dimensional output buckling load where the dimensional buckling load N12,cr is given by

    N12,cr = (N12,cr^*) * pi^2 * 4thRoot(D11 * D22^3) / b^2

nondimCriticalLocalAxialLoad(rho_0, xi, zeta=0.0)

predict the non-dimensional buckling load N11,cr^* for the local axial mode of a panel. Uses the closed-form solution if PanelGPs.axialGP is None or the axialGP ML model if PanelGPs.axialGP is not None Here D11 is for the center of the skin, and the other Dij are for the panel at the panel center.

Parameters:
  • rho_0 (float) -- the affine aspect ratio of the panel, rho_0 = a/s_p * 4thRoot(D22 / D11)

  • xi (float) -- the laminate isotropy, xi = (D12 + 2 * D66) / sqrt(D11 * D22)

  • zeta (float) -- the transverse shear parameters, zeta = A11 / A66 * (h/b)^2

  • Returns --

    N11,cr^* (float)the non-dimensional output buckling load where the dimensional buckling load N11,cr is given by

    N11,cr = (N11,cr^*) * pi^2 * sqrt(D11 * D22) / s_p^2

nondimCriticalLocalShearLoad(rho_0, xi, zeta=0.0)

predict the non-dimensional buckling load N12,cr^* for the local shear mode of a panel. Uses the closed-form solution if PanelGPs.shearGP is None or the shearGP ML model if PanelGPs.shearGP is not None. D11 here is for the center of the panel, and the other Dij are for the panel at the panel center.

Parameters:
  • rho_0 (float) -- the affine aspect ratio of the panel, rho_0 = a/s_p * 4thRoot(D22 / D11)

  • xi (float) -- the laminate isotropy, xi = (D12 + 2 * D66) / sqrt(D11 * D22)

  • zeta (float) -- the transverse shear parameters, zeta = A11 / A66 * (h/b)^2

  • Returns --

    N12,cr^* (float)the non-dimensional output buckling load where the dimensional buckling load N12,cr is given by

    N12,cr = (N12,cr^*) * pi^2 * 4thRoot(D11 * D22^3) / s_p^2

nondimStiffenerCripplingLoad(rho_0, xi, genPoiss, zeta=0.0)

predict the non-dimensional buckling load N11,cr^* for the stiffener crippling mode of the stiffener. Uses the closed-form solution if PanelGPs.cripplingGP is None or the cripplingGP ML model if PanelGPs.cripplingGP is not None

Parameters:
  • rho_0 (float) -- the affine aspect ratio of the stiffener, rho_0 = a/b * 4thRoot(D22 / D11)

  • xi (float) -- the laminate isotropy, xi = (D12 + 2 * D66) / sqrt(D11 * D22)

  • genPoiss (float) -- the generalized Poisson's ratio for weak BCs such as the stiffener with its free edge, eps = D12 / (D12 + 2 * D66)

  • zeta (float) -- the transverse shear parameters, zeta = A11 / A66 * (h/b)^2

  • Returns --

    N11,cr^* (float)the non-dimensional output buckling load where the dimensional buckling load N11,cr is given by

    N11,cr = (N11,cr^*) * pi^2 * sqrt(D11s * D22s) / h_s^2; here D11s, D22s are for the stiffener laminate

setCPTstiffenerCrippling(CPTcripplingMode)

choose whether to use stiffener crippling by: CPT solution from Sean's paper if CPTcripplingMode = True DOD experimental solution from Ali's superclass if CPTcripplingMode = False

Parameters:

CPTcripplingMode (bool) -- whether to use CPT crippling solution or not

setWriteDVMode(newMode)

Set mode for writing DV inputs to the .f5 files 0 - write DVs, 1 - write nondim params, 2 - write failure indexes this is a useful tool to investigate and debug final designs.

Parameters:

newMode (int) -- new mode input for the writeDVMode

test_all_derivative_tests(epsilon, printLevel)

test all the internal derivative tests

class tacs.constitutive.GaussianProcess

Bases: object

Base class for constructing Gaussian Process ML models to predict buckling loads. This is an abstract base class and hence has no constructor, only methods used by its subclasses.

getNparam()

Get the number of input parameters in the Gaussian Process model (4 for all of the current base classes)

Returns:

the number of input parameters in Xtest [4 for the current implemented classes]

Return type:

nparams (int)

getNtrain()

Get the number of input training data points used in the model training weights

Returns:

the number of training data points in alpha_train, the training weights

Return type:

ntrain (int)

getTrainingData()

Get the training dataset Xtrain used for training the Gaussian Process ML model

Returns:

a rank-1 tensor (4*N_train,) of the training dataset non-dim inputs namely [log(1+xi), log(rho0), log(1+gamma), log(1+10^3 * zeta)] for each training point or [log_xi1, log_rho01, log_gamma1, log_zeta1, log_xi2, ...]

Return type:

Xtrain (np.ndarray)

kernel(Xtest, Xtrain)

call the kernel function of the Gaussian Process model on an arbitrary test point and training point pair. this function is for debugging purposes only => to make sure the kernel functions implemented inside the TACS GP models match those of the ml_buckling repo.

Parameters:
  • Xtest (np.ndarray) -- a rank-1 tensor of size (nparams,) currently size (4,) of the log non-dimensional inputs at the test datapoint. the inputs are [log(1+xi), log(rho0), log(1+gamma), log(1+10^3 * zeta)].

  • Xtrain (np.ndarray) -- a rank-1 tensor of size (nparams,) currently size (4,) of the log non-dimensional inputs at the training datapoint. the inputs are [log(1+xi), log(rho0), log(1+gamma), log(1+10^3 * zeta)].

predict_mean_test_data(Xtest)

This method makes buckling predictions at each point of an Xtest dataset and returns the Ytest tensor of buckling load predictions. Note that the Xtest inputs and Ytest outputs are in log-space (more details below).

Parameters:
  • np.ndarray (Xtest -) -- a rank 1 tensor of size (4*N_test) which contains the list of nondim buckling inputs namely [log(1+xi), log(rho0), log(1+gamma), log(1+10^3 * zeta)] concatenated for each point in the test set. specifically: [log_xi1, log_rho01, log_gamma1, log_zeta1, log_xi2, ...]

  • Returns -- Ytest (np.ndarray) : a rank 1 tensor Ytest of size (N_test,) containing the log(N_ij,cr^*) aka log buckling load outputs where N_ij,cr^* is N_11,cr^* for the axial GP, N_12,cr^* for the shear GP.

recompute_alpha(Ytrain)

a helper function that recomputes the training weights alpha_train of size (n_train,) [a rank 1-tensor]. this function solves the linear equation [K(X_train,X_train;theta) + sigma_n^2 * I] * alpha_train = Y_train for the training weights alpha_train with sigma_n from the hyperparameters theta which come from inside the model.

this is called whenever the ksWeight or theta hyperparameters are changed so we update the ML model.

Parameters:

Ytrain (np.ndarray) -- training dataset log(buckling load) outputs in order to retrain the ML model, it's a rank-1 tensor of size (Ntrain,)

setKS(ksWeight, Ytrain)

Set the KS parameter of the Gaussian Process ML model used in the kernel functions. Ytrain is provided since this is needed to retrain the ML model weights for the new KS input,

and Ytrain is not stored in the GP data structure.

Parameters:
  • ksWeight (float) -- the Kresselmeier-Steinhauser aggregation parameter rho_KS used for smooth relu and smooth abs value in the GP kernel functions.

  • Ytrain (np.ndarray) -- training dataset log(buckling load) outputs in order to retrain the ML model, it's a rank-1 tensor of size (Ntrain,)

setTheta(theta, Ytrain)

Set the hyperparameters theta of the Gaussian Process ML model used in the kernel functions. Ytrain is provided since this is needed to retrain the ML model weights for the new KS input,

and Ytrain is not stored in the GP data structure.

Parameters:
  • theta (np.ndarray) -- rank 1-tensor of model hyperparameters (currently a length 13 1-tensor).

  • Ytrain (np.ndarray) -- training dataset log(buckling load) outputs in order to retrain the ML model, it's a rank-1 tensor of size (Ntrain,)

test_all_gp_tests(epsilon, printLevel)

run all the internal Gaussian Process ML model derivative tests

Parameters:
  • epsilon (TacsScalar) -- the step size for each derivative test

  • printLevel (int) -- an input of 0 does not print each derivative test result, an input of 1 does print each test result to terminal.

  • Returns -- max relative error (TacsScalar) : the max relative error among all internal derivative tests of the model

class tacs.constitutive.GeneralMassConstitutive

Bases: Constitutive

This is the base class for the fully general point mass constitutive objects. Assumes 6 dofs (3 translations + 3 rotations). The mass properties of this object are specified using a symmetric 6 x 6 mass matrix, as shown below.

M[ 0] M[ 1] M[ 2] M[ 3] M[ 4] M[ 5]
M[ 1] M[ 6] M[ 7] M[ 8] M[ 9] M[10]
M[ 2] M[ 7] M[11] M[12] M[13] M[14]
M[ 3] M[ 8] M[12] M[15] M[16] M[17]
M[ 4] M[ 9] M[13] M[16] M[18] M[19]
M[ 5] M[10] M[14] M[17] M[19] M[20]
Parameters:

M (array-like[float or complex]) -- Flattened array containing one side of symmetric mass matrix entries.

evalMassMatrix()

Evaluate 6x6 symmetric mass matrix associated with this element.

Returns:

Length 21 flattened array representing unique entries of mass matrix

Return type:

M (numpy.ndarray)

class tacs.constitutive.GeneralSpringConstitutive

Bases: Constitutive

This is the base class for the fully general spring constitutive objects. Assumes 6 dofs (3 translations + 3 rotations). The spring properties of this object are specified using a symmetric 6 x 6 stiffness matrix, as shown below.

K[ 0] K[ 1] K[ 2] K[ 3] K[ 4] K[ 5]
K[ 1] K[ 6] K[ 7] K[ 8] K[ 9] K[10]
K[ 2] K[ 7] K[11] K[12] K[13] K[14]
K[ 3] K[ 8] K[12] K[15] K[16] K[17]
K[ 4] K[ 9] K[13] K[16] K[18] K[19]
K[ 5] K[10] K[14] K[17] K[19] K[20]
Parameters:

K (array-like[float or complex]) -- Flattened array containing one side of symmetric stiffness matrix entries.

class tacs.constitutive.IsoRectangleBeamConstitutive

Bases: BeamConstitutive

Timoshenko theory based constitutive object for a solid rectangular beam.

The thickness dimension is assumed to measured along the beam's local reference axis, the width is perpindicular to the reference axis.

Parameters:
  • props (MaterialProperties) -- The material property.

  • width (float or complex, optional) -- Cross-section width (keyword argument). Defaults to 1.0.

  • wNum (int, optional) -- Design variable number to assign to width (keyword argument). Defaults to -1 (i.e. no design variable).

  • wlb (float or complex, optional) -- Lower bound on width (keyword argument). Defaults to 0.0.

  • wub (float or complex, optional) -- Upper bound on width diameter (keyword argument). Defaults to 10.0.

  • t (float or complex, optional) -- Cross-section thickness (keyword argument). Defaults to 0.1.

  • tNum (int, optional) -- Design variable number to assign to thickness (keyword argument). Defaults to -1 (i.e. no design variable).

  • tlb (float or complex, optional) -- Lower bound on thickness (keyword argument). Defaults to 0.0.

  • tub (float or complex, optional) -- Upper bound on thickness (keyword argument). Defaults to 10.0.

  • wOffset (float or complex, optional) -- Offset distance along width axis of reference axis (where nodes are located) relative to elastic axis. Measured in fraction of cross-section width. A value of 0.5 places the reference axis in the plane z=/2, a value of 0.0 at z=0.0, and a value of -0.5 at z=-w/2 (where z is the local beam axis). Defaults to 0.0.

  • tOffset (float or complex, optional) -- Offset distance along thickness axis of reference axis (where nodes are located) relative to elastic axis. Measured in fraction of cross-section thickness. A value of 0.5 places the reference axis in the plane y=t/2, a value of 0.0 at y=0.0, and a value of -0.5 at y=-t/2 (where y is the local beam axis). Defaults to 0.0.

generateBDFCard()

Generate pyNASTRAN card class based on current design variable values.

Returns:

pyNastran card holding property information

Return type:

card (pyNastran.bdf.cards.properties.bars.PBARL)

class tacs.constitutive.IsoShellConstitutive

Bases: ShellConstitutive

This constitutive class defines the stiffness properties for a isotropic first-order shear deformation theory shell type element.

Parameters:
  • props (MaterialProperties) -- The material property.

  • t (float or complex, optional) -- The material thickness (keyword argument). Defaults to 1.0.

  • tNum (int, optional) -- Design variable number to assign to thickness (keyword argument). Defaults to -1 (i.e. no design variable).

  • tlb (float or complex, optional) -- Thickness variable lower bound (keyword argument). Defaults to 0.0.

  • tub (float or complex, optional) -- Thickness variable upper bound (keyword argument). Defaults to 10.0.

  • tOffset (float or complex, optional) -- Offset distance of reference plane (where nodes are located) relative to thickness mid-plane. Measured in fraction of shell thickness. A value of 0.5 places the reference plane at the top of the plate, a value of 0.0 at the plate mid-plane, and a value of -0.5 at the bottom of the plate. Defaults to 0.0.

generateBDFCard()

Generate pyNASTRAN card class based on current design variable values.

Returns:

pyNastran card holding property information

Return type:

card (pyNastran.bdf.cards.properties.shell.PSHELL)

class tacs.constitutive.IsoTubeBeamConstitutive

Bases: BeamConstitutive

Timoshenko theory based constitutive object for a hollow circular beam.

Parameters:
  • props (MaterialProperties) -- The material property.

  • d (float or complex, optional) -- Tube inner diameter (keyword argument). Defaults to 1.0.

  • dNum (int, optional) -- Design variable number to assign to tube diameter (keyword argument). Defaults to -1 (i.e. no design variable).

  • dlb (float or complex, optional) -- Lower bound on tube diameter (keyword argument). Defaults to 0.0.

  • dub (float or complex, optional) -- Upper bound on tube diameter (keyword argument). Defaults to 10.0.

  • t (float or complex, optional) -- Tube wall thickness (keyword argument). Defaults to 0.1.

  • tNum (int, optional) -- Design variable number to assign to wall thickness (keyword argument). Defaults to -1 (i.e. no design variable).

  • tlb (float or complex, optional) -- Lower bound on wall thickness (keyword argument). Defaults to 0.0.

  • tub (float or complex, optional) -- Upper bound on wall thickness (keyword argument). Defaults to 10.0.

generateBDFCard()

Generate pyNASTRAN card class based on current design variable values.

Returns:

pyNastran card holding property information

Return type:

card (pyNastran.bdf.cards.properties.bars.PBARL)

class tacs.constitutive.LamParamShellConstitutive

Bases: ShellConstitutive

class tacs.constitutive.PanelGPs

Bases: object

This object holds the individual GP objects for a single panel. One PanelGP object is shared among the constitutive object for each panel, and the input and output buckling loads and derivatives are saved and restored so that only one call to each GP object is required per panel (speeds up computation time). The construction of the TACS callback with the panelGPs is such that only one PanelGPs object is made per TACSComponent (assuming each TACSComponent is associated with a different panel).

The typical construction from an example in ml_buckling repo (in file 4_aob_opt/_gp_callback) is: # now build a dictionary of PanelGP objects which manage the GP for each tacs component/panel

def callback_generator(tacs_component_names):

axialGP = constitutive.BucklingGP.from_csv( csv_file=mlb.axialGP_csv, theta_csv=mlb.axial_theta_csv ) shearGP = constitutive.BucklingGP.from_csv( csv_file=mlb.shearGP_csv, theta_csv=mlb.shear_theta_csv ) panelGP_dict = constitutive.PanelGPs.component_dict( tacs_component_names, axialGP=axialGP, shearGP=shearGP )

def gp_callback(dvNum, compID, compDescript, elemDescripts, specialDVs, **kwargs):

# get the panelGPs object associated with this tacs component panelGPs = panelGP_dict[compDescript]

# ... (after this you build the TACSMaterialProperties objects and TACSGPBladeConstitutive objects.

Parameters:
  • BucklingGP (TACSCripplingGaussianProcessModel, optional) -- GP model for axial buckling data, if None uses closed-form instead

  • BucklingGP -- GP model for shear buckling data, if None uses closed-form instead

  • BucklingGP -- GP model for crippling buckling data, if None uses closed-form instead

  • saveData (bool) -- whether to save and restore input and output buckling predictions for more efficient buckling predictions (default True)

classmethod component_dict(tacs_components, axialGP=None, shearGP=None, cripplingGP=None, saveData=True)
constructs a dictionary of PanelGPs objects, one for each tacs component. The dictionary is of the form:

{ tacs_component (str) : PanelGPs object }

Parameters:
  • BucklingGP (TACSCripplingGaussianProcessModel, optional) -- GP model for axial buckling data, if None uses closed-form instead

  • BucklingGP -- GP model for shear buckling data, if None uses closed-form instead

  • BucklingGP -- GP model for crippling buckling data, if None uses closed-form instead

  • saveData (bool) -- whether to save and restore input and output buckling predictions for more efficient buckling predictions (default True)

class tacs.constitutive.PhaseChangeMaterialConstitutive

Bases: Constitutive

This is the base class for the phase change material constitutive objects.

Parameters:
  • solid_props (MaterialProperties) -- The material property of the solid phase.

  • liquid_props (MaterialProperties) -- The material property of the liquid phase.

  • lh (float or complex) -- The specific latent heat of the material.

  • mt (float or complex) -- The melting temperature of the material.

  • t (float or complex, optional) -- The material thickness (keyword argument). Defaults to 1.0.

  • tNum (int, optional) -- Design variable number to assign to thickness (keyword argument). Defaults to -1 (i.e. no design variable).

  • tlb (float or complex, optional) -- Thickness variable lower bound (keyword argument). Defaults to 0.0.

  • tub (float or complex, optional) -- Thickness variable upper bound (keyword argument). Defaults to 10.0.

class tacs.constitutive.PlaneStressConstitutive

Bases: Constitutive

This is the base class for the plane stress constitutive objects. All objects performing plane stress analysis should utilize this class.

Parameters:
  • props (MaterialProperties) -- The material property.

  • t (float or complex, optional) -- The material thickness (keyword argument). Defaults to 1.0.

  • tNum (int, optional) -- Design variable number to assign to thickness (keyword argument). Defaults to -1 (i.e. no design variable).

  • tlb (float or complex, optional) -- Thickness variable lower bound (keyword argument). Defaults to 0.0.

  • tub (float or complex, optional) -- Thickness variable upper bound (keyword argument). Defaults to 10.0.

class tacs.constitutive.PointMassConstitutive

Bases: GeneralMassConstitutive

This is the base class for the traditional point mass constitutive objects with no translation-rotation coupling. Assumes 6 dofs.

Note

TACS uses a negative sign convention in the product of inertia definition, for example:

I12 = -int[x1 * x2 * dm]

The moments of inertia are always positive, as usual:

I11 = int[(x2^2 + x3^2) * dm]

Parameters:
  • m (float or complex, optional) -- Mass value (keyword argument). Defaults to 0.0.

  • I11 (float or complex, optional) -- Moment of inertia in '1' direction (keyword argument). Defaults to 0.0.

  • I22 (float or complex, optional) -- Moment of inertia in '2' direction (keyword argument). Defaults to 0.0.

  • I33 (float or complex, optional) -- Moment of inertia in '3' direction (keyword argument). Defaults to 0.0.

  • I12 (float or complex, optional) -- Moment of inertia in '1-2' plane (keyword argument). Defaults to 0.0.

  • I13 (float or complex, optional) -- Moment of inertia in '1-3' plane (keyword argument). Defaults to 0.0.

  • I23 (float or complex, optional) -- Moment of inertia in '2-3' plane (keyword argument). Defaults to 0.0.

  • mNum (int, optional) -- Design variable number to assign to m (keyword argument). Defaults to -1.

  • mlb (float or complex, optional) -- Lower bound on m (keyword argument). Defaults to 0.0.

  • mub (float or complex, optional) -- Upper bound on wall thickness (keyword argument). Defaults to 1e20.

  • I11Num (int, optional) -- Design variable number to assign to I11 (keyword argument). Defaults to -1.

  • I11lb (float or complex, optional) -- Lower bound on I11 (keyword argument). Defaults to 0.0.

  • I11ub (float or complex, optional) -- Upper bound on I11 (keyword argument). Defaults to 1e20.

  • I22Num (int, optional) -- Design variable number to assign to I22 (keyword argument). Defaults to -1.

  • I22lb (float or complex, optional) -- Lower bound on I22 (keyword argument). Defaults to 0.0.

  • I22ub (float or complex, optional) -- Upper bound on I22 (keyword argument). Defaults to 1e20.

  • I33Num (int, optional) -- Design variable number to assign to I33 (keyword argument). Defaults to -1.

  • I33lb (float or complex, optional) -- Lower bound on I33 (keyword argument). Defaults to 0.0.

  • I33ub (float or complex, optional) -- Upper bound on I33 (keyword argument). Defaults to 1e20.

  • I12Num (int, optional) -- Design variable number to assign to I12 (keyword argument). Defaults to -1.

  • I12lb (float or complex, optional) -- Lower bound on I12 (keyword argument). Defaults to -1e20.

  • I12ub (float or complex, optional) -- Upper bound on I12 (keyword argument). Defaults to 1e20.

  • I13Num (int, optional) -- Design variable number to assign to I13 (keyword argument). Defaults to -1.

  • I13lb (float or complex, optional) -- Lower bound on I13 (keyword argument). Defaults to -1e20.

  • I13ub (float or complex, optional) -- Upper bound on I13 (keyword argument). Defaults to 1e20.

  • I23Num (int, optional) -- Design variable number to assign to I23 (keyword argument). Defaults to -1.

  • I23lb (float or complex, optional) -- Lower bound on I23 (keyword argument). Defaults to -1e20.

  • I23ub (float or complex, optional) -- Upper bound on I23 (keyword argument). Defaults to 1e20.

class tacs.constitutive.ShellConstitutive

Bases: Constitutive

This is the base class for the shell constitutive objects. All objects performing shell elastic analysis should utilize this class.

setDrillingRegularization(kpenalty=10.0)

Update regularization parameter used to stiffen shell in drilling rotation dof.

Parameters:

kpenalty (float) -- Drilling regularization parameter. Defaults to 10.0.

class tacs.constitutive.SmearedCompositeShellConstitutive

Bases: ShellConstitutive

This constitutive class defines the stiffness properties for a composite laminate first-order shear deformation theory (FSDT) shell type element. The stiffness of the laminate is computed by homogenizing the stiffness properties of each ply through the thickness weighted by their relative ply fractions. This class is a good continuous parametrization for laminates used in gradient-based optimizations where stacking sequence effects can be ignored.

Parameters:
  • ply_list (list[OrthotropicPly]) -- List of ply properties in layup.

  • thicknesses (float or complex) -- Total laminate thickness of layup.

  • ply_angles (numpy.ndarray[float or complex]) -- Array of ply angles (in radians) in layup.

  • ply_fractions (numpy.ndarray[float or complex]) -- Fraction of layup contribution of each ply in ply_list.

  • thickness_dv_num (int, optional) -- Design variable number to assign to thickness (keyword argument). Defaults to -1.

  • ply_fraction_dv_nums (numpy.ndarray[int], optional) -- Design variable numbers to assign to ply fractions (keyword argument). Defaults to -1.

  • thickness_lb (float or complex, optional) -- Lower bound for thickness design variable (keyword argument). Defaults to 0.0.

  • thickness_ub (float or complex, optional) -- Upper bound for thickness design variable (keyword argument). Defaults to 1e20.

  • ply_fraction_lb (numpy.ndarray[float or complex], optional) -- Lower bound for ply fraction design variables (keyword argument). Defaults to 0.0.

  • ply_fraction_ub (numpy.ndarray[float or complex], optional) -- Upper bound for ply fraction design variables (keyword argument). Defaults to 1.0.

  • t_offset (float or complex, optional) -- Offset distance of reference plane (where nodes are located) relative to thickness mid-plane. Measured in fraction of shell thickness. A value of 0.5 places the reference plane at the top of the plate, a value of 0.0 at the plate mid-plane, and a value of -0.5 at the bottom of the plate. Defaults to 0.0.

generateBDFCard()

Generate pyNASTRAN card class based on current design variable values.

Returns:

pyNastran card holding property information

Return type:

card (pyNastran.bdf.cards.properties.shell.PCOMP)

class tacs.constitutive.SolidConstitutive

Bases: Constitutive

This is the base class for the solid constitutive objects. All objects performing solid elastic analysis should utilize this class.

Parameters:
  • props (MaterialProperties) -- The material property.

  • t (float or complex, optional) -- The material "thickness" (keyword argument). Weighting factor used for topology optimization. 0.0 corresponds to void, 1.0 corresponds material fully present, values in between are intermediate. Defaults to 1.0.

  • tNum (int, optional) -- Design variable number to assign to "thickness" (keyword argument). Defaults to -1 (i.e. no design variable).

  • tlb (float or complex, optional) -- "Thickness" variable lower bound (keyword argument). Defaults to 0.0.

  • tub (float or complex, optional) -- "Thickness" variable upper bound (keyword argument). Defaults to 10.0.

generateBDFCard()

Generate pyNASTRAN card class based on current design variable values.

Returns:

pyNastran card holding property information

Return type:

card (pyNastran.bdf.cards.properties.solid.PSOLID)

class tacs.constitutive.StiffenedShellConstitutive

Bases: ShellConstitutive

This is a base class for both the BladeStiffenedShellConstitutive and the GPBladeStiffenedShellConstitutive classes

setFailureModes(includePanelMaterialFailure=None, includeStiffenerMaterialFailure=None, includeLocalBuckling=None, includeGlobalBuckling=None, includeStiffenerColumnBuckling=None, includeStiffenerCrippling=None)

Turn on or off each failure mode in the KS failure index computation. If any of the entries are None, they are not modified and will be included in the failure index.

Parameters:
  • includePanelMaterialFailure (bool or None) -- whether to include panel material / panel stress failure in the failure index

  • includeStiffenerMaterialFailure (bool or None) -- whether to include stiffener material / stiffener stress failure in the failure index

  • includeLocalBuckling (bool or None) -- whether to include local buckling (in between stiffeners) in the failure index

  • includeGlobalBuckling (bool or None) -- whether to include global buckling (modes for the full stiffened panel) in the failure index

  • includeStiffenerColumnBuckling (bool or None) -- whether to include column buckling failure in the failure index

  • includeStiffenerCrippling (bool or None) -- whether to include stiffener crippling failure in the failure index

setKSWeight(ksWeight)

Update the ks weight used for aggregating the different failure modes

Parameters:

ksWeight (float) -- KS aggregation weight

setPanelPlyFractionBounds(lowerBound, upperBound)

Set the lower and upper bounds for the panel ply fraction design variables

The default bounds are 0 and 1

Parameters:
Raises:

ValueError -- Raises error if the length of lowerBound or upperBound is not equal to the number of panel plies

setPanelThicknessBounds(lowerBound, upperBound)

Set the lower and upper bounds for the panel thickness design variable

The default bounds are 1e-4 and 1e20

Parameters:
setStiffenerHeightBounds(lowerBound, upperBound)

Set the lower and upper bounds for the stiffener height design variable

The default bounds are 1e-3 and 1e20

Parameters:
setStiffenerPitchBounds(lowerBound, upperBound)

Set the lower and upper bounds for the stiffener pitch design variable

The default bounds are 1e-3 and 1e20

Parameters:
setStiffenerPlyFractionBounds(lowerBound, upperBound)

Set the lower and upper bounds for the stiffener ply fraction design variables

The default bounds are 0 and 1

Parameters:
Raises:

ValueError -- Raises error if the length of lowerBound or upperBound is not equal to the number of stiffener plies

setStiffenerThicknessBounds(lowerBound, upperBound)

Set the lower and upper bounds for the stiffener thickness design variable

The default bounds are 1e-4 and 1e20

Parameters: