OpenVSPAPI  3.20.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 incuded in the Specialized Geometry group.

Click here to return to the main page. More...

Functions

void CutXSec (const string &in geom_id, int index)
 
void CopyXSec (const string &in geom_id, int index)
 
void PasteXSec (const string &in geom_id, int index)
 
void InsertXSec (const string &in geom_id, int index, int type)
 
int GetXSecShape (const string &in xsec_id)
 
double GetXSecWidth (const string &in xsec_id)
 
double GetXSecHeight (const string &in xsec_id)
 
void SetXSecWidthHeight (const string &in xsec_id, double w, double h)
 
string[] GetXSecParmIDs (const string &in xsec_id)
 
string GetXSecParm (const string &in xsec_id, const string &in name)
 
vec3d[] ReadFileXSec (const string &in xsec_id, const string &in file_name)
 
void SetXSecPnts (const string &in xsec_id, vec3d[]@ pnt_arr)
 
vec3d ComputeXSecPnt (const string &in xsec_id, double fract)
 
vec3d ComputeXSecTan (const string &in xsec_id, double fract)
 
void ResetXSecSkinParms (const string &in xsec_id)
 
void SetXSecContinuity (const string &in xsec_id, int cx)
 
void SetXSecTanAngles (const string &in xsec_id, int side, double top, double right=- 1.0e12, double bottom=- 1.0e12, double left=- 1.0e12)
 
void SetXSecTanSlews (const string &in xsec_id, int side, double top, double right=- 1.0e12, double bottom=- 1.0e12, double left=- 1.0e12)
 
void SetXSecTanStrengths (const string &in xsec_id, int side, double top, double right=- 1.0e12, double bottom=- 1.0e12, double left=- 1.0e12)
 
void SetXSecCurvatures (const string &in xsec_id, int side, double top, double right=- 1.0e12, double bottom=- 1.0e12, double left=- 1.0e12)
 
void ReadFileAirfoil (const string &in xsec_id, const string &in file_name)
 
void SetAirfoilPnts (const string &in xsec_id, vec3d[]@ up_pnt_vec, vec3d[]@ low_pnt_vec)
 
vec3d[] GetHersheyBarLiftDist (const int &in npts, const double &in alpha, const double &in Vinf, const double &in span, bool full_span_flag=false)
 
vec3d[] GetHersheyBarDragDist (const int &in npts, const double &in alpha, const double &in Vinf, const double &in span, bool full_span_flag=false)
 
vec3d[] GetVKTAirfoilPnts (const int &in npts, const double &in alpha, const double &in epsilon, const double &in kappa, const double &in tau)
 
double[] GetVKTAirfoilCpDist (const double &in alpha, const double &in epsilon, const double &in kappa, const double &in tau, vec3d[]@ xydata)
 
vec3d[] GetEllipsoidSurfPnts (const vec3d &in center, const vec3d &in abc_rad, int u_npts=20, int w_npts=20)
 
vec3d[] GetFeatureLinePnts (const string &in geom_id)
 
double[] GetEllipsoidCpDist (vec3d[]@ surf_pnt_arr, const vec3d &in abc_rad, const vec3d &in V_inf)
 
vec3d[] GetAirfoilUpperPnts (const string &in xsec_id)
 
vec3d[] GetAirfoilLowerPnts (const string &in xsec_id)
 
double[] GetUpperCSTCoefs (const string &in xsec_id)
 
double[] GetLowerCSTCoefs (const string &in xsec_id)
 
int GetUpperCSTDegree (const string &in xsec_id)
 
int GetLowerCSTDegree (const string &in xsec_id)
 
void SetUpperCST (const string &in xsec_id, int deg, double[]@ coeff_arr)
 
void SetLowerCST (const string &in xsec_id, int deg, double[]@ coeff_arr)
 
void PromoteCSTUpper (const string &in xsec_id)
 
void PromoteCSTLower (const string &in xsec_id)
 
void DemoteCSTUpper (const string &in xsec_id)
 
void DemoteCSTLower (const string &in xsec_id)
 
void FitAfCST (const string &in xsec_surf_id, int xsec_index, int deg)
 
void WriteBezierAirfoil (const string &in file_name, const string &in geom_id, const double &in foilsurf_u)
 
void WriteSeligAirfoil (const string &in file_name, const string &in geom_id, const double &in foilsurf_u)
 
vec3d[] GetAirfoilCoordinates (const string &in geom_id, const double &in foilsurf_u)
 

Detailed Description

Function Documentation

◆ ComputeXSecPnt()

vec3d ComputeXSecPnt ( const string &in  xsec_id,
double  fract 
)

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

//==== Add Spine Surf ====//
string spine_surf = AddXSecSurf();
string spine_xsec = AppendXSec( spine_surf, XS_GENERAL_FUSE );
//==== Get The XSec Surf ====//
string xsec_surf = GetXSecSurf( geom_id, 0 );
//==== Define XSecs ====//
int NUM_XSECS = 11;
for ( int i = 0 ; i < NUM_XSECS ; i++ )
{
string xsec = GetXSec( xsec_surf, i );
vec3d pnt = ComputeXSecPnt( spine_xsec, i*0.10 );
SetCustomXSecLoc( xsec, pnt );
}
Parameters
[in]xsec_idXSec ID
[in]fractCurve parameter value (range: 0 - 1)
Returns
3D coordinate point

◆ ComputeXSecTan()

vec3d ComputeXSecTan ( const string &in  xsec_id,
double  fract 
)

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

//==== Add Spine Surf ====//
string spine_surf = AddXSecSurf();
string spine_xsec = AppendXSec( spine_surf, XS_GENERAL_FUSE );
//==== Get The XSec Surf ====//
string xsec_surf = GetXSecSurf( geom_id, 0 );
//==== Define XSecs ====//
int NUM_XSECS = 11;
for ( int i = 0 ; i < NUM_XSECS ; i++ )
{
string xsec = GetXSec( xsec_surf, i );
vec3d tan = ComputeXSecTan( spine_xsec, i*0.10 );
double ang = signed_angle( tan, vec3d( 1, 0, 0 ), vec3d( 0, -1, 0 ) );
SetCustomXSecRot( xsec, vec3d( 0, Rad2Deg(ang), 0) );
}
Parameters
[in]xsec_idXSec ID
[in]fractCurve parameter value (range: 0 - 1)
Returns
Tangent vector

◆ CopyXSec()

void CopyXSec ( const string &in  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 );
See also
PasteXSec
Parameters
[in]geom_idGeom ID
[in]indexXSec index

◆ CutXSec()

void CutXSec ( const string &in  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
See also
PasteXSec
Parameters
[in]geom_idGeom ID
[in]indexXSec index

◆ DemoteCSTLower()

void DemoteCSTLower ( const string &in  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 DemoteCSTUpper ( const string &in  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 FitAfCST ( const string &in  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()

vec3d [] GetAirfoilCoordinates ( const string &in  geom_id,
const double &in  foilsurf_u 
)

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

See also
WriteSeligAirfoil
Parameters
[in]geom_idGeom 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()

vec3d [] GetAirfoilLowerPnts ( const string &in  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 );
xsec = GetXSec( xsec_surf, 1 );
ReadFileAirfoil( xsec, "airfoil/N0012_VSP.af" );
array< vec3d > @low_array = GetAirfoilLowerPnts( xsec );
See also
SetAirfoilPnts
Parameters
[in]xsec_idXSec ID
Returns
Array of coordinate points for the lower airfoil surface

◆ GetAirfoilUpperPnts()

vec3d [] GetAirfoilUpperPnts ( const string &in  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 );
xsec = GetXSec( xsec_surf, 1 );
ReadFileAirfoil( xsec, "airfoil/N0012_VSP.af" );
array< vec3d > @up_array = GetAirfoilUpperPnts( xsec );
See also
SetAirfoilPnts
Parameters
[in]xsec_idXSec ID
Returns
Array of coordinate points for the upper airfoil surface

◆ GetEllipsoidCpDist()

double [] GetEllipsoidCpDist ( vec3d@[]  surf_pnt_arr,
const vec3d &in  abc_rad,
const vec3d &in  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 );
See also
GetEllipsoidSurfPnts
Parameters
[in]surf_pnt_arrArray of points on the ellipsoid surface to assess
[in]abc_radRadius along the A (X), B (Y), and C (Z) axes
[in]V_inf3D components of freestream velocity
Returns
Array of Cp results corresponding to each point in surf_pnt_arr

◆ GetEllipsoidSurfPnts()

vec3d [] GetEllipsoidSurfPnts ( const vec3d &in  center,
const vec3d &in  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()

vec3d [] GetFeatureLinePnts ( const string &in  geom_id)

Get the points along the feature lines of a particular Geom

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

◆ GetHersheyBarDragDist()

vec3d [] GetHersheyBarDragDist ( const int &in  npts,
const double &in  alpha,
const double &in  Vinf,
const double &in  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 );
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()

vec3d [] GetHersheyBarLiftDist ( const int &in  npts,
const double &in  alpha,
const double &in  Vinf,
const double &in  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()

double [] GetLowerCSTCoefs ( const string &in  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_idXSec ID
Returns
Array of CST coefficients for the lower airfoil surface

◆ GetLowerCSTDegree()

int GetLowerCSTDegree ( const string &in  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
CST Degree for lower airfoil surface

◆ GetUpperCSTCoefs()

double [] GetUpperCSTCoefs ( const string &in  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_idXSec ID
Returns
Array of CST coefficients for the upper airfoil surface

◆ GetUpperCSTDegree()

int GetUpperCSTDegree ( const string &in  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_idXSec ID
Returns
CST Degree for upper airfoil surface

◆ GetVKTAirfoilCpDist()

double [] GetVKTAirfoilCpDist ( const double &in  alpha,
const double &in  epsilon,
const double &in  kappa,
const double &in  tau,
vec3d@[]  xydata 
)

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 );
See also
GetVKTAirfoilPnts
Parameters
[in]alphaAirfoil angle of attack (Radians)
[in]epsilonAirfoil thickness
[in]kappaAirfoil camber
[in]tauAirfoil trailing edge angle (Radians)
[in]xydataArray of points on the airfoil to evaluate
Returns
Array of Cp values for each point in xydata

◆ GetVKTAirfoilPnts()

vec3d [] GetVKTAirfoilPnts ( const int &in  npts,
const double &in  alpha,
const double &in  epsilon,
const double &in  kappa,
const double &in  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 GetXSecHeight ( const string &in  xsec_id)

Get the height of an XSec

// Add Fuselage Geom
string fuseid = AddGeom( "FUSELAGE", "" );
string xsec_surf = GetXSecSurf( fuseid, 0 );
xsec = GetXSec( xsec_surf, GetNumXSec( xsec_surf ) - 1 );
SetXSecWidthHeight( xsec, 3.0, 6.0 );
if ( abs( GetXSecHeight( xsec ) - 6.0 ) > 1e-6 ) { Print( "---> Error: API Get/Set Width " ); }
Parameters
[in]xsec_idXSec ID
Returns
Xsec height

◆ GetXSecParm()

string GetXSecParm ( const string &in  xsec_id,
const string &in  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 );
xsec = GetXSec( xsec_surf, GetNumXSec( xsec_surf ) - 1 );
string wid = GetXSecParm( xsec, "RoundedRect_Width" );
if ( !ValidParm( wid ) ) { Print( "---> Error: API GetXSecParm " ); }
Parameters
[in]xsec_idXSec ID
[in]nameParm name
Returns
Parm ID

◆ GetXSecParmIDs()

string [] GetXSecParmIDs ( const string &in  xsec_id)

Get all Parm IDs for specified XSec Parm Container

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

◆ GetXSecShape()

int GetXSecShape ( const string &in  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" ); }
See also
XSEC_CRV_TYPE
Parameters
[in]xsec_idXSec ID
Returns
XSec type enum (i.e. XS_ELLIPSE)

◆ GetXSecWidth()

double GetXSecWidth ( const string &in  xsec_id)

Get the width of an XSec

// Add Fuselage Geom
string fuseid = AddGeom( "FUSELAGE", "" );
string xsec_surf = GetXSecSurf( fuseid, 0 );
string xsec = GetXSec( xsec_surf, GetNumXSec( xsec_surf ) - 1 );
SetXSecWidthHeight( xsec, 3.0, 6.0 );
if ( abs( GetXSecWidth( xsec ) - 3.0 ) > 1e-6 ) { Print( "---> Error: API Get/Set Width " ); }
Parameters
[in]xsec_idXSec ID
Returns
Xsec width

◆ InsertXSec()

void InsertXSec ( const string &in  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 );
See also
XSEC_CRV_TYPE
Parameters
[in]geom_idGeom ID
[in]indexXSec index
[in]typeXSec type enum (i.e. XS_GENERAL_FUSE)

◆ PasteXSec()

void PasteXSec ( const string &in  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_idGeom ID
[in]indexXSec index

◆ PromoteCSTLower()

void PromoteCSTLower ( const string &in  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 PromoteCSTUpper ( const string &in  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 ReadFileAirfoil ( const string &in  xsec_id,
const string &in  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 );
xsec = GetXSec( xsec_surf, 1 );
ReadFileAirfoil( xsec, "airfoil/N0012_VSP.af" );
Parameters
[in]xsec_idXSec ID
[in]file_nameAirfoil XSec file name

◆ ReadFileXSec()

vec3d [] ReadFileXSec ( const string &in  xsec_id,
const string &in  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 );
xsec = GetXSec( xsec_surf, 2 );
array< vec3d > @vec_array = ReadFileXSec( xsec, "TestXSec.fxs" );
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 ResetXSecSkinParms ( const string &in  xsec_id)

Reset all skining Parms for a specified XSec. Set top, bottom, left, and right strenghts, 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
Parameters
[in]xsec_idXSec ID

◆ SetAirfoilPnts()

void SetAirfoilPnts ( const string &in  xsec_id,
vec3d@[]  up_pnt_vec,
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 );
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 );
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

◆ SetLowerCST()

void SetLowerCST ( const string &in  xsec_id,
int  deg,
double@[]  coeff_arr 
)

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_idXSec ID
[in]degCST degree of lower airfoil surface
[in]coeff_arrArray of CST coefficients for the lower airfoil surface

◆ SetUpperCST()

void SetUpperCST ( const string &in  xsec_id,
int  deg,
double@[]  coeff_arr 
)

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_idXSec ID
[in]degCST degree of upper airfoil surface
[in]coeff_arrArray of CST coefficients for the upper airfoil surface

◆ SetXSecContinuity()

void SetXSecContinuity ( const string &in  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_idXSec ID
[in]cxContinuity level (0, 1, or 2)

◆ SetXSecCurvatures()

void SetXSecCurvatures ( const string &in  xsec_id,
int  side,
double  top,
double  right = - 1.0e12,
double  bottom = - 1.0e12,
double  left = - 1.0e12 
)

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
}
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

◆ SetXSecPnts()

void SetXSecPnts ( const string &in  xsec_id,
vec3d@[]  pnt_arr 
)

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 );
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 );
}
Parameters
[in]xsec_idXSec ID
[in]pnt_arrArray of XSec coordinate points

◆ SetXSecTanAngles()

void SetXSecTanAngles ( const string &in  xsec_id,
int  side,
double  top,
double  right = - 1.0e12,
double  bottom = - 1.0e12,
double  left = - 1.0e12 
)

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_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)

◆ SetXSecTanSlews()

void SetXSecTanSlews ( const string &in  xsec_id,
int  side,
double  top,
double  right = - 1.0e12,
double  bottom = - 1.0e12,
double  left = - 1.0e12 
)

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
}
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 SetXSecTanStrengths ( const string &in  xsec_id,
int  side,
double  top,
double  right = - 1.0e12,
double  bottom = - 1.0e12,
double  left = - 1.0e12 
)

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
}
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

◆ SetXSecWidthHeight()

void SetXSecWidthHeight ( const string &in  xsec_id,
double  w,
double  h 
)

Set the width and height of an XSec

// 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 );
Parameters
[in]xsec_idXSec ID
[in]wXsec width
[in]hXsec height

◆ WriteBezierAirfoil()

void WriteBezierAirfoil ( const string &in  file_name,
const string &in  geom_id,
const double &in  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 );
Parameters
[in]file_nameAirfoil (*.bz) output file name
[in]geom_idGeom 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 WriteSeligAirfoil ( const string &in  file_name,
const string &in  geom_id,
const double &in  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 tesselation.

//==== 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 );
See also
GetAirfoilCoordinates
Parameters
[in]file_nameAirfoil (*.dat) output file name
[in]geom_idGeom 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)
GetEllipsoidCpDist
double[] GetEllipsoidCpDist(vec3d[]@ surf_pnt_arr, const vec3d &in abc_rad, const vec3d &in V_inf)
ValidParm
bool ValidParm(const string &in id)
XS_FILE_AIRFOIL
Definition: openvsp_as.h:2124
SetXSecTanStrengths
void SetXSecTanStrengths(const string &in xsec_id, int side, double top, double right=- 1.0e12, double bottom=- 1.0e12, double left=- 1.0e12)
GetXSecParm
string GetXSecParm(const string &in xsec_id, const string &in name)
XS_ROUNDED_RECTANGLE
Definition: openvsp_as.h:2116
Deg2Rad
double Deg2Rad(double d)
Rad2Deg
double Rad2Deg(double r)
WriteBezierAirfoil
void WriteBezierAirfoil(const string &in file_name, const string &in geom_id, const double &in foilsurf_u)
XS_GENERAL_FUSE
Definition: openvsp_as.h:2117
AppendXSec
string AppendXSec(const string &in xsec_surf_id, int type)
SetXSecWidthHeight
void SetXSecWidthHeight(const string &in xsec_id, double w, double h)
XS_EDIT_CURVE
Definition: openvsp_as.h:2123
WriteSeligAirfoil
void WriteSeligAirfoil(const string &in file_name, const string &in geom_id, const double &in foilsurf_u)
XSEC_BOTH_SIDES
Definition: openvsp_as.h:2140
GetXSecShape
int GetXSecShape(const string &in xsec_id)
PasteXSec
void PasteXSec(const string &in geom_id, int index)
InsertXSec
void InsertXSec(const string &in geom_id, int index, int type)
signed_angle
double signed_angle(const vec3d &in a, const vec3d &in b, const vec3d &in ref)
AddXSecSurf
string AddXSecSurf()
GetXSecParmIDs
string[] GetXSecParmIDs(const string &in xsec_id)
GetXSec
string GetXSec(const string &in xsec_surf_id, int xsec_index)
SetCustomXSecLoc
void SetCustomXSecLoc(const string &in xsec_id, const vec3d &in loc)
SetCustomXSecRot
void SetCustomXSecRot(const string &in xsec_id, const vec3d &in rot)
GetAirfoilLowerPnts
vec3d[] GetAirfoilLowerPnts(const string &in xsec_id)
SetAirfoilPnts
void SetAirfoilPnts(const string &in xsec_id, vec3d[]@ up_pnt_vec, vec3d[]@ low_pnt_vec)
vec3d
A class for representing 3D vectors.
Definition: openvsp_as.h:341
GetXSecWidth
double GetXSecWidth(const string &in xsec_id)
GetHersheyBarDragDist
vec3d[] GetHersheyBarDragDist(const int &in npts, const double &in alpha, const double &in Vinf, const double &in span, bool full_span_flag=false)
ChangeXSecShape
void ChangeXSecShape(const string &in xsec_surf_id, int xsec_index, int type)
array
AngelScript ScriptExtension for representing the C++ std::vector.
Definition: openvsp_as.h:244
ResetXSecSkinParms
void ResetXSecSkinParms(const string &in xsec_id)
GetXSecSurf
string GetXSecSurf(const string &in geom_id, int index)
AddGeom
string AddGeom(const string &in type, const string &in parent=string())
SetXSecContinuity
void SetXSecContinuity(const string &in xsec_id, int cx)
ComputeXSecTan
vec3d ComputeXSecTan(const string &in xsec_id, double fract)
XS_SIX_SERIES
Definition: openvsp_as.h:2120
SetXSecTanAngles
void SetXSecTanAngles(const string &in xsec_id, int side, double top, double right=- 1.0e12, double bottom=- 1.0e12, double left=- 1.0e12)
GetAirfoilUpperPnts
vec3d[] GetAirfoilUpperPnts(const string &in xsec_id)
Print
void Print(const string &in data, bool new_line=true)
GetVKTAirfoilCpDist
double[] GetVKTAirfoilCpDist(const double &in alpha, const double &in epsilon, const double &in kappa, const double &in tau, vec3d[]@ xydata)
ComputeXSecPnt
vec3d ComputeXSecPnt(const string &in xsec_id, double fract)
SetXSecPnts
void SetXSecPnts(const string &in xsec_id, vec3d[]@ pnt_arr)
SetXSecTanSlews
void SetXSecTanSlews(const string &in xsec_id, int side, double top, double right=- 1.0e12, double bottom=- 1.0e12, double left=- 1.0e12)
GetNumXSec
int GetNumXSec(const string &in xsec_surf_id)
XS_FILE_FUSE
Definition: openvsp_as.h:2118
GetXSecHeight
double GetXSecHeight(const string &in xsec_id)
ReadFileXSec
vec3d[] ReadFileXSec(const string &in xsec_id, const string &in file_name)
SetXSecCurvatures
void SetXSecCurvatures(const string &in xsec_id, int side, double top, double right=- 1.0e12, double bottom=- 1.0e12, double left=- 1.0e12)
CopyXSec
void CopyXSec(const string &in geom_id, int index)
CutXSec
void CutXSec(const string &in geom_id, int index)
GetVKTAirfoilPnts
vec3d[] GetVKTAirfoilPnts(const int &in npts, const double &in alpha, const double &in epsilon, const double &in kappa, const double &in tau)
GetHersheyBarLiftDist
vec3d[] GetHersheyBarLiftDist(const int &in npts, const double &in alpha, const double &in Vinf, const double &in span, bool full_span_flag=false)
ReadFileAirfoil
void ReadFileAirfoil(const string &in xsec_id, const string &in file_name)