DVGeometryCST

class pygeo.parameterization.DVGeoCST.DVGeometryCST(datFile, numCST=8, idxChord=0, idxVertical=1, comm=mpi4py.MPI.COMM_WORLD, isComplex=False, debug=False, tolTE=60.0, name=None)[source]

This class implements a 2D geometry parameterisation based on Brenda Kulfan’s CST (Class-Shape Transformation) method. This class can work with 3D coordinates but will only change the point coordinates in one direction.

The CST equation is as follows:

\(y(x) = C(x) * S(x) + y_\text{te}x\)

Where C is the class function:

\(C(x) = (x^{N1} + (1 - x)^{N2})\)

And S is the shape function, in this case a summation of Bernstein polynomials:

\(S(x) = \sum_i^n w_i \binom{n}{i}x^i(1-x)^{n-i}\)

Here x is the normalized chordwise coordinate, ranging from 0 to 1 from front to the rear of the shape.

Assumptions about the point sets being added:

  • Dat file is ordered continuously around the airfoil and the beginning and end of the list is the trailing edge (no jumping around, but CW vs. CCW does not matter)

  • Geometry is exclusively an extruded shape (no spanwise changes allowed)

  • Airfoil’s leading edge is on the left (min x) and trailing edge is on the right (max x)

  • Airfoil is not rotated (trailing edge and leading edge are close to y equals zero)

Parameters:
datFilestr

Filename of dat file that represents the initial airfoil. The coordinates in this file will be used to determine the camber line, which is the dividing line to distinguish upper and lower surface points.

numCSTint or list of two ints

Number of CST parameters to use for the initial fit and the DVs (if DVs with type "upper" or "lower" are added). If numCST is an int, the value will be used for both upper and lower. If it is a two-item list, the first value defines the number of upper CST coefficients and the second is the number of lower coefficients, by default 8.

idxChordint, optional

Index of the column in the point set to use as the chordwise (x in CST) coordinates, by default 0

idxVerticalint, optional

Index of the column in the point set to use as the vertical (y in CST) airfoil coordinates, by default 1

commMPI communicator, optional

Communicator for DVGeometryCST instance, by default MPI.COMM_WORLD

isComplexbool, optional

Initialize variables to complex types where necessary, by default False

debugbool, optional

Show plots when addPointSet is called to visually verify that it is correctly splitting the upper and lower surfaces of the airfoil points, by default False

tolTEfloat, optional

Tolerance used to detect trailing edge corners on the airfoil. The value represents the angle difference in degrees between adjacent edges of the airfoil, by default 60 deg.

namestring, optional

Name of this DVGeo object, not necessary unless multiple DVGeos are used in one optimization.

addDV(dvName, dvType, lowerBound=None, upperBound=None, scale=1.0, default=None)[source]

Add design variables to the DVGeometryCST object. For upper and lower CST coefficient DVs, the number of design variables is defined using the numCST parameter in DVGeoCST’s init function.

Parameters:
dvNamestr

A unique name to be given to this design variable group

dvTypestr

Define the type of CST design variable being added. Either the upper/lower surface class shape parameter DV can be defined (e.g., "N1_upper"), or the DV for both the upper and lower surfaces’ class shape parameter can be defined (e.g., "N1"), but not both. The options (not case sensitive) are

  • "upper": upper surface CST coefficients (specify dvNum to define how many)

  • "lower": lower surface CST coefficients (specify dvNum to define how many)

  • "N1": first class shape parameter for both upper and lower surfaces (adds a single DV)

  • "N2": second class shape parameter for both upper and lower surfaces (adds a single DV)

  • "N1_upper": first class shape parameters for upper surface (adds a single DV)

  • "N1_lower": first class shape parameters for lower surface (adds a single DV)

  • "N2_upper": second class shape parameters for upper surface (adds a single DV)

  • "N2_lower": second class shape parameters for lower surface (adds a single DV)

  • "chord": chord length in whatever units the point set length is defined and scaled to keep the leading edge at the same position (adds a single DV)

lowerBoundfloat or ndarray, optional

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

upperBoundfloat or ndarray, optional

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

scalefloat, optional

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.

defaultndarray, optional

Default value for design variable (must be same length as number of DVs added).

Returns:
Nint

The number of design variables added.

addPointSet(points, ptName, boundTol=1e-10, **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.

Note

Even if isComplex=True, the imaginary portion of coordinates passed in here is ignored when determining if a given point is on the upper or lower surface.

Parameters:
pointsarray, size (N,3)

The coordinates to embed.

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.

boundTolfloat, optional

Small absolute deviation by which the airfoil coordinates can exceed the initial minimum and maximum x coordinates, by default 1e-10.

**kwargs

Any other parameters are ignored.

addVariablesPyOpt(optProb)[source]

Add the current set of variables to the optProb object.

Parameters:
optProbpyOpt_optimization class

Optimization problem definition to which variables are added

static computeCSTCoordinates(x, N1, N2, w, yte, dtype=<class 'float'>)[source]

Compute the vertical coordinates of a CST curve.

This function assumes x has been normalized to the range [0,1].

Parameters:
xndarray (# pts,)

x coordinates at which to compute the CST curve height

N1float

First class shape parameter

N2float

Second class shape parameter

wndarray (# coeff,)

CST coefficient array

ytefloat

y coordinate of the trailing edge (used to define trailing edge thickness). Note that the trailing edge will be twice this thick, assuming the same yte value is used for both the upper and lower surfaces.

dtypetype, optional

Type for instantiated arrays, by default float

Returns:
ndarray (# pts,)

y coordinates of the CST curve

static computeCSTdydN1(x, N1, N2, w, dtype=<class 'float'>)[source]

Compute the derivatives of the height of a CST curve with respect to N1

Given \(y = C(x, N1, N2) * S(x)\)

\(\frac{dy}{dN1} = S(x) * \frac{dC}{dN1} = S(x) * C(x, N1, N2) * \ln{x}\)

This function assumes x has been normalised to the range [0,1].

Parameters:
xndarray (# pts,)

x coordinates at which to compute the CST curve height

N1float

First class shape parameter

N2float

Second class shape parameter

wndarray (# coeff,)

CST coefficient array

dtypetype, optional

Type for instantiated arrays, by default float

Returns:
ndarray (# pts,)

Derivative of the y coordinates with respect to the first class shape parameter

static computeCSTdydN2(x, N1, N2, w, dtype=<class 'float'>)[source]

Compute the derivatives of the height of a CST curve with respect to N2

Given \(y = C(x, N1, N2) * S(x)\)

\(\frac{dy}{dN2} = S(x) * \frac{dC}{dN2} = S(x) * C(x, N1, N2) * \ln(1-x)\)

This function assumes x has been normalised to the range [0,1].

Parameters:
xndarray (# pts,)

x coordinates at which to compute the CST curve height

N1float

First class shape parameter

N2float

Second class shape parameter

wndarray (# coeff,)

CST coefficient array

dtypetype, optional

Type for instantiated arrays, by default float

Returns:
ndarray (# pts,)

Derivative of the y coordinates with respect to the second class shape parameter

static computeCSTdydw(x, N1, N2, w, dtype=<class 'float'>)[source]

Compute the derivatives of the height of a CST curve with respect to the shape function coefficients

Given \(y = C(x) * sum [w_i * p_i(x)]\)

\(\frac{dy}{dw_i} = C(x) * p_i(x)\)

This function assumes x has been normalized to the range [0,1].

Only the shape of w is used, not the values.

Parameters:
xndarray (# pts,)

x coordinates at which to compute the CST curve height

N1float

First class shape parameter

N2float

Second class shape parameter

wndarray (# coeff,)

CST coefficient array

dtypetype, optional

Type for instantiated arrays, by default float

Returns:
ndarray (# coeff, # pts)

Derivatives of the y coordinates with respect to the CST coefficients

static computeCSTfromCoords(xCoord, yCoord, nCST, N1=0.5, N2=1.0, dtype=<class 'float'>)[source]

Compute the CST coefficients that fit a set of airfoil coordinates (either for the upper or lower surface, not both).

This function internally normalizes the x and y-coordinates.

Parameters:
xCoordndarray

Upper or lower surface airfoil x-coordinates (same length as yCoord vector).

yCoordndarray

Upper or lower surface airfoil y-coordinates (same length as xCoord vector).

nCSTint

Number of CST coefficients to fit.

N1float, optional

First class shape parameter to assume in fitting, by default 0.5

N2float, optional

Second class shape parameter to assume in fitting, by default 1.0

dtypetype, optional

Type for instantiated arrays, by default float

Returns:
np.ndarray (nCST,)

CST coefficients fit to the airfoil surface.

static computeClassShape(x, N1, N2, dtype=<class 'float'>)[source]

Compute the class shape of a CST curve

Parameters:
xndarray (# pts,)

x coordinates at which to compute the CST curve height

N1float

First class shape parameter

N2float

Second class shape parameter

dtypetype, optional

Type for instantiated arrays, by default float

Returns:
ndarray (# pts,)

y coordinates of the class shape

static computeShapeFunctions(x, w, dtype=<class 'float'>)[source]

Compute the Bernstein polynomial shape function of a CST curve

This function assumes x has been normalized to the range [0,1].

Parameters:
xndarray (# pts,)

x coordinates at which to compute the CST curve height

wndarray (# coeff,)

CST coefficient array

dtypetype, optional

Type for instantiated arrays, by default float

Returns:
ndarray (# coeff, # pts)

Bernstein polynomials for each CST coefficient

getNDV()[source]

Return the total number of design variables this object has.

Returns:
nDVint

Total number of design variables

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 setValues()

Returns:
dvDictdict

Dictionary of design variables

getVarNames(**kwargs)[source]

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

Examples

>>> optProb.addCon(.....wrt=DVGeo.getVarNames())
static plotCST(upperCoeff, lowerCoeff, N1=0.5, N2=1.0, nPts=100, ax=None, **kwargs)[source]

Simple utility to generate a plot from CST coefficients.

Parameters:
upperCoeffndarray

One dimensional array of CST coefficients for the upper surface.

lowerCoeffndarray

One dimensional array of CST coefficients for the lower surface.

N1float

First class shape parameter.

N2float

Second class shape parameter.

nPtsint, optional

Number of coordinates to compute on each surface.

axmatplotlib Axes, optional

Axes on which to plot airfoil.

**kwargs

Keyword arguments passed to matplotlib.pyplot.plot

Returns:
matplotlib Axes

Axes with airfoil plotted

printDesignVariables()[source]

Print a formatted list of design variables to the screen

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, **kwargs)[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 (N, Npt, 3). 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

\*\*kwargs

Any other parameters ignored, but this is maintained to allow the same interface as other DVGeo implementations.

Returns:
dIdxDictdict

The dictionary containing the derivatives, suitable for pyOptSparse

totalSensitivityProd(vec, ptSetName, **kwargs)[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

\*\*kwargs

Any other parameters ignored, but this is maintained to allow the same interface as other DVGeo implementations.

Returns:
xsdotarray (Nx3)

Array with derivative seeds of the surface nodes.

update(ptSetName, **kwargs)[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.

**kwargs

Any other parameters ignored, but this is maintained to allow the same interface as other DVGeo implementations.

Returns:
pointsndarray (N x 3)

Updated point set coordinates.