OpenVSP API Documentation  3.38.0
Functions
XSec and Airfoil Functions

This group of functions provides API control of cross-sections (XSecs). Airfoils are a type of XSec included in this group as well. API functions for Body of Revolution XSecs are included in the Specialized Geometry group. Click here to return to the main page. More...

Functions

void vsp::CutXSec (const std::string &geom_id, int index)
 
void vsp::CopyXSec (const std::string &geom_id, int index)
 
void vsp::PasteXSec (const std::string &geom_id, int index)
 
void vsp::InsertXSec (const std::string &geom_id, int index, int type)
 
int vsp::GetXSecShape (const std::string &xsec_id)
 
double vsp::GetXSecWidth (const std::string &xsec_id)
 
double vsp::GetXSecHeight (const std::string &xsec_id)
 
void vsp::SetXSecWidthHeight (const std::string &xsec_id, double w, double h)
 
void vsp::SetXSecWidth (const std::string &xsec_id, double w)
 
void vsp::SetXSecHeight (const std::string &xsec_id, double h)
 
std::vector< std::string > vsp::GetXSecParmIDs (const std::string &xsec_id)
 
std::string vsp::GetXSecParm (const std::string &xsec_id, const std::string &name)
 
std::vector< vec3dvsp::ReadFileXSec (const std::string &xsec_id, const std::string &file_name)
 
void vsp::SetXSecPnts (const std::string &xsec_id, std::vector< vec3d > &pnt_vec)
 
vec3d vsp::ComputeXSecPnt (const std::string &xsec_id, double fract)
 
vec3d vsp::ComputeXSecTan (const std::string &xsec_id, double fract)
 
void vsp::ResetXSecSkinParms (const std::string &xsec_id)
 
void vsp::SetXSecContinuity (const std::string &xsec_id, int cx)
 
void vsp::SetXSecTanAngles (const std::string &xsec_id, int side, double top, double right, double bottom, double left)
 
void vsp::SetXSecTanSlews (const std::string &xsec_id, int side, double top, double right, double bottom, double left)
 
void vsp::SetXSecTanStrengths (const std::string &xsec_id, int side, double top, double right, double bottom, double left)
 
void vsp::SetXSecCurvatures (const std::string &xsec_id, int side, double top, double right, double bottom, double left)
 
void vsp::ReadFileAirfoil (const std::string &xsec_id, const std::string &file_name)
 
void vsp::SetAirfoilUpperPnts (const std::string &xsec_id, const std::vector< vec3d > &up_pnt_vec)
 
void vsp::SetAirfoilLowerPnts (const std::string &xsec_id, const std::vector< vec3d > &low_pnt_vec)
 
void vsp::SetAirfoilPnts (const std::string &xsec_id, const std::vector< vec3d > &up_pnt_vec, const std::vector< vec3d > &low_pnt_vec)
 
std::vector< vec3dvsp::GetHersheyBarLiftDist (const int &npts, const double &alpha, const double &Vinf, const double &span, bool full_span_flag=false)
 
std::vector< vec3dvsp::GetHersheyBarDragDist (const int &npts, const double &alpha, const double &Vinf, const double &span, bool full_span_flag=false)
 
std::vector< vec3dvsp::GetVKTAirfoilPnts (const int &npts, const double &alpha, const double &epsilon, const double &kappa, const double &tau)
 
std::vector< double > vsp::GetVKTAirfoilCpDist (const double &alpha, const double &epsilon, const double &kappa, const double &tau, std::vector< vec3d > xyz_data)
 
std::vector< vec3dvsp::GetEllipsoidSurfPnts (const vec3d &center, const vec3d &abc_rad, int u_npts=20, int w_npts=20)
 
std::vector< vec3dvsp::GetFeatureLinePnts (const string &geom_id)
 
std::vector< double > vsp::GetEllipsoidCpDist (const std::vector< vec3d > &surf_pnt_vec, const vec3d &abc_rad, const vec3d &V_inf)
 
std::vector< vec3dvsp::GetAirfoilUpperPnts (const std::string &xsec_id)
 
std::vector< vec3dvsp::GetAirfoilLowerPnts (const std::string &xsec_id)
 
std::vector< double > vsp::GetUpperCSTCoefs (const std::string &xsec_id)
 
std::vector< double > vsp::GetLowerCSTCoefs (const std::string &xsec_id)
 
int vsp::GetUpperCSTDegree (const std::string &xsec_id)
 
int vsp::GetLowerCSTDegree (const std::string &xsec_id)
 
void vsp::SetUpperCST (const std::string &xsec_id, int deg, const std::vector< double > &coefs)
 
void vsp::SetLowerCST (const std::string &xsec_id, int deg, const std::vector< double > &coefs)
 
void vsp::PromoteCSTUpper (const std::string &xsec_id)
 
void vsp::PromoteCSTLower (const std::string &xsec_id)
 
void vsp::DemoteCSTUpper (const std::string &xsec_id)
 
void vsp::DemoteCSTLower (const std::string &xsec_id)
 
void vsp::FitAfCST (const std::string &xsec_surf_id, int xsec_index, int deg)
 
void vsp::WriteBezierAirfoil (const std::string &file_name, const std::string &geom_id, const double &foilsurf_u)
 
void vsp::WriteSeligAirfoil (const std::string &file_name, const std::string &geom_id, const double &foilsurf_u)
 
std::vector< vec3dvsp::GetAirfoilCoordinates (const std::string &geom_id, const double &foilsurf_u)
 

Detailed Description

Function Documentation

◆ ComputeXSecPnt()

vec3d vsp::ComputeXSecPnt ( const std::string &  xsec_id,
double  fract 
)

Compute 3D coordinate for a point on an XSec curve given the parameter value (U) along the curve

//==== Add Geom ====//
string stack_id = AddGeom( "STACK" );
//==== Get The XSec Surf ====//
string xsec_surf = GetXSecSurf( stack_id, 0 );
string xsec = GetXSec( xsec_surf, 2 );
double u_fract = 0.25;
vec3d pnt = ComputeXSecPnt( xsec, u_fract );
Definition: Vec3d.h:235
std::string AddGeom(const std::string &type, const std::string &parent=std::string())
vec3d ComputeXSecPnt(const std::string &xsec_id, double fract)
std::string GetXSecSurf(const std::string &geom_id, int index)
std::string GetXSec(const std::string &xsec_surf_id, int xsec_index)
Parameters
[in]xsec_idstring XSec ID
[in]fractdouble Curve parameter value (range: 0 - 1)
Returns
vec3d 3D coordinate point

◆ ComputeXSecTan()

vec3d vsp::ComputeXSecTan ( const std::string &  xsec_id,
double  fract 
)

Compute the tangent vector of a point on an XSec curve given the parameter value (U) along the curve

//==== Add Geom ====//
string stack_id = AddGeom( "STACK" );
//==== Get The XSec Surf ====//
string xsec_surf = GetXSecSurf( stack_id, 0 );
string xsec = GetXSec( xsec_surf, 2 );
double u_fract = 0.25;
vec3d tan = ComputeXSecTan( xsec, u_fract );
vec3d ComputeXSecTan(const std::string &xsec_id, double fract)
Parameters
[in]xsec_idstring XSec ID
[in]fractdouble Curve parameter value (range: 0 - 1)
Returns
vec3d Tangent vector

◆ CopyXSec()

void vsp::CopyXSec ( const std::string &  geom_id,
int  index 
)

Copy a cross-section from the specified geometry and maintain it in memory

// Add Stack
string sid = AddGeom( "STACK", "" );
// Copy XSec To Clipboard
CopyXSec( sid, 1 );
// Paste To XSec 3
PasteXSec( sid, 3 );
void PasteXSec(const std::string &geom_id, int index)
void CopyXSec(const std::string &geom_id, int index)
See also
PasteXSec
Parameters
[in]geom_idstring Geom ID
[in]indexXSec index

◆ CutXSec()

void vsp::CutXSec ( const std::string &  geom_id,
int  index 
)

Cut a cross-section from the specified geometry and maintain it in memory

string fid = AddGeom( "FUSELAGE", "" ); // Add Fuselage
//==== Insert, Cut, Paste Example ====//
InsertXSec( fid, 1, XS_ROUNDED_RECTANGLE ); // Insert A Cross-Section
CopyXSec( fid, 2 ); // Copy Just Created XSec To Clipboard
PasteXSec( fid, 1 ); // Paste Clipboard
CutXSec( fid, 2 ); // Cut Created XSec
@ XS_ROUNDED_RECTANGLE
Definition: APIDefines.h:1383
void InsertXSec(const std::string &geom_id, int index, int type)
void CutXSec(const std::string &geom_id, int index)
See also
PasteXSec
Parameters
[in]geom_idstring Geom ID
[in]indexXSec index

◆ DemoteCSTLower()

void vsp::DemoteCSTLower ( const std::string &  xsec_id)

Demote the CST for the lower airfoil surface. The XSec must be of type XS_CST_AIRFOIL

See also
GetLowerCSTDegree
Parameters
[in]xsec_idXSec ID

◆ DemoteCSTUpper()

void vsp::DemoteCSTUpper ( const std::string &  xsec_id)

Demote the CST for the upper airfoil surface. The XSec must be of type XS_CST_AIRFOIL

See also
GetUpperCSTDegree
Parameters
[in]xsec_idXSec ID

◆ FitAfCST()

void vsp::FitAfCST ( const std::string &  xsec_surf_id,
int  xsec_index,
int  deg 
)

Fit a CST airfoil for an existing airfoil of type XS_FOUR_SERIES, XS_SIX_SERIES, XS_FOUR_DIGIT_MOD, XS_FIVE_DIGIT, XS_FIVE_DIGIT_MOD, XS_ONE_SIX_SERIES, or XS_FILE_AIRFOIL.

Parameters
[in]xsec_surf_idXsecSurf ID
[in]xsec_indexXSec index
[in]degCST degree

◆ GetAirfoilCoordinates()

std::vector< vec3d > vsp::GetAirfoilCoordinates ( const std::string &  geom_id,
const double &  foilsurf_u 
)

Get the untwisted unit-length 2D coordinate points for the specified airfoil

See also
WriteSeligAirfoil
Parameters
[in]geom_idstring Geom ID
[in]foilsurf_uU location (range: 0 - 1) along the surface. The foil surface does not include root and tip caps (i.e. 2 section wing -> XSec0 @ u=0, XSec1 @ u=0.5, XSec2 @ u=1.0)

◆ GetAirfoilLowerPnts()

std::vector<vec3d> vsp::GetAirfoilLowerPnts ( const std::string &  xsec_id)

Get the coordinate points for the lower surface of an airfoil. The XSec must be of type XS_FILE_AIRFOIL

// Add Fuselage Geom
string fuseid = AddGeom( "FUSELAGE", "" );
string xsec_surf = GetXSecSurf( fuseid, 0 );
string xsec = GetXSec( xsec_surf, 1 );
ReadFileAirfoil( xsec, "airfoil/N0012_VSP.af" );
array< vec3d > @low_array = GetAirfoilLowerPnts( xsec );
@ XS_FILE_AIRFOIL
Definition: APIDefines.h:1391
std::vector< vec3d > GetAirfoilLowerPnts(const std::string &xsec_id)
void ReadFileAirfoil(const std::string &xsec_id, const std::string &file_name)
void ChangeXSecShape(const std::string &xsec_surf_id, int xsec_index, int type)
See also
SetAirfoilPnts
Parameters
[in]xsec_idstring XSec ID
Returns
vector<vec3d> Vector of coordinate points for the lower airfoil surface

◆ GetAirfoilUpperPnts()

std::vector<vec3d> vsp::GetAirfoilUpperPnts ( const std::string &  xsec_id)

Get the coordinate points for the upper surface of an airfoil. The XSec must be of type XS_FILE_AIRFOIL

// Add Fuselage Geom
string fuseid = AddGeom( "FUSELAGE", "" );
string xsec_surf = GetXSecSurf( fuseid, 0 );
string xsec = GetXSec( xsec_surf, 1 );
ReadFileAirfoil( xsec, "airfoil/N0012_VSP.af" );
array< vec3d > @up_array = GetAirfoilUpperPnts( xsec );
std::vector< vec3d > GetAirfoilUpperPnts(const std::string &xsec_id)
See also
SetAirfoilPnts
Parameters
[in]xsec_idstring XSec ID
Returns
vector<vec3d> VectorArray of coordinate points for the upper airfoil surface

◆ GetEllipsoidCpDist()

std::vector<double> vsp::GetEllipsoidCpDist ( const std::vector< vec3d > &  surf_pnt_vec,
const vec3d abc_rad,
const vec3d V_inf 
)

Generate Analytical Solution for Potential Flow for specified ellipsoid shape at input surface points for input velocity vector. Based on Munk, M. M., 'Remarks on the Pressure Distribution over the Surface of an Ellipsoid, Moving Translationally Through a Perfect Fluid,' NACA TN-196, June 1924. Function initially created to compare VSPAERO results to theory.

const double pi = 3.14159265358979323846;
const int npts = 101;
const vec3d abc_rad = vec3d(1.0, 2.0, 3.0);
const double alpha = 5; // deg
const double beta = 5; // deg
const double V_inf = 100.0;
array < vec3d > x_slice_pnt_vec(npts);
array<double> theta_vec(npts);
theta_vec[0] = 0;
for ( int i = 1; i < npts; i++ )
{
theta_vec[i] = theta_vec[i-1] + (2 * pi / ( npts - 1) );
}
for ( int i = 0; i < npts; i++ )
{
x_slice_pnt_vec[i] = vec3d( 0, abc_rad[1] * cos( theta_vec[i] ), abc_rad[2] *sin( theta_vec[i] ) );
}
vec3d V_vec = vec3d( ( V_inf * cos( Deg2Rad( alpha ) ) * cos( Deg2Rad( beta ) ) ), ( V_inf * sin( Deg2Rad( beta ) ) ), ( V_inf * sin( Deg2Rad( alpha ) ) * cos( Deg2Rad( beta ) ) ) );
array < double > cp_dist = GetEllipsoidCpDist( x_slice_pnt_vec, abc_rad, V_vec );
std::vector< double > GetEllipsoidCpDist(const std::vector< vec3d > &surf_pnt_vec, const vec3d &abc_rad, const vec3d &V_inf)
See also
GetEllipsoidSurfPnts
Parameters
[in]surf_pnt_vecvector<vec3d> Vector of points on the ellipsoid surface to assess
[in]abc_radvec3d Radius along the A (X), B (Y), and C (Z) axes
[in]V_infvec3d 3D components of freestream velocity
Returns
vector<double> Vector of Cp results corresponding to each point in surf_pnt_arr

◆ GetEllipsoidSurfPnts()

std::vector<vec3d> vsp::GetEllipsoidSurfPnts ( const vec3d center,
const vec3d abc_rad,
int  u_npts = 20,
int  w_npts = 20 
)

Generate the surface coordinate points for a ellipsoid at specified center of input radius along each axis. Based on the MATLAB function ellipsoid (https://in.mathworks.com/help/matlab/ref/ellipsoid.html).

See also
GetVKTAirfoilPnts
Parameters
[in]center3D location of the ellipsoid center
[in]abc_radRadius along the A (X), B (Y), and C (Z) axes
[in]u_nptsNumber of points in the U direction
[in]w_nptsNumber of points in the W direction
Returns
Array of coordinates describing the ellipsoid surface

◆ GetFeatureLinePnts()

std::vector<vec3d> vsp::GetFeatureLinePnts ( const string &  geom_id)

Get the points along the feature lines of a particular Geom

Parameters
[in]geom_idstring Geom ID
Returns
Array of points along the Geom's feature lines

◆ GetHersheyBarDragDist()

std::vector<vec3d> vsp::GetHersheyBarDragDist ( const int &  npts,
const double &  alpha,
const double &  Vinf,
const double &  span,
bool  full_span_flag = false 
)

Get the theoretical drag (Cd) distribution for a Hershey Bar wing with unit chord length using Glauert's Method. This function was initially created to compare VSPAERO results to Lifting Line Theory. If full_span_flag is set to true symmetry is applied to the results.

// Compute theoretical lift and drag distributions using 100 points
double Vinf = 100;
double halfAR = 20;
double alpha_deg = 10;
int n_pts = 100;
array<vec3d> cl_dist_theo = GetHersheyBarLiftDist( int( n_pts ), Deg2Rad( alpha_deg ), Vinf, ( 2 * halfAR ), false );
array<vec3d> cd_dist_theo = GetHersheyBarDragDist( int( n_pts ), Deg2Rad( alpha_deg ), Vinf, ( 2 * halfAR ), false );
std::vector< vec3d > GetHersheyBarLiftDist(const int &npts, const double &alpha, const double &Vinf, const double &span, bool full_span_flag=false)
std::vector< vec3d > GetHersheyBarDragDist(const int &npts, const double &alpha, const double &Vinf, const double &span, bool full_span_flag=false)
Parameters
[in]nptsNumber of points along the span to assess
[in]alphaWing angle of attack (Radians)
[in]VinfFreestream velocity
[in]spanHershey Bar full-span. Note, only half is used in the calculation
[in]full_span_flagFlag to apply symmetry to results (default: false)
Returns
Theoretical coefficient of drag distribution array (size = 2*npts if full_span_flag = true)

◆ GetHersheyBarLiftDist()

std::vector<vec3d> vsp::GetHersheyBarLiftDist ( const int &  npts,
const double &  alpha,
const double &  Vinf,
const double &  span,
bool  full_span_flag = false 
)

Get the theoretical lift (Cl) distribution for a Hershey Bar wing with unit chord length using Glauert's Method. This function was initially created to compare VSPAERO results to Lifting Line Theory. If full_span_flag is set to true symmetry is applied to the results.

// Compute theoretical lift and drag distributions using 100 points
double Vinf = 100;
double halfAR = 20;
double alpha_deg = 10;
int n_pts = 100;
array<vec3d> cl_dist_theo = GetHersheyBarLiftDist( int( n_pts ), Deg2Rad( alpha_deg ), Vinf, ( 2 * halfAR ), false );
array<vec3d> cd_dist_theo = GetHersheyBarDragDist( int( n_pts ), Deg2Rad( alpha_deg ), Vinf, ( 2 * halfAR ), false );
Parameters
[in]nptsNumber of points along the span to assess
[in]alphaWing angle of attack (Radians)
[in]VinfFreestream velocity
[in]spanHershey Bar full-span. Note, only half is used in the calculation
[in]full_span_flagFlag to apply symmetry to results
Returns
Theoretical coefficient of lift distribution array (size = 2*npts if full_span_flag = true)

◆ GetLowerCSTCoefs()

std::vector<double> vsp::GetLowerCSTCoefs ( const std::string &  xsec_id)

Get the CST coefficients for the lower surface of an airfoil. The XSec must be of type XS_CST_AIRFOIL

See also
SetLowerCST
Parameters
[in]xsec_idstring XSec ID
Returns
vector<double> Vector of CST coefficients for the lower airfoil surface

◆ GetLowerCSTDegree()

int vsp::GetLowerCSTDegree ( const std::string &  xsec_id)

Get the CST degree for the lower surface of an airfoil. The XSec must be of type XS_CST_AIRFOIL

See also
SetLowerCST
Parameters
[in]xsec_idXSec ID
Returns
int CST Degree for lower airfoil surface

◆ GetUpperCSTCoefs()

std::vector<double> vsp::GetUpperCSTCoefs ( const std::string &  xsec_id)

Get the CST coefficients for the upper surface of an airfoil. The XSec must be of type XS_CST_AIRFOIL

See also
SetUpperCST
Parameters
[in]xsec_idstring XSec ID
Returns
vector<double> Vector of CST coefficients for the upper airfoil surface

◆ GetUpperCSTDegree()

int vsp::GetUpperCSTDegree ( const std::string &  xsec_id)

Get the CST degree for the upper surface of an airfoil. The XSec must be of type XS_CST_AIRFOIL

See also
SetUpperCST
Parameters
[in]xsec_idstring XSec ID
Returns
int CST Degree for upper airfoil surface

◆ GetVKTAirfoilCpDist()

std::vector<double> vsp::GetVKTAirfoilCpDist ( const double &  alpha,
const double &  epsilon,
const double &  kappa,
const double &  tau,
std::vector< vec3d xyz_data 
)

Get the pressure coefficient (Cp) along a Von Kármán-Trefftz airfoil of specified shape at specified points along the airfoil

const double pi = 3.14159265358979323846;
const int npts = 122;
const double alpha = 0.0;
const double epsilon = 0.1;
const double kappa = 0.1;
const double tau = 10;
array<vec3d> xyz_airfoil = GetVKTAirfoilPnts(npts, alpha, epsilon, kappa, tau*(pi/180) );
array<double> cp_dist = GetVKTAirfoilCpDist( alpha, epsilon, kappa, tau*(pi/180), xyz_airfoil );
std::vector< double > GetVKTAirfoilCpDist(const double &alpha, const double &epsilon, const double &kappa, const double &tau, std::vector< vec3d > xyz_data)
std::vector< vec3d > GetVKTAirfoilPnts(const int &npts, const double &alpha, const double &epsilon, const double &kappa, const double &tau)
See also
GetVKTAirfoilPnts
Parameters
[in]alphadouble Airfoil angle of attack (Radians)
[in]epsilondouble Airfoil thickness
[in]kappadouble Airfoil camber
[in]taudouble Airfoil trailing edge angle (Radians)
[in]xyz_datavector<vec3d> Vector of points on the airfoil to evaluate
Returns
vector<double> Vector of Cp values for each point in xydata

◆ GetVKTAirfoilPnts()

std::vector<vec3d> vsp::GetVKTAirfoilPnts ( const int &  npts,
const double &  alpha,
const double &  epsilon,
const double &  kappa,
const double &  tau 
)

Get the 2D coordinates an input number of points along a Von K�rm�n-Trefftz airfoil of specified shape

const double pi = 3.14159265358979323846;
const int npts = 122;
const double alpha = 0.0;
const double epsilon = 0.1;
const double kappa = 0.1;
const double tau = 10;
array<vec3d> xyz_airfoil = GetVKTAirfoilPnts(npts, alpha, epsilon, kappa, tau*(pi/180) );
array<double> cp_dist = GetVKTAirfoilCpDist( alpha, epsilon, kappa, tau*(pi/180), xyz_airfoil );
Parameters
[in]nptsNumber of points along the airfoil to return
[in]alphaAirfoil angle of attack (Radians)
[in]epsilonAirfoil thickness
[in]kappaAirfoil camber
[in]tauAirfoil trailing edge angle (Radians)
Returns
Array of points on the VKT airfoil (size = npts)

◆ GetXSecHeight()

double vsp::GetXSecHeight ( const std::string &  xsec_id)

Get the height of an XSec. Note that POINT type XSecs have a width and height of 0, regardless of what width and height it is set to.

// Add Fuselage Geom
string fuseid = AddGeom( "FUSELAGE", "" );
string xsec_surf = GetXSecSurf( fuseid, 0 );
string xsec = GetXSec( xsec_surf, GetNumXSec( xsec_surf ) - 2 ); // Get 2nd to last XSec
SetXSecWidthHeight( xsec, 3.0, 6.0 );
if ( abs( GetXSecHeight( xsec ) - 6.0 ) > 1e-6 ) { Print( "---> Error: API Get/Set Width " ); }
void SetXSecWidthHeight(const std::string &xsec_id, double w, double h)
double GetXSecHeight(const std::string &xsec_id)
int GetNumXSec(const std::string &xsec_surf_id)
See also
SetXSecHeight
Parameters
[in]xsec_idXSec ID
Returns
Xsec height

◆ GetXSecParm()

std::string vsp::GetXSecParm ( const std::string &  xsec_id,
const std::string &  name 
)

Get a specific Parm ID from an Xsec

// Add Fuselage Geom
string fuseid = AddGeom( "FUSELAGE", "" );
string xsec_surf = GetXSecSurf( fuseid, 0 );
ChangeXSecShape( xsec_surf, GetNumXSec( xsec_surf ) - 1, XS_ROUNDED_RECTANGLE );
string xsec = GetXSec( xsec_surf, GetNumXSec( xsec_surf ) - 1 );
string wid = GetXSecParm( xsec, "RoundedRect_Width" );
if ( !ValidParm( wid ) ) { Print( "---> Error: API GetXSecParm " ); }
bool ValidParm(const std::string &id)
std::string GetXSecParm(const std::string &xsec_id, const std::string &name)
Parameters
[in]xsec_idXSec ID
[in]nameParm name
Returns
Parm ID

◆ GetXSecParmIDs()

std::vector<std::string> vsp::GetXSecParmIDs ( const std::string &  xsec_id)

Get all Parm IDs for specified XSec Parm Container

// Add Fuselage Geom
string fuseid = AddGeom( "FUSELAGE", "" );
string xsec_surf = GetXSecSurf( fuseid, 0 );
string xsec = GetXSec( xsec_surf, GetNumXSec( xsec_surf ) - 1 );
array< string > @parm_array = GetXSecParmIDs( xsec );
if ( parm_array.size() < 1 ) { Print( "---> Error: API GetXSecParmIDs " ); }
std::vector< std::string > GetXSecParmIDs(const std::string &xsec_id)
Parameters
[in]xsec_idXSec ID
Returns
Array of Parm IDs

◆ GetXSecShape()

int vsp::GetXSecShape ( const std::string &  xsec_id)

Get the shape of an XSec

// Add Stack
string sid = AddGeom( "STACK", "" );
// Get First (and Only) XSec Surf
string xsec_surf = GetXSecSurf( sid, 0 );
ChangeXSecShape( xsec_surf, 1, XS_EDIT_CURVE );
string xsec = GetXSec( xsec_surf, 1 );
if ( GetXSecShape( xsec ) != XS_EDIT_CURVE ) { Print( "ERROR: GetXSecShape" ); }
@ XS_EDIT_CURVE
Definition: APIDefines.h:1390
int GetXSecShape(const std::string &xsec_id)
See also
XSEC_CRV_TYPE
Parameters
[in]xsec_idXSec ID
Returns
XSec type enum (i.e. XS_ELLIPSE)

◆ GetXSecWidth()

double vsp::GetXSecWidth ( const std::string &  xsec_id)

Get the width of an XSec. Note that POINT type XSecs have a width and height of 0, regardless of what width and height it is set to.

// Add Fuselage Geom
string fuseid = AddGeom( "FUSELAGE", "" );
string xsec_surf = GetXSecSurf( fuseid, 0 );
string xsec = GetXSec( xsec_surf, GetNumXSec( xsec_surf ) - 2 ); // Get 2nd to last XSec
SetXSecWidthHeight( xsec, 3.0, 6.0 );
if ( abs( GetXSecWidth( xsec ) - 3.0 ) > 1e-6 ) { Print( "---> Error: API Get/Set Width " ); }
double GetXSecWidth(const std::string &xsec_id)
See also
SetXSecWidth
Parameters
[in]xsec_idXSec ID
Returns
Xsec width

◆ InsertXSec()

void vsp::InsertXSec ( const std::string &  geom_id,
int  index,
int  type 
)

Insert a cross-section of particular type to the specified geometry after the given index

string wing_id = AddGeom( "WING" );
//===== Add XSec ====//
InsertXSec( wing_id, 1, XS_SIX_SERIES );
@ XS_SIX_SERIES
Definition: APIDefines.h:1387
See also
XSEC_CRV_TYPE
Parameters
[in]geom_idstring Geom ID
[in]indexXSec index
[in]typeXSec type enum (i.e. XS_GENERAL_FUSE)

◆ PasteXSec()

void vsp::PasteXSec ( const std::string &  geom_id,
int  index 
)

Paste the cross-section currently held in memory to the specified geometry

// Add Stack
string sid = AddGeom( "STACK", "" );
// Copy XSec To Clipboard
CopyXSec( sid, 1 );
// Paste To XSec 3
PasteXSec( sid, 3 );
See also
CutXSec, CopyXSec
Parameters
[in]geom_idstring Geom ID
[in]indexXSec index

◆ PromoteCSTLower()

void vsp::PromoteCSTLower ( const std::string &  xsec_id)

Promote the CST for the lower airfoil surface. The XSec must be of type XS_CST_AIRFOIL

See also
GetLowerCSTDegree
Parameters
[in]xsec_idXSec ID

◆ PromoteCSTUpper()

void vsp::PromoteCSTUpper ( const std::string &  xsec_id)

Promote the CST for the upper airfoil surface. The XSec must be of type XS_CST_AIRFOIL

See also
GetUpperCSTDegree
Parameters
[in]xsec_idXSec ID

◆ ReadFileAirfoil()

void vsp::ReadFileAirfoil ( const std::string &  xsec_id,
const std::string &  file_name 
)

Read in XSec shape from airfoil file and set to the specified XSec. The XSec must be of type XS_FILE_AIRFOIL. Airfoil files may be in Lednicer or Selig format with *.af or *.dat extensions.

// Add Fuselage Geom
string fuseid = AddGeom( "FUSELAGE", "" );
string xsec_surf = GetXSecSurf( fuseid, 0 );
string xsec = GetXSec( xsec_surf, 1 );
ReadFileAirfoil( xsec, "airfoil/N0012_VSP.af" );
Parameters
[in]xsec_idXSec ID
[in]file_nameAirfoil XSec file name

◆ ReadFileXSec()

std::vector<vec3d> vsp::ReadFileXSec ( const std::string &  xsec_id,
const std::string &  file_name 
)

Read in XSec shape from fuselage (*.fsx) file and set to the specified XSec. The XSec must be of type XS_FILE_FUSE.

// Add Fuselage Geom
string fuseid = AddGeom( "FUSELAGE", "" );
string xsec_surf = GetXSecSurf( fuseid, 0 );
ChangeXSecShape( xsec_surf, 2, XS_FILE_FUSE );
string xsec = GetXSec( xsec_surf, 2 );
array< vec3d > @vec_array = ReadFileXSec( xsec, "TestXSec.fxs" );
@ XS_FILE_FUSE
Definition: APIDefines.h:1385
std::vector< vec3d > ReadFileXSec(const std::string &xsec_id, const std::string &file_name)
Parameters
[in]xsec_idXSec ID
[in]file_nameFuselage XSec file name
Returns
Array of coordinate points read from the file and set to the XSec

◆ ResetXSecSkinParms()

void vsp::ResetXSecSkinParms ( const std::string &  xsec_id)

Reset all skinning Parms for a specified XSec. Set top, bottom, left, and right strengths, slew, angle, and curvature to 0. Set all symmetry and equality conditions to false.

string fid = AddGeom( "FUSELAGE", "" ); // Add Fuselage
string xsec_surf = GetXSecSurf( fid, 0 ); // Get First (and Only) XSec Surf
int num_xsecs = GetNumXSec( xsec_surf );
string xsec = GetXSec( xsec_surf, 1 );
SetXSecTanAngles( xsec, XSEC_BOTH_SIDES, 0.0 ); // Set Tangent Angles At Cross Section
SetXSecContinuity( xsec, 1 ); // Set Continuity At Cross Section
@ XSEC_BOTH_SIDES
Definition: APIDefines.h:1417
void SetXSecContinuity(const std::string &xsec_id, int cx)
void SetXSecTanAngles(const std::string &xsec_id, int side, double top, double right, double bottom, double left)
void ResetXSecSkinParms(const std::string &xsec_id)
Parameters
[in]xsec_idstring XSec ID

◆ SetAirfoilLowerPnts()

void vsp::SetAirfoilLowerPnts ( const std::string &  xsec_id,
const std::vector< vec3d > &  low_pnt_vec 
)

Set the lower points for an airfoil. The XSec must be of type XS_FILE_AIRFOIL.

// Add Fuselage Geom
string fuseid = AddGeom( "FUSELAGE", "" );
string xsec_surf = GetXSecSurf( fuseid, 0 );
string xsec = GetXSec( xsec_surf, 1 );
ReadFileAirfoil( xsec, "airfoil/N0012_VSP.af" );
array< vec3d > @low_array = GetAirfoilLowerPnts( xsec );
for ( int i = 0 ; i < int( low_array.size() ) ; i++ )
{
low_array[i].scale_y( 0.5 );
}
SetAirfoilUpperPnts( xsec, low_array );
void SetAirfoilUpperPnts(const std::string &xsec_id, const std::vector< vec3d > &up_pnt_vec)
Parameters
[in]xsec_idXSec ID
[in]low_pnt_vecArray of points defining the lower surface of the airfoil

◆ SetAirfoilPnts()

void vsp::SetAirfoilPnts ( const std::string &  xsec_id,
const std::vector< vec3d > &  up_pnt_vec,
const std::vector< vec3d > &  low_pnt_vec 
)

Set the upper and lower points for an airfoil. The XSec must be of type XS_FILE_AIRFOIL.

// Add Fuselage Geom
string fuseid = AddGeom( "FUSELAGE", "" );
string xsec_surf = GetXSecSurf( fuseid, 0 );
string xsec = GetXSec( xsec_surf, 1 );
ReadFileAirfoil( xsec, "airfoil/N0012_VSP.af" );
array< vec3d > @up_array = GetAirfoilUpperPnts( xsec );
array< vec3d > @low_array = GetAirfoilLowerPnts( xsec );
for ( int i = 0 ; i < int( up_array.size() ) ; i++ )
{
up_array[i].scale_y( 2.0 );
low_array[i].scale_y( 0.5 );
}
SetAirfoilPnts( xsec, up_array, low_array );
void SetAirfoilPnts(const std::string &xsec_id, const std::vector< vec3d > &up_pnt_vec, const std::vector< vec3d > &low_pnt_vec)
Parameters
[in]xsec_idXSec ID
[in]up_pnt_vecArray of points defining the upper surface of the airfoil
[in]low_pnt_vecArray of points defining the lower surface of the airfoil

◆ SetAirfoilUpperPnts()

void vsp::SetAirfoilUpperPnts ( const std::string &  xsec_id,
const std::vector< vec3d > &  up_pnt_vec 
)

Set the upper points for an airfoil. The XSec must be of type XS_FILE_AIRFOIL.

// Add Fuselage Geom
string fuseid = AddGeom( "FUSELAGE", "" );
string xsec_surf = GetXSecSurf( fuseid, 0 );
string xsec = GetXSec( xsec_surf, 1 );
ReadFileAirfoil( xsec, "airfoil/N0012_VSP.af" );
array< vec3d > @up_array = GetAirfoilUpperPnts( xsec );
for ( int i = 0 ; i < int( up_array.size() ) ; i++ )
{
up_array[i].scale_y( 2.0 );
}
SetAirfoilUpperPnts( xsec, up_array );
Parameters
[in]xsec_idXSec ID
[in]up_pnt_vecArray of points defining the upper surface of the airfoil

◆ SetLowerCST()

void vsp::SetLowerCST ( const std::string &  xsec_id,
int  deg,
const std::vector< double > &  coefs 
)

Set the CST degree and coefficients for the lower surface of an airfoil. The number of coefficients should be one more than the CST degree. The XSec must be of type XS_CST_AIRFOIL

See also
GetLowerCSTDegree, GetLowerCSTCoefs
Parameters
[in]xsec_idstring XSec ID
[in]degint CST degree of lower airfoil surface
[in]coefsvector<double> Vector of CST coefficients for the lower airfoil surface

◆ SetUpperCST()

void vsp::SetUpperCST ( const std::string &  xsec_id,
int  deg,
const std::vector< double > &  coefs 
)

Set the CST degree and coefficients for the upper surface of an airfoil. The number of coefficients should be one more than the CST degree. The XSec must be of type XS_CST_AIRFOIL

See also
GetUpperCSTDegree, GetUpperCSTCoefs
Parameters
[in]xsec_idstring XSec ID
[in]degint CST degree of upper airfoil surface
[in]coefsvector<double> Vector of CST coefficients for the upper airfoil surface

◆ SetXSecContinuity()

void vsp::SetXSecContinuity ( const std::string &  xsec_id,
int  cx 
)

Set C-type continuity enforcement for a particular XSec

string fid = AddGeom( "FUSELAGE", "" ); // Add Fuselage
string xsec_surf = GetXSecSurf( fid, 0 ); // Get First (and Only) XSec Surf
int num_xsecs = GetNumXSec( xsec_surf );
for ( int i = 0 ; i < num_xsecs ; i++ )
{
string xsec = GetXSec( xsec_surf, i );
SetXSecContinuity( xsec, 1 ); // Set Continuity At Cross Section
}
Parameters
[in]xsec_idstring XSec ID
[in]cxint Continuity level (0, 1, or 2)

◆ SetXSecCurvatures()

void vsp::SetXSecCurvatures ( const std::string &  xsec_id,
int  side,
double  top,
double  right,
double  bottom,
double  left 
)

Set curvatures for the specified XSec

// Add Stack
string sid = AddGeom( "STACK", "" );
// Get First (and Only) XSec Surf
string xsec_surf = GetXSecSurf( sid, 0 );
// Flatten ends
int num_xsecs = GetNumXSec( xsec_surf );
for ( int i = 0 ; i < num_xsecs ; i++ )
{
string xsec = GetXSec( xsec_surf, i );
SetXSecCurvatures( xsec, XSEC_BOTH_SIDES, 0.2 ); // Set Tangent Strengths At Cross Section
}
void SetXSecCurvatures(const std::string &xsec_id, int side, double top, double right, double bottom, double left)
See also
XSEC_SIDES_TYPE
Parameters
[in]xsec_idXSec ID
[in]sideSide type enum (i.e. XSEC_BOTH_SIDES)
[in]topTop curvature
[in]rightRight curvature
[in]bottomBottom curvature
[in]leftLeft curvature

◆ SetXSecHeight()

void vsp::SetXSecHeight ( const std::string &  xsec_id,
double  h 
)

Set the height of an XSec. Note that POINT type XSecs have a width and height of 0, regardless of what is input to SetXSecHeight.

// Add Stack
string sid = AddGeom( "STACK", "" );
// Get First (and Only) XSec Surf
string xsec_surf = GetXSecSurf( sid, 0 );
// Identify XSec 2
string xsec_2 = GetXSec( xsec_surf, 2 );
SetXSecHeight( xsec_2, 1.5 );
void SetXSecHeight(const std::string &xsec_id, double h)
See also
GetXSecHeight
Parameters
[in]xsec_idXSec ID
[in]hXsec height

◆ SetXSecPnts()

void vsp::SetXSecPnts ( const std::string &  xsec_id,
std::vector< vec3d > &  pnt_vec 
)

Set the coordinate points for a specific XSec. The XSec must be of type XS_FILE_FUSE.

// Add Fuselage Geom
string fuseid = AddGeom( "FUSELAGE", "" );
string xsec_surf = GetXSecSurf( fuseid, 0 );
ChangeXSecShape( xsec_surf, 2, XS_FILE_FUSE );
string xsec = GetXSec( xsec_surf, 2 );
array< vec3d > @vec_array = ReadFileXSec( xsec, "TestXSec.fxs" );
if ( vec_array.size() > 0 )
{
vec_array[1] = vec_array[1] * 2.0;
vec_array[3] = vec_array[3] * 2.0;
SetXSecPnts( xsec, vec_array );
}
void SetXSecPnts(const std::string &xsec_id, std::vector< vec3d > &pnt_vec)
Parameters
[in]xsec_idstring XSec ID
[in]pnt_vecvector<vec3d> Vector of XSec coordinate points

◆ SetXSecTanAngles()

void vsp::SetXSecTanAngles ( const std::string &  xsec_id,
int  side,
double  top,
double  right,
double  bottom,
double  left 
)

Set the tangent angles for the specified XSec

// Add Stack
string sid = AddGeom( "STACK", "" );
// Get First (and Only) XSec Surf
string xsec_surf = GetXSecSurf( sid, 0 );
int num_xsecs = GetNumXSec( xsec_surf );
for ( int i = 0 ; i < num_xsecs ; i++ )
{
string xsec = GetXSec( xsec_surf, i );
SetXSecTanAngles( xsec, XSEC_BOTH_SIDES, 10.0 ); // Set Tangent Angles At Cross Section
}
See also
XSEC_SIDES_TYPE
Parameters
[in]xsec_idstring XSec ID
[in]sideint Side type enum (i.e. XSEC_BOTH_SIDES)
[in]topdouble Top angle (degrees)
[in]rightdouble Right angle (degrees)
[in]bottomdouble Bottom angle (degrees)
[in]leftdouble Left angle (degrees)

◆ SetXSecTanSlews()

void vsp::SetXSecTanSlews ( const std::string &  xsec_id,
int  side,
double  top,
double  right,
double  bottom,
double  left 
)

Set the tangent slew angles for the specified XSec

// Add Stack
string sid = AddGeom( "STACK", "" );
// Get First (and Only) XSec Surf
string xsec_surf = GetXSecSurf( sid, 0 );
int num_xsecs = GetNumXSec( xsec_surf );
for ( int i = 0 ; i < num_xsecs ; i++ )
{
string xsec = GetXSec( xsec_surf, i );
SetXSecTanSlews( xsec, XSEC_BOTH_SIDES, 5.0 ); // Set Tangent Slews At Cross Section
}
void SetXSecTanSlews(const std::string &xsec_id, int side, double top, double right, double bottom, double left)
See also
XSEC_SIDES_TYPE
Parameters
[in]xsec_idXSec ID
[in]sideSide type enum (i.e. XSEC_BOTH_SIDES)
[in]topTop angle (degrees)
[in]rightRight angle (degrees)
[in]bottomBottom angle (degrees)
[in]leftLeft angle (degrees)

◆ SetXSecTanStrengths()

void vsp::SetXSecTanStrengths ( const std::string &  xsec_id,
int  side,
double  top,
double  right,
double  bottom,
double  left 
)

Set the tangent strengths for the specified XSec

// Add Stack
string sid = AddGeom( "STACK", "" );
// Get First (and Only) XSec Surf
string xsec_surf = GetXSecSurf( sid, 0 );
// Flatten ends
int num_xsecs = GetNumXSec( xsec_surf );
for ( int i = 0 ; i < num_xsecs ; i++ )
{
string xsec = GetXSec( xsec_surf, i );
SetXSecTanStrengths( xsec, XSEC_BOTH_SIDES, 0.8 ); // Set Tangent Strengths At Cross Section
}
void SetXSecTanStrengths(const std::string &xsec_id, int side, double top, double right, double bottom, double left)
See also
XSEC_SIDES_TYPE
Parameters
[in]xsec_idXSec ID
[in]sideSide type enum (i.e. XSEC_BOTH_SIDES)
[in]topTop strength
[in]rightRight strength
[in]bottomBottom strength
[in]leftLeft strength

◆ SetXSecWidth()

void vsp::SetXSecWidth ( const std::string &  xsec_id,
double  w 
)

Set the width of an XSec. Note that POINT type XSecs have a width and height of 0, regardless of what is input to SetXSecWidth.

// Add Stack
string sid = AddGeom( "STACK", "" );
// Get First (and Only) XSec Surf
string xsec_surf = GetXSecSurf( sid, 0 );
// Identify XSec 2
string xsec_2 = GetXSec( xsec_surf, 2 );
SetXSecWidth( xsec_2, 1.5 );
void SetXSecWidth(const std::string &xsec_id, double w)
See also
GetXSecWidth
Parameters
[in]xsec_idXSec ID
[in]wXsec width

◆ SetXSecWidthHeight()

void vsp::SetXSecWidthHeight ( const std::string &  xsec_id,
double  w,
double  h 
)

Set the width and height of an XSec. Note, if the XSec is an EDIT_CURVE type and PreserveARFlag is true, the input width value will be ignored and instead set from on the input height and aspect ratio. Use SetXSecWidth and SetXSecHeight directly to avoid this.

// Add Stack
string sid = AddGeom( "STACK", "" );
// Get First (and Only) XSec Surf
string xsec_surf = GetXSecSurf( sid, 0 );
// Identify XSec 2
string xsec_2 = GetXSec( xsec_surf, 2 );
SetXSecWidthHeight( xsec_2, 1.5, 1.5 );
See also
SetXSecWidth, SetXSecHeight
Parameters
[in]xsec_idXSec ID
[in]wXsec width
[in]hXsec height

◆ WriteBezierAirfoil()

void vsp::WriteBezierAirfoil ( const std::string &  file_name,
const std::string &  geom_id,
const double &  foilsurf_u 
)

Write out the untwisted unit-length 2D Bezier curve for the specified airfoil in custom *.bz format. The output will describe the analytical shape of the airfoil. See BezierAirfoilExample.m and BezierCtrlToCoordPnts.m for examples of discretizing the Bezier curve and generating a Selig airfoil file.

//==== Add Wing Geometry and Set Parms ====//
string wing_id = AddGeom( "WING", "" );
const double u = 0.5; // export airfoil at mid span location
//==== Write Bezier Airfoil File ====//
WriteBezierAirfoil( "Example_Bezier.bz", wing_id, u );
void WriteBezierAirfoil(const std::string &file_name, const std::string &geom_id, const double &foilsurf_u)
Parameters
[in]file_nameAirfoil (*.bz) output file name
[in]geom_idstring Geom ID
[in]foilsurf_uU location (range: 0 - 1) along the surface. The foil surface does not include root and tip caps (i.e. 2 section wing -> XSec0 @ u=0, XSec1 @ u=0.5, XSec2 @ u=1.0)

◆ WriteSeligAirfoil()

void vsp::WriteSeligAirfoil ( const std::string &  file_name,
const std::string &  geom_id,
const double &  foilsurf_u 
)

Write out the untwisted unit-length 2D coordinate points for the specified airfoil in Selig format. Coordinate points follow the on-screen wire frame W tessellation.

//==== Add Wing Geometry and Set Parms ====//
string wing_id = AddGeom( "WING", "" );
const double u = 0.5; // export airfoil at mid span location
//==== Write Selig Airfoil File ====//
WriteSeligAirfoil( "Example_Selig.dat", wing_id, u );
void WriteSeligAirfoil(const std::string &file_name, const std::string &geom_id, const double &foilsurf_u)
See also
GetAirfoilCoordinates
Parameters
[in]file_nameAirfoil (*.dat) output file name
[in]geom_idstring Geom ID
[in]foilsurf_uU location (range: 0 - 1) along the surface. The foil surface does not include root and tip caps (i.e. 2 section wing -> XSec0 @ u=0, XSec1 @ u=0.5, XSec2 @ u=1.0)