DVGeometry

class pygeo.parameterization.DVGeo.DVGeometry(fileName, *args, isComplex=False, child=False, faceFreeze=None, name=None, kmax=4, **kwargs)[source]

A class for manipulating geometry.

The purpose of the DVGeometry class is to provide a mapping from user-supplied design variables to an arbitrary set of discrete, three-dimensional coordinates. These three-dimensional coordinates can in general represent anything, but will typically be the surface of an aerodynamic mesh, the nodes of a FE mesh or the nodes of another geometric construct.

In a very general sense, DVGeometry performs two primary functions:

  1. Given a new set of design variables, update the three-dimensional coordinates: \(X_{DV}\rightarrow X_{pt}\) where \(X_{pt}\) are the coordinates and \(X_{DV}\) are the user variables.

  2. Determine the derivative of the coordinates with respect to the design variables. That is the derivative \(\frac{dX_{pt}}{dX_{DV}}\)

DVGeometry uses the Free-Form Deformation approach for geometry manipulation. The basic idea is the coordinates are embedded in a clear-flexible jelly-like block. Then by stretching moving and ‘poking’ the volume, the coordinates that are embedded inside move along with overall deformation of the volume.

Parameters
fileNamestr

filename of FFD file. This must be a ascii formatted plot3D file in fortran ordering.

isComplexbool

Make the entire object complex. This should only be used when debugging the entire tool-chain with the complex step method.

childbool

Flag to indicate that this object is a child of parent DVGeo object

faceFreezedict

A dictionary of lists of strings specifying which faces should be ‘frozen’. Each dictionary represents one block in the FFD. For example if faceFreeze =[‘0’:[‘iLow’],’1’:[]], then the plane of control points corresponding to i=0, and i=1, in block ‘0’ will not be able to move in DVGeometry.

namestr

This is prepended to every DV name for ensuring design variables names are unique to pyOptSparse. Only useful when using multiple DVGeos with addTriangulatedSurfaceConstraint()

kmaxint

Maximum order of the splines used for the underlying formulation. Default is a 4th order spline in each direction if the dimensions allow.

Examples

The general sequence of operations for using DVGeometry is as follows::
>>> from pygeo import *
>>> DVGeo = DVGeometry('FFD_file.fmt')
>>> # Embed a set of coordinates Xpt into the object
>>> DVGeo.addPointSet(Xpt, 'myPoints')
>>> # Associate a 'reference axis' for large-scale manipulation
>>> DVGeo.addRefAxis('wing_axis', axis_curve)
>>> # Define a global design variable function:
>>> def twist(val, geo):
>>>    geo.rot_z['wing_axis'].coef[:] = val[:]
>>> # Now add this as a global variable:
>>> DVGeo.addGlobalDV('wing_twist', 0.0, twist, lower=-10, upper=10)
>>> # Now add local (shape) variables
>>> DVGeo.addLocalDV('shape', lower=-0.5, upper=0.5, axis='y')
>>>
addChild(childDVGeo)[source]

Embed a child FFD into this object.

An FFD child is a ‘sub’ FFD that is fully contained within another, parent FFD. A child FFD is also an instance of DVGeometry which may have its own global and/or local design variables. Coordinates do not need to be added to the children. The parent object will take care of that in a call to addPointSet().

See https://github.com/mdolab/pygeo/issues/7 for a description of an issue with Child FFDs that you should be aware of if you are combining shape changes of a parent FFD with rotation or shape changes of a child FFD.

Parameters
childDVGeoinstance of DVGeometry

DVGeo object to use as a sub-FFD

addCompositeDV(dvName, ptSetName=None, u=None, scale=None)[source]

Add composite DVs. Note that this is essentially a preprocessing call which only works in serial at the moment.

Parameters
dvNamestr

The name of the composite DVs

ptSetNamestr, optional

If the matrices need to be computed, then a point set must be specified, by default None

undarray, optional

The u matrix used for the composite DV, by default None

scalefloat or ndarray, optional

The scaling applied to this DV, by default None

addGlobalDV(dvName, value, func, lower=None, upper=None, scale=1.0, config=None)[source]

Add a global design variable to the DVGeometry object. This type of design variable acts on one or more reference axis.

Parameters
dvNamestr

A unique name to be given to this design variable group

valuefloat, or iterable list of floats

The starting value(s) for the design variable. This parameter may be a single variable or a numpy array (or list) if the function requires more than one variable. The number of variables is determined by the rank (and if rank ==1, the length) of this parameter.

lowerfloat, or iterable list of floats

The lower bound(s) for the variable(s). A single variable is permissable even if an array is given for value. However, if an array is given for lower, it must be the same length as value

funcpython function

The python function handle that will be used to apply the design variable

upperfloat, or iterable list of floats

The upper bound(s) for the variable(s). Same restrictions as lower

scalefloat, or iterable list of floats

The scaling of the variables. A good approximate scale to start with is approximately 1.0/(upper-lower). This gives variables that are of order ~1.0.

configstr or list

Define what configurations this design variable will be applied to Use a string for a single configuration or a list for multiple configurations. The default value of None implies that the design variable applies to ALL configurations.

addLocalDV(dvName, lower=None, upper=None, scale=1.0, axis='y', volList=None, pointSelect=None, config=None)[source]

Add one or more local design variables ot the DVGeometry object. Local variables are used for small shape modifications.

Parameters
dvNamestr

A unique name to be given to this design variable group

lowerfloat

The lower bound for the variable(s). This will be applied to all shape variables

upperfloat

The upper bound for the variable(s). This will be applied to all shape variables

scalefloat

The scaling of the variables. A good approximate scale to start with is approximately 1.0/(upper-lower). This gives variables that are of order ~1.0.

axisstr. Default is y

The coordinate directions to move. Permissible values are x, y and z. If more than one direction is required, use multiple calls to addLocalDV() with different axis values.

volListlist

Use the control points on the volume indices given in volList. You should use pointSelect = None, otherwise this will not work.

pointSelectpointSelect object. Default is None Use a

pointSelect object to select a subset of the total number of control points. See the documentation for the pointSelect class in geo_utils. Using pointSelect discards everything in volList.

configstr or list

Define what configurations this design variable will be applied to Use a string for a single configuration or a list for multiple configurations. The default value of None implies that the design variable applies to ALL configurations.

Returns
Nint

The number of design variables added.

Examples

>>> # Add all variables in FFD as local shape variables
>>> # moving in the y direction, within +/- 1.0 units
>>> DVGeo.addLocalDV('shape_vars', lower=-1.0, upper= 1.0, axis='y')
>>> # As above, but moving in the x and y directions.
>>> nVar = DVGeo.addLocalDV('shape_vars_x', lower=-1.0, upper= 1.0, axis='x')
>>> nVar = DVGeo.addLocalDV('shape_vars_y', lower=-1.0, upper= 1.0, axis='y')
>>> # Create a point select to use: (box from (0,0,0) to (10,0,10) with
>>> # any point projecting into the point along 'y' axis will be selected.
>>> PS = geo_utils.PointSelect(type = 'y', pt1=[0,0,0], pt2=[10, 0, 10])
>>> nVar = DVGeo.addLocalDV('shape_vars', lower=-1.0, upper=1.0, pointSelect=PS)
addLocalSectionDV(dvName, secIndex, lower=None, upper=None, scale=1.0, axis=1, pointSelect=None, volList=None, orient0=None, orient2='svd', config=None)[source]

Add one or more section local design variables to the DVGeometry object. Section local variables are used as an alternative to local variables when it is desirable to deform a cross-section shape within a plane that is consistent with the original cross-section orientation. This is helpful in at least two common scenarios:

  1. The original geometry has cross-sections that are not aligned with

    the global coordinate axes. For instance, with a winglet, we want the shape variables to deform normal to the winglet surface instead of in the x, y, or z directions.

  2. The global design variables cause changes in the geometry that

    rotate the orientation of the original cross-section planes. In this case, we want the shape variables to deform in directions aligned with the rotated cross-section plane, which may not be the x, y, or z directions.

** Warnings **
  • Rotations in an upper level (parent) FFD will not propagate down

    to the lower level FFDs due to limitations of the current implementation.

  • Section local design variables should not be specified at the same

    time as local design variables. This will most likely not result in the desired behavior.

Parameters
dvNamestr

A unique name to be given to this design variable group

lowerfloat

The lower bound for the variable(s). This will be applied to all shape variables

upperfloat

The upper bound for the variable(s). This will be applied to all shape variables

scalefloat

The scaling of the variables. A good approximate scale to start with is approximately 1.0/(upper-lower). This gives variables that are of order ~1.0.

axisint
The coordinate directions to move. Permissible values are

0: longitudinal direction (in section plane) 1: latitudinal direction (in section plane) 2: transverse direction (out of section plane)

If more than one direction is required, use multiple calls to addLocalSectionDV with different axis values.

                    1
                    ^
                    |
o-----o--------o----|----o--------o--------o-----o
|                   |                            |  j
|                   x---------> 0                |  ^
|                  /                             |  |
o-----o--------o--/------o--------o--------o-----o
                 /      ----> i
                /
               2
pointSelectpointSelect object. Default is None

Use a pointSelect object to select a subset of the total number of control points. See the documentation for the pointSelect class in geo_utils. Using pointSelect discards everything in volList. You can create a PointSelect object by using, for instance: >>> PS = geo_utils.PointSelect(type = y, pt1=[0,0,0], pt2=[10, 0, 10]) Check the other PointSelect options in geo_utils.py

volListlist

Use the control points on the volume indices given in volList. If None, all volumes will be included. PointSelect has priority over volList. So if you use PointSelect, the values defined in volList will have no effect.

secIndexchar or list of chars

For each volume, we need to specify along which index we would like to subdivide the volume into sections. Entries in list can be i, j, or k. This index will be designated as the transverse (2) direction in terms of the direction of perturbation for the ‘axis’ parameter.

orient0None, i, j, k, or numpy vector. Default is None.

Although secIndex defines the 2 axis, the 0 and 1 axes are still free to rotate within the section plane. We will choose the orientation of the 0 axis and let 1 be orthogonal. We have three options:

  1. <None> (default) If nothing is prescribed, the 0 direction will

    be the best fit line through the section points. In the case of an airfoil, this would roughly align with the chord.

  2. <i,`j` or k> In this case, the 0 axis will be aligned

    with the mean vector between the FFD edges corresponding to this index. In the ascii art above, if j were given for this option, we would average the vectors between the points on the top and bottom surfaces and project this vector on to the section plane as the 0 axis. If a list is given, each index will be applied to its corresponding volume in volList.

  3. <[x, y, z]> If a numpy vector is given, the 0 axis

    will be aligned with a projection of this vector onto the section plane. If a numpy array of len(volList) x 3 is given, each vector will apply to its corresponding volume.

orient2svd or ffd. Default is svd

How to compute the orientation 2 axis. SVD is the default behaviour and is taken from the svd of the plane points. ffd Uses the vector along the FFD direction of secIndex. This is required to get consistent normals if you have a circular-type FFD when the SVD will swap the normals.

configstr or list

Define what configurations this design variable will be applied to Use a string for a single configuration or a list for multiple configurations. The default value of None implies that the design variable applies to ALL configurations.

Returns
Nint

The number of design variables added.

Examples

>>> # Add all control points in FFD as local shape variables
>>> # moving in the 1 direction, within +/- 1.0 units
>>> DVGeo.addLocalSectionDV('shape_vars', secIndex='k', lower=-1, upper=1, axis=1)
addPointSet(points, ptName, origConfig=True, coord_xfer=None, **kwargs)[source]

Add a set of coordinates to DVGeometry

The is the main way that geometry, in the form of a coordinate list is given to DVGeometry to be manipulated.

Parameters
pointsarray, size (N,3)

The coordinates to embed. These coordinates should all project into the interior of the FFD volume.

ptNamestr

A user supplied name to associate with the set of coordinates. This name will need to be provided when updating the coordinates or when getting the derivatives of the coordinates.

origConfigbool

Flag determine if the coordinates are projected into the undeformed or deformed configuration. This should almost always be True except in circumstances when the user knows exactly what they are doing.

coord_xferfunction

A callback function that performs a coordinate transformation between the DVGeo reference frame and any other reference frame. The DVGeo object uses this routine to apply the coordinate transformation in “forward” and “reverse” directions to go between the two reference frames. Derivatives are also adjusted since they are vectors coming into DVGeo (in the reverse AD mode) and need to be rotated. We have a callback function here that lets the user to do whatever they want with the coordinate transformation. The function must have the first positional argument as the array that is (npt, 3) and the two keyword arguments that must be available are “mode” (“fwd” or “bwd”) and “apply_displacement” (True or False). This function can then be passed to DVGeo through something like ADflow, where the set DVGeo call can be modified as: CFDSolver.setDVGeo(DVGeo, pointSetKwargs={“coord_xfer”: coord_xfer})

An example function is as follows:

def coord_xfer(coords, mode="fwd", apply_displacement=True, **kwargs):
    # given the (npt by 3) array "coords" apply the coordinate transformation.
    # The "fwd" mode implies we go from DVGeo reference frame to the
    # application, e.g. CFD, the "bwd" mode is the opposite;
    # goes from the CFD reference frame back to the DVGeo reference frame.
    # the apply_displacement flag needs to be correctly implemented
    # by the user; the derivatives are also passed through this routine
    # and they only need to be rotated when going between reference frames,
    # and they should NOT be displaced.

    # In summary, all the displacements MUST be within the if apply_displacement == True
    # checks, otherwise the derivatives will be wrong.

    #  Example transfer: The CFD mesh
    # is rotated about the x-axis by 90 degrees with the right hand rule
    # and moved 5 units below (in z) the DVGeo reference.
    # Note that the order of these operations is important.

    # a different rotation matrix can be created during the creation of
    # this function. This is a simple rotation about x-axis.
    # Multiple rotation matrices can be used; the user is completely free
    # with whatever transformations they want to apply here.
    rot_mat = np.array([
        [1, 0, 0],
        [0, 0, -1],
        [0, 1, 0],
    ])

    if mode == "fwd":
        # apply the rotation first
        coords_new = np.dot(coords, rot_mat)

        # then the translation
        if apply_displacement:
            coords_new[:, 2] -= 5
    elif mode == "bwd":
        # apply the operations in reverse
        coords_new = coords.copy()
        if apply_displacement:
            coords_new[:, 2] += 5

        # and the rotation. note the rotation matrix is transposed
        # for switching the direction of rotation
        coords_new = np.dot(coords_new, rot_mat.T)

    return coords_new
addRefAxis(name, curve=None, xFraction=None, yFraction=None, zFraction=None, volumes=None, rotType=5, axis='x', alignIndex=None, rotAxisVar=None, rot0ang=None, rot0axis=[1, 0, 0], includeVols=[], ignoreInd=[], raySize=1.5)[source]

This function is used to add a ‘reference’ axis to the DVGeometry object. Adding a reference axis is only required when ‘global’ design variables are to be used, i.e. variables like span, sweep, chord etc — variables that affect many FFD control points.

There are two different ways that a reference can be specified:

  1. The first is explicitly a pySpline curve object using the keyword argument curve=<curve>.

  2. The second is to specify the xFraction variable. There are a few caveats with the use of this method. First, DVGeometry will try to determine automatically the orientation of the FFD volume. Then, a reference axis will consist of the same number of control points as the number of span-wise sections in the FFD volume and will be oriented in the streamwise (x-direction) according to the xPercent keyword argument.

Parameters
namestr

Name of the reference axis. This name is used in the user-supplied design variable functions to determine what axis operations occur on.

curvepySpline curve object

Supply exactly the desired reference axis

xFractionfloat

Specify the parametric stream-wise (axis: 0) location of the reference axis node relative to front and rear control points location. Constant for every spanwise section.

yFractionfloat

Specify the parametric location of the reference axis node along axis: 1 relative to top and bottom control points location. Constant for every spanwise section. NOTE: if this is the spanwise axis of the FFD box, the refAxis node will remain in-plane and the option will not have any effect.

zFractionfloat

Specify the parametric location of the reference axis node along axis: 2 relative to top and bottom control points location. Constant for every spanwise section. NOTE: if this is the spanwise axis of the FFD box, the refAxis node will remain in-plane and the option will not have any effect.

volumeslist or array or integers

List of the volume indices, in 0-based ordering that this reference axis should manipulate. If xFraction is specified, the volumes argument must contain at most 1 volume. If the volumes is not given, then all volumes are taken.

rotTypeint

Integer in range 0->6 (inclusive) to determine the order that the rotations are made.

  1. Intrinsic rotation, rot_theta is rotation about axis

  2. x-y-z

  3. x-z-y

  4. y-z-x

  5. y-x-z

  6. z-x-y Default (x-streamwise y-up z-out wing)

  7. z-y-x

  8. z-x-y + rot_theta

  9. z-x-y + rotation about section axis (to allow for winglet rotation)

axis: str

Axis along which to project points/control points onto the ref axis. Default is x which will project rays.

alignIndex: str

FFD axis along which the reference axis will lie. Can be i, j, or k. Only necessary when using xFraction.

rotAxisVar: str

If rotType == 8, then you must specify the name of the section local variable which should be used to compute the orientation of the theta rotation.

rot0ang: float

If rotType == 0, defines the offset angle of the (child) FFD with respect to the main system of reference. This is necessary to use the scaling functions scale_x, scale_y, and scale_z with rotType == 0. The axis of rotation is defined by rot0axis.

rot0axis: list

If rotType == 0, defines the rotation axis for the rotation offset of the FFD grid given by rot0ang. The variable has to be a list of 3 floats defining the [x,y,z] components of the axis direction. This is necessary to use the scaling functions scale_x, scale_y, and scale_z with rotType == 0.

includeVolslist

List of additional volumes to add to reference axis after the automatic generation of the ref axis based on the volumes list using xFraction.

ignoreIndlist

List of indices that should be ignored from the volumes that were added to this reference axis. This can be handy if you have a single volume but you want to link different sets of indices to different reference axes.

raySizefloat

Used in projection to find attachment point on reference axis. See full description in pyNetwork.projectRays function doc string. In most cases the default value is sufficient. In the case of highly swept wings its sometimes necessary to increase this value.

Returns
nAxisint

The number of control points on the reference axis.

Notes

One of curve or xFraction must be specified.

Examples

>>> # Simple wing with single volume FFD, reference axis at 1/4 chord:
>>> DVGeo.addRefAxis('wing', xFraction=0.25)
>>> # Multiblock FFD, wing is volume 6.
>>> DVGeo.addRefAxis('wing', xFraction=0.25, volumes=[6])
>>> # Multiblock FFD, multiple volumes attached refAxis
>>> DVGeo.addRefAxis('wing', myCurve, volumes=[2,3,4])
addSpanwiseLocalDV(dvName, spanIndex, axis='y', lower=None, upper=None, scale=1.0, pointSelect=None, volList=None, config=None)[source]

Add one or more spanwise local design variables to the DVGeometry object. Spanwise local variables are alternative form of local shape variables used to apply equal DV changes in a chosen direction. Some scenarios were this could be useful are:

  1. 2D airfoil shape optimization. Because adflow works with 3D meshes, 2D problems are represented my a mesh a single cell wide. Therefor, to change the 2D representation of the airfoil both sides of the mesh must be moved equally. This can be done with the addition of linear constraints on a set of local shape variables, however this approach requires more DVs than necessary (which complicates DV sweeps) and the constaints are only enforced to a tolerance. Using spanwise local design variables insures the airfoil is always correctly represented in the 3D mesh using the correct amount of design variables.

  2. 3D wing optimization with constant airfoil shape. If the initial wing geometry has a constant airfoil shape and constant chord, then spanwise local dvs can be used to change the airfoil shape of the wing while still keeping it constant along the span of the wing.

Parameters
dvNamestr

A unique name to be given to this design variable group

spanIndexstr, (‘i’, ‘j’, ‘k’)

the axis of the FFD along which the DVs are constant all shape variables

axisstr. Default is y

The coordinate directions to move. Permissible values are x, y and z. If more than one direction is required, use multiple calls to addLocalDV with different axis values.

lowerfloat

The lower bound for the variable(s). This will be applied to all shape variables

upperfloat

The upper bound for the variable(s). This will be applied to all shape variables

scalefloat

The scaling of the variables. A good approximate scale to start with is approximately 1.0/(upper-lower). This gives variables that are of order ~1.0.

pointSelectpointSelect object. Default is None Use a

pointSelect object to select a subset of the total number of control points. See the documentation for the pointSelect class in geo_utils. Using pointSelect discards everything in volList.

volListlist

Use the control points on the volume indices given in volList. You should use pointSelect = None, otherwise this will not work.

configstr or list

Define what configurations this design variable will be applied to Use a string for a single configuration or a list for multiple configurations. The default value of None implies that the design variable applies to ALL configurations.

Returns
Nint

The number of design variables added.

Examples

>>> # Add all spanwise local variables
>>> # moving in the y direction, within +/- 0.5 units
>>> DVGeo.addSpanwiseLocalDV("shape", 'k', lower=-0.5, upper=0.5, axis="z", scale=1.0)
addVariablesPyOpt(optProb, globalVars=True, localVars=True, sectionlocalVars=True, spanwiselocalVars=True, ignoreVars=None, freezeVars=None)[source]

Add the current set of variables to the optProb object.

Parameters
optProbpyOpt_optimization class

Optimization problem definition to which variables are added

globalVarsbool

Flag specifying whether global variables are to be added

localVarsbool

Flag specifying whether local variables are to be added

sectionlocalVarsbool

Flag specifying whether section local variables are to be added

spanwiselocalVarsbool

Flag specifying whether spanwiselocal variables are to be added

ignoreVarslist of strings

List of design variables the user DOESN’T want to use as optimization variables.

freezeVarslist of string

List of design variables the user WANTS to add as optimization variables, but to have the lower and upper bounds set at the current variable. This effectively eliminates the variable, but it the variable is still part of the optimization.

applyToChild(iChild)[source]

This function is used to apply the changes in the parent FFD to the child FFD points and child reference axis points.

checkDerivatives(ptSetName)[source]

Run a brute force FD check on ALL design variables

Parameters
ptSetNamestr

name of the point set to check

computeDVJacobian(config=None)[source]

return J_temp for a given config

computeTotalJacobian(ptSetName, config=None)[source]

Return the total point jacobian in CSR format since we need this for TACS

computeTotalJacobianCS(ptSetName, config=None)[source]

Return the total point jacobian in CSR format since we need this for TACS

computeTotalJacobianFD(ptSetName, config=None)[source]

This function takes the total derivative of an objective, I, with respect the points controlled on this processor using FD. We take the transpose prodducts and mpi_allreduce them to get the resulting value on each processor. Note that this function is slow and should eventually be replaced by an analytic version.

convertDictToSensitivity(dIdxDict)[source]

This function performs the reverse operation of convertSensitivityToDict(); it transforms the dictionary back into an array. This function is important for the matrix-free interface.

Parameters
dIdxDictdictionary

Dictionary of information keyed by this object’s design variables

Returns
dIdxarray

Flattened array of length getNDV().

convertSensitivityToDict(dIdx, out1D=False, useCompositeNames=False)[source]

This function takes the result of totalSensitivity and converts it to a dict for use in pyOptSparse

Parameters
dIdxarray

Flattened array of length getNDV(). Generally it comes from a call to totalSensitivity()

out1Dboolean

If true, creates a 1D array in the dictionary instead of 2D. This function is used in the matrix-vector product calculation.

useCompositeNamesboolean

Whether the sensitivity dIdx is with respect to the composite DVs or the original DVGeo DVs. If False, the returned dictionary will have keys corresponding to the original set of geometric DVs. If True, the returned dictionary will have replace those with a single key corresponding to the composite DV name.

Returns
dIdxDictdictionary

Dictionary of the same information keyed by this object’s design variables

demoDesignVars(directory, includeLocal=True, includeGlobal=True, pointSet=None, CFDSolver=None, callBack=None, freq=2)[source]

This function can be used to “test” the design variable parametrization for a given optimization problem. It should be called in the script after DVGeo has been set up. The function will loop through all the design variables and write out a deformed FFD volume for the upper and lower bound of every design variable. It can also write out deformed point sets and surface meshes.

Parameters
directorystr

The directory where the files should be written.

includeLocalboolean

False if you don’t want to include the shape variables.

includeGlobalboolean

False if you don’t want to include global variables.

pointSetstr

Name of the point set to write out. If None, no point set output is generated.

CFDSolverstr

An ADflow instance that will be used to write out deformed surface meshes. In addition to having a DVGeo object, CFDSolver must have an AeroProblem set, for example with CFDSolver.setAeroProblem(ap). If CFDSolver is None, no surface mesh output is generated.

callBackfunction

This allows the user to perform an additional task at each new design variable iteration. The callback function must take two inputs:

  1. the output directory name (str) and

  2. the iteration count (int).

freqint

Number of snapshots to take between the upper and lower bounds of a given variable. If greater than 2, will do a sinusoidal sweep.

extractCoef(axisID)[source]

Extract the coefficients for the selected reference axis. This should be used only inside design variable functions

extractS(axisID)[source]

Extract the parametric positions of the control points. This is usually used in conjunction with extractCoef()

getFlattenedChildren()[source]

Return a flattened list of all DVGeo objects in the family hierarchy.

getLocalIndex(iVol, comp=None)[source]

Return the local index mapping that points to the global coefficient list for a given volume

getNDV()[source]

Return the total number of design variables this object has.

Returns
nDVint

Total number of design variables

getSymmetricCoefList(volList=None, pointSelect=None, tol=1e-08, getSymmPlane=False)[source]

Determine the pairs of coefs that need to be constrained for symmetry.

Parameters
volListlist

Use the control points on the volume indices given in volList

pointSelectpointSelect object. Default is None Use a

pointSelect object to select a subset of the total number of control points. See the documentation for the pointSelect class in geo_utils.

tolfloat

Tolerance for ignoring nodes around the symmetry plane. These should be merged by the network/connectivity anyway

getSymmPlanebool

If this flag is set to True, we also return the points on the symmetry plane for all volumes. e.g. a reduced point on the symmetry plane with the same indices on both volumes will show up as the same value in both arrays. This is useful when determining the indices to ignore when adding pointsets. The default behavior will not include the points exactly on the symmetry plane. this is more useful for adding them as linear constraints

Returns
indSetAlist of ints

One half of the coefs to be constrained

indSetBlist of ints

Other half of the coefs to be constrained

getValues()[source]

Generic routine to return the current set of design variables. Values are returned in a dictionary format that would be suitable for a subsequent call to setDesignVars()

Returns
dvDictdict

Dictionary of design variables

getVarNames(pyOptSparse=False)[source]

Return a list of the design variable names. This is typically used when specifying a wrt= argument for pyOptSparse.

Parameters
pyOptSparsebool

Flag to specify whether the DVs returned should be those in the optProb or those internal to DVGeo. Only relevant if using composite DVs.

Examples

optProb.addCon(…..wrt=DVGeo.getVarNames())

mapSensToComp(inVec)[source]

Maps the sensitivity matrix to the composite design space

Parameters
inVecndarray

The sensitivities to be mapped

Returns
ndarray

The mapped sensitivity matrix

mapVecToComp(inVec)[source]

This is the vector version of mapXDictToComp(), where the actual mapping is done

Parameters
inVecndarray

The DVs in a single 1D array

Returns
ndarray

The mapped DVs in a single 1D array

mapVecToDVGeo(inVec)[source]

This is the vector version of mapXDictToDVGeo(), where the actual mapping is done

Parameters
inVecndarray

The DVs in a single 1D array

Returns
ndarray

The mapped DVs in a single 1D array

mapXDictToComp(inDict)[source]

The inverse of mapXDictToDVGeo(), where we map the DVs to the composite space

Parameters
inDictdict

The DVs to be mapped

Returns
dict

The mapped DVs

mapXDictToDVGeo(inDict)[source]

Map a dictionary of DVs to the ‘DVGeo’ design, while keeping non-DVGeo DVs in place without modifying them

Parameters
inDictdict

The dictionary of DVs to be mapped

Returns
dict

The mapped DVs in the same dictionary format

printDesignVariables()[source]

Print a formatted list of design variables to the screen

restoreCoef(coef, axisID)[source]

Restore the coefficients for the selected reference axis. This should be used inside design variable functions

sectionFrame(sectionIndex, sectionTransform, sectionLink, ivol=0, orient0=None, orient2='svd')[source]

This function computes a unique reference coordinate frame for each section of an FFD volume. You can choose which axis of the FFD you would like these sections to be defined by. For example, if we have a wing with a winglet, the airfoil sections which make up the wing will not all lie in parallel planes. We want to find a reference frame for each of these airfoil sections so that we can constrain local control points to deform within the sectional plane. Let’s say the wing FFD is oriented with indices:

i

along chord

j

normal to wing surface

k

along span

If we choose sectionIndex=’k’, this function will compute a frame which has two axes aligned with the k-planes of the FFD volume. This is useful because in some cases (as with a winglet), we want to perturb sectional control points within the section plane instead of in the global coordinate directions.

Assumptions:

  • the normal direction is computed along the block index with size 2

  • all point for a given sectionIndex lie within a plane

Parameters
sectionIndexi, j, or k

This the index of the FFD which defines a section plane.

orient0None, i, j, k, or numpy vector. Default is None.

Although secIndex defines the ‘2’ axis, the ‘0’ and ‘1’ axes are still free to rotate within the section plane. We will choose the orientation of the ‘0’ axis and let ‘1’ be orthogonal. See addLocalSectionDV for a more detailed description.

ivolinteger

Volume ID for the volume in which section normals will be computed.

alignStreamwisex, y, or z (optional)

If given, section frames are rotated about the k-plane normal so that the longitudinal axis is parallel with the given streamwise direction.

rootGloballist

List of sections along specified axis that will be fixed to the global coordinate frame.

Returns
sectionTransformlist of 3x3 arrays

List of transformation matrices for the sections of a given volume. Transformations are set up from local section frame to global frame.

setDesignVars(dvDict)[source]

Standard routine for setting design variables from a design variable dictionary.

Parameters
dvDictdict

Dictionary of design variables. The keys of the dictionary must correspond to the design variable names. Any additional keys in the dictionary are simply ignored.

totalSensitivity(dIdpt, ptSetName, comm=None, config=None)[source]

This function computes sensitivity information.

Specifically, it computes the following: \(\frac{dX_{pt}}{dX_{DV}}^T \frac{dI}{d_{pt}}\)

Parameters
dIdptarray of size (Npt, 3) or (N, Npt, 3)

This is the total derivative of the objective or function of interest with respect to the coordinates in ‘ptSetName’. This can be a single array of size (Npt, 3) or a group of N vectors of size (Npt, 3, N). If you have many to do, it is faster to do many at once.

ptSetNamestr

The name of set of points we are dealing with

commMPI.IntraComm

The communicator to use to reduce the final derivative. If comm is None, no reduction takes place.

configstr or list

Define what configurations this design variable will be applied to Use a string for a single configuration or a list for multiple configurations. The default value of None implies that the design variable applies to ALL configurations.

Returns
dIdxDictdic

The dictionary containing the derivatives, suitable for pyOptSparse

Notes

The child and nDVStore options are only used internally and should not be changed by the user.

totalSensitivityProd(vec, ptSetName, config=None)[source]

This function computes sensitivity information.

Specifically, it computes the following: \(\frac{dX_{pt}}{dX_{DV}} \times\mathrm{vec}\)

This is useful for forward AD mode.

Parameters
vecdictionary whose keys are the design variable names, and whose

values are the derivative seeds of the corresponding design variable.

ptSetNamestr

The name of set of points we are dealing with

configstr or list

Define what configurations this design variable will be applied to Use a string for a single configuration or a list for multiple configurations. The default value of None implies that the design variable applies to ALL configurations.

Returns
xsdotarray (Nx3) -> Array with derivative seeds of the surface nodes.
totalSensitivityTransProd(vec, ptSetName, config=None)[source]

This function computes sensitivity information.

Specifically, it computes the following: \(\frac{dX_{pt}}{dX_{DV}}^T \times\mathrm{vec}\)

This is useful for reverse AD mode.

Parameters
dIdptarray of size (Npt, 3) or (N, Npt, 3)

This is the total derivative of the objective or function of interest with respect to the coordinates in ‘ptSetName’. This can be a single array of size (Npt, 3) or a group of N vectors of size (Npt, 3, N). If you have many to do, it is faster to do many at once.

ptSetNamestr

The name of set of points we are dealing with

commMPI.IntraComm

The communicator to use to reduce the final derivative. If comm is None, no reduction takes place.

configstr or list

Define what configurations this design variable will be applied to Use a string for a single configuration or a list for multiple configurations. The default value of None implies that the design variable applies to ALL configurations.

Returns
dIdxDictdic

The dictionary containing the derivatives, suitable for pyOptSparse

Notes

The child and nDVStore options are only used internally and should not be changed by the user.

update(ptSetName, childDelta=True, config=None)[source]

This is the main routine for returning coordinates that have been updated by design variables.

Parameters
ptSetNamestr

Name of point-set to return. This must match ones of the given in an addPointSet() call.

childDeltabool

Return updates on child as a delta. The user should not need to ever change this parameter.

configstr or list

Define what configurations this design variable will be applied to Use a string for a single configuration or a list for multiple configurations. The default value of None implies that the design variable applies to ALL configurations.

updateCalculations(new_pts, isComplex, config)[source]

The core update routine. pulled out here to eliminate duplication between update and update_deriv.

updatePyGeo(geo, outputType, fileName, nRefU=0, nRefV=0)[source]

Deform a pyGeo object and write to a file of specified type given the (deformed) current state of the FFD object.

Parameters
geopyGeo object

A pyGeo object containing an initialized object

outputType: str

Type of output file to be written. Can be iges or tecplot

fileName: str

Filename for the output file. Should have no extension, an extension will be added

nRefU: int or list of ints

Number of spline refinement points to add in the surface B-Spline u-direction. If scalar, it is applied across each surface. If list, the length must match the number of surfaces in the object and corresponding entries are matched with surfaces.

nRefV: int or list of ints

Number of spline refinement points to add in the surface B-Spline v-direction. If scalar, it is applied across each surface. If list, the length must match the number of surfaces in the object and corresponding entries are matched with surfaces

Write the links attaching the control points to the reference axes

Parameters
fileNamestr

Filename for tecplot file. Should have .dat extension

writePlot3d(fileName)[source]

Write the (deformed) current state of the FFD object into a plot3D file. This file could then be used as the base-line FFD for a subsequent optimization. This function is not typically used in a regular basis, but may be useful in certain situaions, i.e. a sequence of optimizations

Parameters
fileNamestr

Filename of the plot3D file to write. Should have a .fmt file extension.

writePointSet(name, fileName, solutionTime=None)[source]

Write a given point set to a tecplot file

Parameters
namestr

The name of the point set to write to a file

fileNamestr

Filename for tecplot file. Should have no extension, an extension will be added

SolutionTimefloat

Solution time to write to the file. This could be a fictitious time to make visualization easier in tecplot.

writeRefAxes(fileName)[source]

Write the (deformed) current state of the RefAxes to a tecplot file, including the children

Parameters
fileNamestr

Filename for tecplot file. Should have a no extension,an extension will be added.

writeTecplot(fileName, solutionTime=None)[source]

Write the (deformed) current state of the FFD’s to a tecplot file, including the children

Parameters
fileNamestr

Filename for tecplot file. Should have a .dat extension

SolutionTimefloat

Solution time to write to the file. This could be a fictitious time to make visualization easier in tecplot.

zeroJacobians(ptSetNames)[source]

set stored jacobians to None for ptSetNames

Parameters
ptSetNameslist

list of ptSetNames to zero the jacobians.