#include <Point2D.h>
Public Member Functions | |
| void | average (const Point2D &p) |
| double | azimuthTo (const Point2D &p) const |
| short | compare (const Point2D &p) const |
| void | degreesToDM () |
| void | degreesToDMS () |
| double | distanceTo (const Point2D &p) const |
| double | distanceToSegment (const Point2D &p1, const Point2D &p2) const |
| void | DMSToDegrees () |
| void | DMToDegrees () |
| bool | fromString (const StringSection &str) |
| double | geographicalAzimuthTo (const Point2D &p) const |
| double | geographicalDistanceTo (const Point2D &p) const |
| void | geographicalToRectangular (const Point2D &ref, double rotation=0.0) |
| void | geographicalToUtm (const QPoint *zone=0) |
| void | interpole (double valX, const Point2D &p1, const Point2D &p2) |
| bool | isHit (const Point2D &p, double tolX, double tolY) const |
| bool | isSimilar (const Point2D &p, double tolerance) const |
| bool | isValid () const |
| void | move (double distance, const Angle &azimuth) |
| bool | operator!= (const Point2D &obj) const |
| Point2D | operator* (double mul) const |
| Point2D | operator* (const Point2D &p) const |
| void | operator*= (const Point2D &p) |
| void | operator*= (double mul) |
| void | operator*= (const Matrix2x2 &transformation) |
| Point2D | operator+ (const Point2D &p) const |
| void | operator+= (const Point2D &p) |
| Point2D | operator- (const Point2D &p) const |
| void | operator-= (const Point2D &p) |
| Point2D | operator/ (double div) const |
| Point2D | operator/ (const Point2D &p) const |
| void | operator/= (const Point2D &p) |
| void | operator/= (double div) |
| bool | operator< (const Point2D &p) const |
| bool | operator<= (const Point2D &obj) const |
| Point2D & | operator= (const Point2D &p) |
| bool | operator== (const Point2D &p) const |
| bool | operator> (const Point2D &obj) const |
| bool | operator>= (const Point2D &obj) const |
| Point2D () | |
| Point2D (double x, double y) | |
| Point2D (const QPoint &p) | |
| Point2D (const QPointF &p) | |
| Point2D (const Point2D &p) | |
| Point2D (const QVector< double > &p) | |
| void | rectangularToGeographical (const Point2D &ref, double rotation=0.0) |
| void | rotate (const Angle &a) |
| double | scalarProduct (const Point2D &p) const |
| void | set (double xi, double yi) |
| void | setValid (bool) |
| void | setX (double v) |
| void | setY (double v) |
| void | setY (double v, const CurvePointOptions *) |
| QString | toString (int precision=6, char format='g') const |
| void | translate (const Point2D &p) |
| void | utmToGeographical (const QPoint &zone) |
| QPoint | utmZoneIndex () const |
| double | x () const |
| double | y () const |
| double | y (const CurvePointOptions *) const |
Static Public Member Functions | |
| static double | utmLatitudeReference (QPoint z) |
| static double | utmLongitudeReference (QPoint z) |
| static QString | utmZone (QPoint z, bool *ok=0) |
| static QPoint | utmZone (const QString &z, bool *ok=0) |
Basic properties of a 2D point (coordinates in double precision)
| QGpCoreTools::Point2D::Point2D | ( | ) | [inline] |
Default constructor
Referenced by distanceToSegment(), operator*(), operator+(), operator-(), and operator/().
{
_x=0.;
_y=0.;
}
| QGpCoreTools::Point2D::Point2D | ( | double | x, |
| double | y | ||
| ) | [inline] |
Constructor setting x and y
{
_x=xi;
_y=yi;
}
| QGpCoreTools::Point2D::Point2D | ( | const QPoint & | p | ) | [inline] |
{
_x=p.x();
_y=p.y();
}
| QGpCoreTools::Point2D::Point2D | ( | const QPointF & | p | ) | [inline] |
{
_x=p.x();
_y=p.y();
}
| QGpCoreTools::Point2D::Point2D | ( | const Point2D & | p | ) | [inline] |
| QGpCoreTools::Point2D::Point2D | ( | const QVector< double > & | p | ) |
Copy constructor from a vector
{
ASSERT(p.count()>1);
_x=p.at(0);
_y=p.at(1);
}
| void QGpCoreTools::Point2D::average | ( | const Point2D & | p | ) | [inline] |
{_y= 0.5*(_y+p._y);}
| double QGpCoreTools::Point2D::azimuthTo | ( | const Point2D & | p | ) | const |
Returns the azimuth from this point to p in radians in the mathematical sense (from 0 to 2 PI, 0 is East).
References TRACE, x(), and y().
Referenced by SourceParameters::azimuth(), SourceParameters::azimuthMath(), FKTimeWindows::currentVelocitySlowness(), FKLoopTask::exportResults(), ToolRefra::initStations(), main(), SciFigs::PolarGridPlot::mouseReleaseEvent(), SciFigs::SlopeEstimator::paintText(), CoordReader::parse(), SignalReader::parse(), SciFigs::CircleViewer::Limits::polarLimits(), SourceParameters::setDistance(), FKArrayMap::setSourceSignals(), and GeopsyGui::RelativePositions::updateCurrentDistAz().
{
TRACE;
double a=atan2(p.y()-_y, p.x()-_x);
if(a<0.0) a+=2*M_PI;
return a;
}
| short QGpCoreTools::Point2D::compare | ( | const Point2D & | p | ) | const [inline] |
{
if(_x<p._x) return(-1);
else if(_x>p._x) return(1);
else if(_y<p._y) return(-1);
else if(_y>p._y) return(1);
else return(0);
}
| void QGpCoreTools::Point2D::degreesToDM | ( | ) |
Converts point from geographical coordinates in degrees to D.MMmm... (degrees, minutes).
Referenced by CoordReader::parse().
{
_x=Angle::degreesToDM(_x);
_y=Angle::degreesToDM(_y);
}
| void QGpCoreTools::Point2D::degreesToDMS | ( | ) |
Converts point from geographical coordinates in degrees to D.MMSSsss... (degrees, minutes, seconds).
Referenced by CoordReader::parse().
{
_x=Angle::degreesToDMS(_x);
_y=Angle::degreesToDMS(_y);
}
| double QGpCoreTools::Point2D::distanceTo | ( | const Point2D & | p | ) | const |
Calculates distance to p.
References QGpCoreTools::sqrt(), TRACE, x(), and y().
Referenced by QGpCoreTools::IrregularGrid2DData::crossSection(), FKTimeWindows::currentVelocitySlowness(), SourceParameters::distance(), QGpCoreTools::Circle::distanceTo(), QGpCoreTools::Segment2D::distanceTo(), distanceToSegment(), FKLoopTask::exportResults(), SciFigs::PolarGridPlot::mouseReleaseEvent(), SourceParameters::setAzimuth(), TapePositioningSystem::Triangulator::setPriorCoordinates(), and FKArrayMap::setSourceSignals().
{
TRACE;
double tmp, dist;
dist=0.;
tmp=_x-p.x();
dist+=tmp*tmp;
tmp=_y-p.y();
dist+=tmp*tmp;
return ::sqrt(dist);
}
| double QGpCoreTools::Point2D::distanceToSegment | ( | const Point2D & | p1, |
| const Point2D & | p2 | ||
| ) | const |
Calculate the distance to a segment defined by p1 and p2.
References distanceTo(), Point2D(), TRACE, x(), and y().
{
TRACE;
// Find the intersection between segment and its perpendicular passing at this
double dys=p2.y()-p1.y();
if(dys==0) {
if(p1.x()<p2.x()) {
if(_x>p2.x()) return distanceTo(p2);
if(_x<p1.x()) return distanceTo(p1);
} else {
if(_x>p1.x()) return distanceTo(p1);
if(_x<p2.x()) return distanceTo(p2);
}
return fabs(_y-p1.y());
}
double dxs=p2.x()-p1.x();
if(dxs==0) {
if(p1.y()<p2.y()) {
if(_y>p2.y()) return distanceTo(p2);
if(_y<p1.y()) return distanceTo(p1);
} else {
if(_y>p1.y()) return distanceTo(p1);
if(_y<p2.y()) return distanceTo(p2);
}
return fabs(_x-p1.x());
}
double coeff=dys/dxs;
double coeffp=dxs/dys;
double interx=(coeff*p1.x()-coeffp*_x+_y-p1.y())/(coeff-coeffp);
if(p1.x()<p2.x()) {
if(interx>p2.x()) return distanceTo(p2);
if(interx<p1.x()) return distanceTo(p1);
} else {
if(interx>p1.x()) return distanceTo(p1);
if(interx<p2.x()) return distanceTo(p2);
}
double intery=coeffp*(interx-_x)+_y;
if(p1.y()<p2.y()) {
if(intery>p2.y()) return distanceTo(p2);
if(intery<p1.y()) return distanceTo(p1);
} else {
if(intery>p1.y()) return distanceTo(p1);
if(intery<p2.y()) return distanceTo(p2);
}
return distanceTo(Point2D(interx,intery));
}
| void QGpCoreTools::Point2D::DMSToDegrees | ( | ) |
Converts point from geographical coordinates in D.MMSSsss... (degrees, minutes, seconds) to degrees.
Referenced by CoordReader::parse().
{
_x=Angle::DMSToDegrees(_x);
_y=Angle::DMSToDegrees(_y);
}
| void QGpCoreTools::Point2D::DMToDegrees | ( | ) |
Converts point from geographical coordinates in D.MMmm... (degrees, minutes) to degrees.
Referenced by CoordReader::parse().
{
_x=Angle::DMToDegrees(_x);
_y=Angle::DMToDegrees(_y);
}
| bool QGpCoreTools::Point2D::fromString | ( | const StringSection & | str | ) |
Extracts coordinates from string str.
Reimplemented in QGpCoreTools::Point.
References QGpCoreTools::StringSection::isValid(), QGpCoreTools::StringSection::nextField(), QGpCoreTools::StringSection::toDouble(), and TRACE.
Referenced by Group2PhaseReader::parse(), SciFigs::ParallelBand::xml_setProperty(), and SciFigs::ImageLayer::ReferencePoint::xml_setProperty().
| double QGpCoreTools::Point2D::geographicalAzimuthTo | ( | const Point2D & | p | ) | const |
Returns the azimuth from this point to p in radians from north, assuming that points have geographical coordinates. Returned values are between 0 and 2 PI.
The azimuth is defined by the great circle according to Snyder, J. P. (1987). Map Projections - A Working Manual, USGS Professional Paper 1395.
References QGpCoreTools::cos(), QGpCoreTools::Angle::degreesToRadians(), QGpCoreTools::sin(), and TRACE.
Referenced by CoordReader::parse().
{
TRACE;
double phi1=Angle::degreesToRadians(_y);
double phi2=Angle::degreesToRadians(p._y);
double cphi1=::cos(phi1);
double cphi2=::cos(phi2);
double sphi1=::sin(phi1);
double sphi2=::sin(phi2);
double deltaL=Angle::degreesToRadians(p._x-_x);
double a=::atan2(cphi2*::sin(deltaL), cphi1*sphi2-sphi1*cphi2*::cos(deltaL));
if(a<0.0) a+=2*M_PI;
return a;
}
| double QGpCoreTools::Point2D::geographicalDistanceTo | ( | const Point2D & | p | ) | const |
Returns the distance from this point to p in meters, assuming that points have geographical coordinates.
The distance is defined by the great circle according to Snyder, J. P. (1987). Map Projections - A Working Manual, USGS Professional Paper 1395. Same formula can be found at http://en.wikipedia.org/wiki/Haversine_formula.
Radius of earth taken at equator from WSG84 (http://en.wikipedia.org/wiki/World_Geodetic_System).
References QGpCoreTools::asin(), QGpCoreTools::cos(), QGpCoreTools::Angle::degreesToRadians(), QGpCoreTools::sin(), QGpCoreTools::sqrt(), and TRACE.
Referenced by CoordReader::parse().
{
TRACE;
double phi1=Angle::degreesToRadians(_y);
double phi2=Angle::degreesToRadians(p._y);
double deltaL=Angle::degreesToRadians(p._x-_x);
double sdeltaphi=::sin(0.5*(phi2-phi1));
sdeltaphi*=sdeltaphi;
double sdeltalamda=::sin(0.5*deltaL);
sdeltalamda*=sdeltalamda;
double sc=sdeltaphi+::cos(phi1)*::cos(phi2)*sdeltalamda;
if(sc<0.0) sc=0.0;
sc=::sqrt(sc);
if(sc>1.0) sc=1.0;
return 6378137.0*2.0*::asin(sc);
}
| void QGpCoreTools::Point2D::geographicalToRectangular | ( | const Point2D & | ref, |
| double | rotation = 0.0 |
||
| ) |
Converts point from geographical coordinates (long and lat, WGS84) to rectangular coordinates in meters, using reference point ref (in geographical coordinates). Longitude is along X axis, latitude is along Y axis.
rotation is an angle to rotate north reference (clockwize, in degrees).
Kindly provided by Matthias Ohrnberger.
References QGpCoreTools::cos(), DRAD, DRLT, E_FLAT, E_RAD, QGpCoreTools::sin(), QGpCoreTools::tan(), x(), and y().
Referenced by CoordReader::parse().
{
double olon;
double olat;
double lat_fac; // conversion factor for latitude in km
double lon_fac; // conversion factor for longitude in km
double snr; // sin of rotation angle
double csr; // cos of rotation angle
double dlt1;
double dlt2;
double del;
double radius;
double tmp;
double tmp_x, tmp_y;
/* convert everything to minutes */
olon=60.0 * ref.x();
olat=60.0 * ref.y();
/* latitude */
dlt1=::atan(DRLT * ::tan((double)olat * DRAD/60.0));
dlt2=::atan(DRLT * ::tan(((double)olat +1.0) * DRAD/60.0));
del =dlt2 - dlt1;
radius=E_RAD * (1.0 - (::sin(dlt1)*::sin(dlt1)/E_FLAT));
lat_fac=radius * del;
/* longitude */
del=acos(1.0 - (1.0 - ::cos(DRAD/60.0)) * ::cos(dlt1) * ::cos(dlt1));
dlt2=radius * del;
lon_fac=dlt2/::cos(dlt1);
/* rotation */
snr=::sin((double)rotation * DRAD);
csr=::cos((double)rotation * DRAD);
_y *= 60.0;
_x *= 60.0;
tmp_x=_x - olon;
tmp_y=_y - olat;
tmp =::atan(DRLT * ::tan(DRAD * (_y+olat)/120.0));
tmp_x=(double)tmp_x * lon_fac * ::cos(tmp);
tmp_y=(double)tmp_y * lat_fac;
_x=(csr*tmp_x - snr*tmp_y)*1000.0;
_y=(csr*tmp_y + snr*tmp_x)*1000.0;
}
| void QGpCoreTools::Point2D::geographicalToUtm | ( | const QPoint * | zone = 0 | ) |
Converts from latitude, longitude to Universal Transverse Mercator (UTM). The implemented method delivers a cm precision. Results are in agreement with Google Earth coordinates with a 1-cm error.
Equations from USGS Bulletin 1532 (p 63 to 69), by Snyder (1982).
If zone is not null, the zone is forced to zone.
References A, QGpCoreTools::cos(), QGpCoreTools::Angle::degreesToRadians(), setX(), setY(), QGpCoreTools::sin(), QGpCoreTools::sqrt(), TRACE, utmLongitudeReference(), utmZoneIndex(), x(), and y().
Referenced by CoordReader::parse().
{
TRACE;
double phi=Angle::degreesToRadians(y()); // Latitude (radians)
double lamda=Angle::degreesToRadians(x()); // Longitude (radians)
// Find reference longitude according to UTM band
QPoint zoneIndex;
if(zone) {
zoneIndex=*zone;
} else {
zoneIndex=utmZoneIndex();
}
double lamda0=Angle::degreesToRadians(utmLongitudeReference(zoneIndex)); // Reference longitude (radians)
//double phi0=utmLatitudeReference(zoneIndex);
double e=0.0818192; // Eccentricity
double a=6378137; // Equatorial radius in m for reference ellipsoide
double sphi=::sin(phi);
double cphi=::cos(phi);
double e2=e*e;
double ep2=e2/(1.0-e2);
double N=a/::sqrt(1.0-e2*sphi*sphi);
double T=sphi/cphi;
double T2=T*T;
double T4=T2*T2;
double C=ep2*cphi*cphi;
double A=(lamda-lamda0)*cphi;
double M=utmM(phi);
//double M0=utmM(phi0);
double N0;
if(phi>0.0) {
N0=0.0; // Northern hemisphere in m
} else {
N0=10000000.0; // Southern hemisphere in m
}
double k0=0.9996; // Scale on central meridian, standard value for UTM
double E0=500000.0; // in m
double A2=A*A;
double A3=A2*A;
double A4=A3*A;
double A5=A4*A;
double A6=A5*A;
setX(E0+k0*N*(A+(1.0-T2+C)*A3/6.0+(5.0-18.0*T2+T4+72.0*C-58.0*ep2)*A5/120.0));
// M-M0 is replaced by M, mistake in Snyder?
setY(N0+k0*(M+N*T*(A2*0.5+(5.0-T2+9.0*C+4.0*C*C)*A4/24.0+(61.0-58.0*T2+T4)*A6/720.0)));
}
| void QGpCoreTools::Point2D::interpole | ( | double | valX, |
| const Point2D & | p1, | ||
| const Point2D & | p2 | ||
| ) | [inline] |
| bool QGpCoreTools::Point2D::isHit | ( | const Point2D & | p, |
| double | tolX, | ||
| double | tolY | ||
| ) | const |
| bool QGpCoreTools::Point2D::isSimilar | ( | const Point2D & | p, |
| double | tolerance | ||
| ) | const [inline] |
| bool QGpCoreTools::Point2D::isValid | ( | ) | const [inline] |
| void QGpCoreTools::Point2D::move | ( | double | distance, |
| const Angle & | azimuth | ||
| ) |
Moves this point by distance in direction azimuth (radians, mathematical sense, 0 is east).
References QGpCoreTools::Angle::cos(), and QGpCoreTools::Angle::sin().
{
_x+=distance*azimuth.cos();
_y+=distance*azimuth.sin();
}
| bool QGpCoreTools::Point2D::operator!= | ( | const Point2D & | obj | ) | const [inline] |
{return (! (*this==obj));}
| Point2D QGpCoreTools::Point2D::operator* | ( | double | mul | ) | const [inline] |
Reimplemented in QGpCoreTools::Point, and QGpCoreTools::NamedPoint.
References Point2D().
{
return Point2D (_x*mul, _y*mul);
}
| void QGpCoreTools::Point2D::operator*= | ( | const Point2D & | p | ) | [inline] |
{
_x*=p._x;
_y*=p._y;
}
| void QGpCoreTools::Point2D::operator*= | ( | double | mul | ) | [inline] |
Reimplemented in QGpCoreTools::Point.
{
_x*=mul;
_y*=mul;
}
| void QGpCoreTools::Point2D::operator*= | ( | const Matrix2x2 & | transformation | ) |
{
*this=transformation*(*this);
}
| void QGpCoreTools::Point2D::operator+= | ( | const Point2D & | p | ) | [inline] |
{
_x+=p._x;
_y+=p._y;
}
| void QGpCoreTools::Point2D::operator-= | ( | const Point2D & | p | ) | [inline] |
{
_x-=p._x;
_y-=p._y;
}
| Point2D QGpCoreTools::Point2D::operator/ | ( | double | div | ) | const [inline] |
Reimplemented in QGpCoreTools::Point, and QGpCoreTools::NamedPoint.
References Point2D().
{
double f=1.0/div;
return Point2D (_x*f, _y*f);
}
| void QGpCoreTools::Point2D::operator/= | ( | const Point2D & | p | ) | [inline] |
{
_x/=p._x;
_y/=p._y;
}
| void QGpCoreTools::Point2D::operator/= | ( | double | div | ) | [inline] |
Reimplemented in QGpCoreTools::Point.
{
double f=1.0/div;
_x*=f;
_y*=f;
}
| bool QGpCoreTools::Point2D::operator< | ( | const Point2D & | p | ) | const [inline] |
{
if(_x < p._x) return true;
else if(_x==p._x && _y<p._y) return true;
return false;
}
| bool QGpCoreTools::Point2D::operator<= | ( | const Point2D & | obj | ) | const [inline] |
{return (*this==obj || *this < obj);}
Reimplemented in QGpCoreTools::Point, QGpCoreTools::NamedPoint, and TapePoint.
{
_x=p._x;
_y=p._y;
return *this;
}
| bool QGpCoreTools::Point2D::operator== | ( | const Point2D & | p | ) | const [inline] |
{
if((_x==p._x) &&
(_y==p._y)) return true;
return false;
}
| bool QGpCoreTools::Point2D::operator> | ( | const Point2D & | obj | ) | const [inline] |
{return (! (*this <= obj));}
| bool QGpCoreTools::Point2D::operator>= | ( | const Point2D & | obj | ) | const [inline] |
{return (! (*this < obj));}
| void QGpCoreTools::Point2D::rectangularToGeographical | ( | const Point2D & | ref, |
| double | rotation = 0.0 |
||
| ) |
Converts point from rectangular coordinates in meters to geographical coordinates (lat and long, WGS84), using reference point ref (in geographical coordinates). Longitude is along X axis, latitude is along Y axis.
rotation is an angle to rotate north reference (clockwize in degrees).
Kindly provided by Matthias Ohrnberger.
References QGpCoreTools::cos(), DRAD, DRLT, E_FLAT, E_RAD, QGpCoreTools::sin(), QGpCoreTools::tan(), x(), and y().
Referenced by CoordReader::parse(), and QGpGuiTools::CoordinateFile::write().
{
double olon;
double olat;
double lat_fac; // conversion factor for latitude in km
double lon_fac; // conversion factor for longitude in km
double snr; // sin of rotation angle
double csr; // cos of rotation angle
double dlt1;
double dlt2;
double del;
double radius;
double tmp;
double tmp_x, tmp_y;
// Convert to km
_x*=0.001;
_y*=0.001;
/* convert everything to minutes */
olon=60.0 * ref.x();
olat=60.0 * ref.y();
/* latitude */
dlt1=atan(DRLT * ::tan((double)olat * DRAD/60.0));
dlt2=atan(DRLT * ::tan(((double)olat +1.0) * DRAD/60.0));
del =dlt2 - dlt1;
radius=E_RAD * (1.0 - (::sin(dlt1)*::sin(dlt1)/E_FLAT));
lat_fac=radius * del;
/* longitude */
del=acos(1.0 - (1.0 - ::cos(DRAD/60.0)) * ::cos(dlt1) * ::cos(dlt1));
dlt2=radius * del;
lon_fac=dlt2/::cos(dlt1);
/* rotation */
snr=::sin((double)rotation * DRAD);
csr=::cos((double)rotation * DRAD);
tmp_x=snr* _y + csr* _x;
tmp_y=csr* _y - snr* _x;
tmp_y=tmp_y/lat_fac;
tmp_y += olat;
tmp=::atan(DRLT * ::tan(DRAD * (tmp_y+olat)/120.0));
tmp_x=tmp_x/(lon_fac * ::cos(tmp));
tmp_x += olon;
_x=tmp_x/60.0;
_y=tmp_y/60.0;
}
| void QGpCoreTools::Point2D::rotate | ( | const Angle & | a | ) |
Rotates point by angle degrees clockwize around Z axis.
References QGpCoreTools::Angle::cos(), QGpCoreTools::Angle::sin(), and x().
Referenced by SciFigs::SlopeEstimator::paintData(), CoordReader::parse(), and QGpCoreTools::Segment2D::rotated().
| double QGpCoreTools::Point2D::scalarProduct | ( | const Point2D & | p | ) | const [inline] |
{
return _x*p._x+_y*p._y;
}
| void QGpCoreTools::Point2D::set | ( | double | xi, |
| double | yi | ||
| ) | [inline] |
{
_x=xi;
_y=yi;
}
| void QGpCoreTools::Point2D::setValid | ( | bool | ) | [inline] |
Reimplemented in QGpCoreTools::Point, and Sample.
{}
| void QGpCoreTools::Point2D::setX | ( | double | v | ) | [inline] |
Referenced by SciFigs::AxisWindow::alignHScales(), QGpCoreTools::Histogram::boxes(), SciFigs::GraphContent::coordinateTipText(), createDots(), SpacSelector::createObjects(), QGpCoreTools::IrregularGrid2DData::crossSection(), SciFigs::GraphicSheet::currentOrigin(), QGpCoreWave::ShAmplification::curve(), QGpCoreTools::Segment2D::distanceTo(), ArrayCore::FKHorizontal::FKHorizontal(), QGpCoreTools::NamedPoint::fromString(), QGpCoreTools::Point::fromString(), geographicalToUtm(), SciFigs::ImageLayer::ImageLayer(), QGpCoreTools::IrregularGrid2DData::integralCrossSection(), QGpCoreTools::Point::interpole(), interpole(), QGpCoreTools::Line2D::intersect(), QGpCoreTools::Rect::intersect(), QGpCoreTools::Segment2D::intersects(), SciFigs::CircleViewer::Limits::Limits(), main(), DinverDCGui::GroundModelViewer::minMaxProfiles(), QGpCoreTools::Histogram::normalize(), QGpCoreTools::Matrix2x2::operator*(), QGpCoreTools::Matrix3x3::operator*(), QGpCoreTools::Matrix4x4::operator*(), QGpCoreTools::operator>>(), SciFigs::XUniqueYColorLines::paintData(), CoordReader::parse(), SciFigs::GraphicSheet::paste(), SciFigs::ComplexStatisticalLine::point(), SciFigs::RealStatisticalLine::point(), QGpGuiWave::MagnetoTelluricLine::point(), QGpGuiWave::ModalLine::point(), QGpGuiWave::RefractionLine::point(), QGpCoreTools::PointLocate::position(), Sample::read(), GeopsyCore::Signal::read(), QGpGuiTools::GeographicalReference::reference(), DinverDCGui::RefractionViewer::report2plot(), DinverDCGui::MagnetoTelluricViewer::report2plot(), QGpCoreWave::Profile::report2plot(), ArrayCore::FKRadial::rotationFactors(), ArrayCore::FKTransverse::rotationFactors(), SciFigs::ImageLayer::scaling(), Measurement::set(), QGpCoreTools::Segment2D::set(), SourceParameters::setAzimuth(), ToolSPAC::setCoArrayMap(), SciFigs::PlotLine2D::setCurve(), TapeCoordinateItem::setData(), GeopsyGui::StationCoordinatesItem::setData(), SourceParameters::setDistance(), QGpGuiWave::DispersionLimitLayer::setFrequencySampling(), GeopsyCore::Signal::setHeader(), QGpCoreTools::Plane::setNormalVector(), QGpCoreTools::Plane::setNormalVectorXY(), CoordReader::setOptions(), SignalReader::setOptions(), SciFigs::GraphicObject::setPrintAnchor(), SciFigs::GraphicObject::setPrintLeft(), SciFigs::GraphicObject::setPrintRight(), SciFigs::GraphicObject::setPrintXAnchor(), SciFigs::ImageLayer::setProperty(), GeopsyCore::SignalHeaderObject::setReceiverX(), SourceParameters::setSourceX(), GeopsyCore::SignalHeaderObject::setSourceX(), TapePositioningSystem::Triangulator::solutions(), DinverDCCore::ModalStorageReader::toPlot(), QGpCoreWave::MagnetoTelluricFactory::toStream(), QGpCoreWave::RefractionFactory::toStream(), utmToGeographical(), and QGpCoreTools::Point::vectorialProduct().
{_x=v;}
| void QGpCoreTools::Point2D::setY | ( | double | v | ) | [inline] |
Referenced by QGpCoreTools::Histogram::addValue(), SciFigs::AxisWindow::alignVScales(), SciFigs::GraphContent::coordinateTipText(), createDots(), QGpCoreTools::IrregularGrid2DData::crossSection(), SciFigs::GraphicSheet::currentOrigin(), NewNoiseModel::curve(), QGpCoreWave::ShAmplification::curve(), QGpCoreTools::Segment2D::distanceTo(), SciFigs::XYColorLines::divYbyX(), ArrayCore::FKHorizontal::FKHorizontal(), QGpCoreTools::NamedPoint::fromString(), QGpCoreTools::Point::fromString(), geographicalToUtm(), QGpCoreTools::GridSearch::globalMax(), SciFigs::ImageLayer::ImageLayer(), QGpCoreTools::IrregularGrid2DData::integralCrossSection(), QGpCoreTools::Point::interpole(), interpole(), QGpCoreTools::Line2D::intersect(), QGpCoreTools::Rect::intersect(), QGpCoreTools::Segment2D::intersects(), SciFigs::CircleViewer::Limits::Limits(), main(), SciFigs::XYColorLines::mulYbyX(), QGpCoreTools::Histogram::normalize(), StatGridAnalyser::on_freqScroll_valueChanged(), QGpCoreTools::Matrix2x2::operator*(), QGpCoreTools::Matrix3x3::operator*(), QGpCoreTools::Matrix4x4::operator*(), QGpCoreTools::operator>>(), SciFigs::XUniqueYColorLines::paintData(), CoordReader::parse(), SciFigs::GraphicSheet::paste(), SciFigs::ComplexStatisticalLine::point(), SciFigs::RealStatisticalLine::point(), QGpGuiWave::MagnetoTelluricLine::point(), QGpGuiWave::ModalLine::point(), QGpGuiWave::RefractionLine::point(), QGpCoreTools::PointLocate::position(), Sample::read(), GeopsyCore::Signal::read(), QGpGuiTools::GeographicalReference::reference(), DinverDCGui::RefractionViewer::report2plot(), DinverDCGui::MagnetoTelluricViewer::report2plot(), QGpCoreWave::Profile::report2plot(), ArrayCore::FKTransverse::rotationFactors(), ArrayCore::FKRadial::rotationFactors(), SciFigs::ImageLayer::scaling(), Measurement::set(), QGpCoreTools::Segment2D::set(), SourceParameters::setAzimuth(), ToolSPAC::setCoArrayMap(), SciFigs::PlotLine2D::setCurve(), TapeCoordinateItem::setData(), GeopsyGui::StationCoordinatesItem::setData(), SourceParameters::setDistance(), GeopsyCore::Signal::setHeader(), QGpCoreTools::Plane::setNormalVector(), QGpCoreTools::Plane::setNormalVectorXY(), CoordReader::setOptions(), SignalReader::setOptions(), SciFigs::GraphicObject::setPrintAnchor(), SciFigs::GraphicObject::setPrintBottom(), SciFigs::GraphicObject::setPrintTop(), SciFigs::GraphicObject::setPrintYAnchor(), SciFigs::ImageLayer::setProperty(), GeopsyCore::SignalHeaderObject::setReceiverY(), SourceParameters::setSourceY(), GeopsyCore::SignalHeaderObject::setSourceY(), TapePositioningSystem::Triangulator::solutions(), DinverDCCore::ModalStorageReader::toPlot(), QGpCoreWave::MagnetoTelluricFactory::toStream(), QGpCoreWave::RefractionFactory::toStream(), utmToGeographical(), and QGpCoreTools::Point::vectorialProduct().
{_y=v;}
| void QGpCoreTools::Point2D::setY | ( | double | v, |
| const CurvePointOptions * | |||
| ) | [inline] |
{_y=v;}
| QString QGpCoreTools::Point2D::toString | ( | int | precision = 6, |
| char | format = 'g' |
||
| ) | const |
Returns the point as a string with space separation between coordinates. precision is the number of significant digits. format is 'g' or 'f'.
Reimplemented in QGpCoreTools::Point, and QGpCoreTools::NamedPoint.
References TRACE.
Referenced by QGpGuiTools::GeographicalReference::save(), QGpCoreTools::Ellipse::toString(), and SciFigs::ImageLayer::ReferencePoint::xml_writeProperties().
| void QGpCoreTools::Point2D::translate | ( | const Point2D & | p | ) | [inline] |
Referenced by TapePositioningSystem::Cluster::setCoordinateSystem().
{operator+=(p);}
| double QGpCoreTools::Point2D::utmLatitudeReference | ( | QPoint | z | ) | [static] |
Returns the latitude reference of UTM zone with index z. Returns 1000.0 for unvalid index.
Referenced by utmToGeographical().
{
return (z.y() << 3)-76.0;
}
| double QGpCoreTools::Point2D::utmLongitudeReference | ( | QPoint | z | ) | [static] |
Returns the longitude reference of UTM zone with index z. Returns 1000.0 for invalid index.
Referenced by geographicalToUtm(), and utmToGeographical().
{
if(z.y()<0 || z.y()>19) {
return 1000.0;
} else {
// Norvegian exceptions...
switch(z.y()) {
case 17:
switch(z.x()) {
case 31: return 1.5;
case 32: return 7.5;
default: break;
}
break;
case 19:
switch(z.x()) {
case 31: return 4.5;
case 32: return 1000.0;
case 33: return 15.0;
case 34: return 1000.0;
case 35: return 27.0;
case 36: return 1000.0;
case 37: return 37.5;
default: break;
}
break;
default:
break;
}
if(z.x()<1 || z.x()>60) {
return 1000.0;
} else {
return z.x()*6.0-183.0;
}
}
}
| void QGpCoreTools::Point2D::utmToGeographical | ( | const QPoint & | zone | ) |
Convert UTM coordinates into geographical coordinates (WSG84)
Equations from USGS Bulletin 1532 (p 63 to 69), by Snyder (1982).
References QGpCoreTools::cos(), QGpCoreTools::Angle::degreesToRadians(), QGpCoreTools::Angle::radiansToDegrees(), setX(), setY(), QGpCoreTools::sin(), QGpCoreTools::sqrt(), TRACE, utmLatitudeReference(), utmLongitudeReference(), x(), and y().
Referenced by CoordReader::parse().
{
TRACE;
double lamda0=Angle::degreesToRadians(utmLongitudeReference(zone)); // Reference longitude (radians)
double phi0=Angle::degreesToRadians(utmLatitudeReference(zone)); // Reference latitude (radians)
//double M0=utmM(phi0);
double N0;
if(phi0>0.0) {
N0=0.0; // Northern hemisphere in m
} else {
N0=10000000.0; // Southern hemisphere in m
}
double E0=500000.0; // in m
double e=0.0818192; // Eccentricity
double a=6378137; // Equatorial radius in m for reference ellipsoide
double k0=0.9996; // Scale on central meridian, standard value for UTM
double e2=e*e;
double e4=e2*e2;
double e6=e4*e2;
double ep2=e2/(1.0-e2);
// M-M0 is replaced by M, mistake in Snyder? Original formula: M=M0+(y()-N0)/k0
double M=(y()-N0)/k0;
double mu=M/(a*(1.0-e2*0.25-3.0/64.0*e4-5.0/256.0*e6));
double e2m1=1.0-e2;
double sqrte2=::sqrt(e2m1);
double e1=(1.0-sqrte2)/(1.0+sqrte2);
double e12=e1*e1;
double e13=e12*e1;
double e14=e13*e1;
double phi1=mu+(3.0/2.0*e1-27.0/32.0*e13)*::sin(2.0*mu)+(21.0/16.0*e12-55.0/32.0*e14)*::sin(4.0*mu)+(151.0/96.0*e13)*::sin(6.0*mu)+
(1097.0/512.0*e14)*::sin(8.0*mu);
double sphi1=::sin(phi1);
double cphi1=::cos(phi1);
double C1=ep2*cphi1*cphi1;
double C12=C1*C1;
double T1=sphi1/cphi1;
double T12=T1*T1;
double T14=T12*T12;
double esphim1=1.0-e2*sphi1*sphi1;
double N1=a/::sqrt(esphim1);
double R1=a*e2m1/::pow(esphim1, 1.5);
double D=(x()-E0)/(N1*k0);
double D2=D*D;
double D3=D2*D;
double D4=D2*D2;
double D5=D4*D;
double D6=D3*D3;
setY(Angle::radiansToDegrees(phi1-(N1*T1/R1)*(D2*0.5-(5.0+3*T12+10*C1-4.0*C12-9.0*ep2)*D4/24.0+
(61.0+90.0*T12+298.0*C1+45.0*T14+252.0*ep2-3.0*C12)*D6/720.0)));
setX(Angle::radiansToDegrees(lamda0+(D-(1.0+2.0*T12+C1)*D3/6.0+(5.0-2.0*C1+28.0*T12-3.0*C12+8.0*ep2+24.0*T14)*D5/120.0)/cphi1));
}
| QString QGpCoreTools::Point2D::utmZone | ( | QPoint | z, |
| bool * | ok = 0 |
||
| ) | [static] |
Returns UTM zone name (3 characters) from a UTM zone index. ok is set to false if z is not a valid zone.
{
QString zone;
if(z.x()<1 || z.x()>60) {
if(ok) *ok=false;
} else {
if(ok) *ok=true;
zone=QString::number(z.x());
switch(z.y()) {
case 0:
case 1:
case 2:
case 3:
case 4:
case 5:
zone+=QChar(z.y()+'C');
break;
case 6:
case 7:
case 8:
case 9:
case 10:
zone+=QChar(z.y()-6+'J');
break;
case 11:
case 12:
case 13:
case 14:
case 15:
case 16:
case 17:
case 18:
case 19:
zone+=QChar(z.y()-11+'P');
break;
default:
if(ok) *ok=false;
zone=QString::null;
break;
}
}
return zone;
}
| QPoint QGpCoreTools::Point2D::utmZone | ( | const QString & | z, |
| bool * | ok = 0 |
||
| ) | [static] |
Returns UTM zone indexes from a UTM zone name (3 characters). ok is set to false if z is not a valid zone.
{
QPoint p;
bool intOk;
if(z.isEmpty()) {
if(ok) *ok=false;
return p;
}
char l=z[z.count()-1].toUpper().toAscii();
switch(l) {
case 'C': // 0
case 'D':
case 'E':
case 'F':
case 'G':
case 'H': // 5
p.setY(l-'C');
break;
case 'J': // 6
case 'K':
case 'L':
case 'M':
case 'N': // 10
p.setY(l-'J'+6);
break;
case 'P': // 11
case 'Q':
case 'R':
case 'S':
case 'T':
case 'U':
case 'V':
case 'W':
case 'X': // 19
p.setY(l-'P'+11);
break;
default:
if(ok) *ok=false;
return p;
}
p.setX(z.left(z.count()-1).toInt(&intOk));
if(!intOk || p.y()<0 || p.y()>19 || p.x()<1 || p.x()>60) {
if(ok) *ok=false;
} else {
if(ok) *ok=true;
}
return p;
}
| QPoint QGpCoreTools::Point2D::utmZoneIndex | ( | ) | const |
Returns the Universal Transverse Mercator zone indexes of this point. x is the longitude from 1 to 60, y is the latitude band from 0 to 19.
Referenced by geographicalToUtm(), and CoordReader::parse().
{
QPoint z;
if(y()<-80.0 || y()>84.0) {
z.setX(-1);
z.setY(-1);
} else {
z.setX(((int)floor(x()+180.0)/6.0)+1);
z.setY(((int)floor(y()+80.0) >> 3));
// Norvegian exceptions...
switch(z.y()) {
case 17:
switch(z.x()) {
case 31:
if(x()>3.0) {
z.setX(32); // From 31V to 32V wich is 9 degree wide instead of 6 degrees
}
break;
default:
break;
}
break;
case 19:
switch(z.x()) {
case 32:
if(x()<9.0) {
z.setX(31);
} else {
z.setX(33);
}
break;
case 34:
if(x()<21.0) {
z.setX(33);
} else {
z.setX(35);
}
break;
case 36:
if(x()<33.0) {
z.setX(35);
} else {
z.setX(37);
}
break;
default:
break;
}
break;
default:
break;
}
}
return z;
}
| double QGpCoreTools::Point2D::x | ( | ) | const [inline] |
Referenced by RealTimeHistogram::add(), QGpCoreTools::Parallelepiped::add(), QGpCoreTools::Rect::add(), addDotPlot(), SciFigs::AxisWindow::alignHScales(), azimuthTo(), T0GridSearch::bestT0(), SciFigs::CircleMask::boundingRect(), SciFigs::CircleViewer::boundingRect(), SciFigs::ImageLayer::boundingRect(), SciFigs::XYColorLines::boundingRect(), SciFigs::LineLayer::boundingRect(), QGpCoreTools::Histogram::boxes(), QGpCoreTools::Circle::Circle(), SciFigs::ImageLayer::colorHistogram(), DampingResults::compute(), SciFigs::ColorPaletteLayer::coordinateTipInfo(), SciFigs::IrregularGrid2DPlot::coordinateTipInfo(), SciFigs::GraphContentLayer::coordinateTipInfo(), SciFigs::GraphContent::coordinateTipText(), Histogram2D::countSamples(), TapePositioningSystem::Triangulator::covariance(), QGpCoreTools::IrregularGrid2DData::crossSection(), NewNoiseModel::curve(), SourceItemModel::data(), TapeCoordinateItem::data(), GeopsyGui::StationCoordinatesItem::data(), SciFigs::LineItem::data(), QGpCoreTools::Segment2D::distanceTo(), QGpCoreTools::Point::distanceTo(), distanceTo(), QGpCoreTools::Point::distanceToSegment(), distanceToSegment(), SciFigs::XYColorLines::divYbyX(), SciFigs::LineLayer::divYbyX(), SciFigs::GridPlot::drawGrid2DBlock(), SciFigs::GridPlot::drawGrid2DYSmooth(), QGpCoreTools::Point::elevationTo(), geographicalToRectangular(), geographicalToUtm(), geopsyToFile(), HRFKLoopTask::getPower(), FKLoopTask::getPower(), GeopsyCore::Signal::header(), QGpCoreTools::Parallelepiped::includes(), QGpCoreTools::Segment::includes(), QGpCoreTools::Segment2D::includes(), QGpCoreTools::Rect::includes(), PhaseShifter::initGrid(), QGpCoreTools::Point::interpole(), interpole(), QGpCoreTools::Line2D::intersect(), QGpCoreTools::Rect::intersect(), QGpCoreTools::Segment2D::intersects(), QGpCoreTools::Segment::intersects(), isHit(), isSimilar(), SciFigs::CircleViewer::Limits::Limits(), main(), DinverDCGui::GroundModelViewer::minMaxProfiles(), SciFigs::LineEditor::mouseMoveEvent(), SciFigs::GridViewer::mouseReleaseEvent(), SciFigs::IrregularGrid2DPlot::mouseReleaseEvent(), SciFigs::LineLayer::mouseReleaseEvent(), SciFigs::XYColorLines::mulYbyX(), SciFigs::LineLayer::mulYbyX(), LinearFKActiveArrayStations::name(), QGpCoreTools::Histogram::normalize(), StatGridAnalyser::on_freqScroll_valueChanged(), QGpCoreTools::NamedPoint::operator*(), QGpCoreTools::Matrix2x2::operator*(), QGpCoreTools::Matrix3x3::operator*(), QGpCoreTools::Point::operator*(), QGpCoreTools::Matrix4x4::operator*(), QGpCoreTools::NamedPoint::operator+(), QGpCoreTools::Point::operator+(), QGpCoreTools::NamedPoint::operator-(), QGpCoreTools::Point::operator-(), QGpCoreTools::NamedPoint::operator/(), QGpCoreTools::Point::operator/(), QGpCoreTools::Point::operator<(), QGpCoreTools::operator<<(), QGpCoreTools::operator>>(), SciFigs::CircleMask::paintData(), SciFigs::CircleViewer::paintData(), SciFigs::ImageLayer::paintData(), SciFigs::SlopeEstimator::paintText(), QGpCoreTools::Parallelepiped::Parallelepiped(), CoordReader::parse(), SciFigs::GraphicSheet::paste(), Point2D(), QGpCoreTools::PointLocate::position(), SciFigs::GraphicObject::printLeft(), SciFigs::GraphicObject::printRight(), SciFigs::ImageLayer::properties(), QGpCoreTools::qHash(), SciFigs::GraphContentOptions::r2s(), SciFigs::GraphContentOptions::r2sF(), GeopsyCore::SignalHeaderObject::receiverX(), QGpCoreTools::Rect::Rect(), rectangularToGeographical(), QGpCompatibility::CompatMultiModalCurves::refinesToReport(), DinverDCGui::MagnetoTelluricViewer::report2plot(), rotate(), GeopsyCore::SubSignalPool::rotateComponents(), QGpCoreTools::Point::scalarProduct(), SciFigs::ImageLayer::scaling(), QGpCoreTools::Segment::set(), QGpCoreTools::Segment2D::set(), QGpCompatibility::CompatFunction::set(), MonoStation::AbstractSummary::setBubbleValues(), ToolSPAC::setCoArrayMap(), GeopsyCore::SEGYTraceHeader::setCoordinateFactor(), ArrayCore::FKStationSignals::setCurrentShift(), ArrayCore::FKStationSignals::setPhaseShift(), SciFigs::GraphicObject::setPrintAnchor(), GeopsyCore::SEGYTraceHeader::setReceiver(), QGpGuiTools::GeographicalReference::setReference(), GeopsyCore::SEGYTraceHeader::setSource(), ArrayCore::StationCouple::setStations(), TapePositioningSystem::Triangulator::solutions(), GeopsyCore::SignalHeaderObject::sourceX(), CoordReader::terminate(), QGpCompatibility::CompatMultiModalCurves::toPointVector(), QGpCoreWave::MagnetoTelluricFactory::toStream(), QGpCoreTools::Point::toString(), utmToGeographical(), utmZoneIndex(), QGpCoreWave::TheoreticalFK::value(), ArrayCore::FKHorizontal::value(), PhaseShifter::value(), QGpCoreTools::Point::vectorialProduct(), Acquisition::write(), GeopsyCore::Signal::writeSac(), GeopsyCore::Signal::writeSeg2(), SciFigs::ParallelBand::xml_setProperty(), SciFigs::ImageLayer::ReferencePoint::xml_setProperty(), SciFigs::CircleViewer::xml_writeProperties(), SciFigs::GraphContentOptions::yr2s(), and QGpCoreTools::Plane::z().
{return _x;}
| double QGpCoreTools::Point2D::y | ( | ) | const [inline] |
Referenced by RealTimeHistogram::add(), QGpCoreTools::Parallelepiped::add(), QGpCoreTools::Rect::add(), addDotPlot(), QGpCoreTools::Histogram::addValue(), SciFigs::AxisWindow::alignVScales(), azimuthTo(), SciFigs::CircleMask::boundingRect(), SciFigs::CircleViewer::boundingRect(), SciFigs::ImageLayer::boundingRect(), SciFigs::XYColorLines::boundingRect(), SciFigs::LineLayer::boundingRect(), QGpCoreTools::Circle::Circle(), SciFigs::ImageLayer::colorHistogram(), DampingResults::compute(), SciFigs::IrregularGrid2DPlot::coordinateTipInfo(), SciFigs::GraphContentLayer::coordinateTipInfo(), SciFigs::GraphContent::coordinateTipText(), Histogram2D::countSamples(), TapePositioningSystem::Triangulator::covariance(), QGpCoreTools::IrregularGrid2DData::crossSection(), SourceItemModel::data(), TapeCoordinateItem::data(), GeopsyGui::StationCoordinatesItem::data(), SciFigs::LineItem::data(), QGpCoreTools::Segment2D::distanceTo(), QGpCoreTools::Point::distanceTo(), distanceTo(), QGpCoreTools::Point::distanceToSegment(), distanceToSegment(), SciFigs::XYColorLines::divYbyX(), SciFigs::LineLayer::divYbyX(), SciFigs::GridPlot::drawGrid2DBlock(), SciFigs::GridPlot::drawGrid2DYSmooth(), QGpCoreTools::Point::elevationTo(), geographicalToRectangular(), geographicalToUtm(), geopsyToFile(), HRFKLoopTask::getPower(), FKLoopTask::getPower(), GeopsyCore::Signal::header(), MonoStation::StationResults::highestPeak(), QGpCoreTools::Parallelepiped::includes(), QGpCoreTools::Segment::includes(), QGpCoreTools::Segment2D::includes(), QGpCoreTools::Rect::includes(), PhaseShifter::initGrid(), QGpCoreTools::Point::interpole(), interpole(), QGpCoreTools::Rect::intersect(), QGpCoreTools::Segment2D::intersects(), QGpCoreTools::Segment::intersects(), isHit(), isSimilar(), QGpCoreTools::Histogram::limits(), SciFigs::CircleViewer::Limits::Limits(), MonoStation::StationResults::lowestPeak(), main(), MonoStation::StationResults::maximumAmplitudePeak(), SciFigs::LineEditor::mouseMoveEvent(), SciFigs::GridViewer::mouseReleaseEvent(), SciFigs::IrregularGrid2DPlot::mouseReleaseEvent(), SciFigs::XYColorLines::mulYbyX(), SciFigs::LineLayer::mulYbyX(), LinearFKActiveArrayStations::name(), QGpCoreTools::Histogram::normalize(), StatGridAnalyser::on_freqScroll_valueChanged(), QGpCoreTools::NamedPoint::operator*(), QGpCoreTools::Matrix2x2::operator*(), QGpCoreTools::Matrix3x3::operator*(), QGpCoreTools::Point::operator*(), QGpCoreTools::Matrix4x4::operator*(), QGpCoreTools::NamedPoint::operator+(), QGpCoreTools::Point::operator+(), QGpCoreTools::NamedPoint::operator-(), QGpCoreTools::Point::operator-(), QGpCoreTools::NamedPoint::operator/(), QGpCoreTools::Point::operator/(), QGpCoreTools::Point::operator<(), QGpCoreTools::operator<<(), QGpCoreTools::operator>>(), SciFigs::CircleMask::paintData(), SciFigs::CircleViewer::paintData(), SciFigs::ImageLayer::paintData(), SciFigs::SlopeEstimator::paintText(), QGpCoreTools::Parallelepiped::Parallelepiped(), CoordReader::parse(), SciFigs::GraphicSheet::paste(), MonoStation::StationResults::peak(), Point2D(), QGpCoreTools::PointLocate::position(), SciFigs::GraphicObject::printBottom(), SciFigs::GraphicObject::printTop(), SciFigs::ImageLayer::properties(), QGpCoreTools::qHash(), SciFigs::GraphContentOptions::r2s(), SciFigs::GraphContentOptions::r2sF(), GeopsyCore::SignalHeaderObject::receiverY(), QGpCoreTools::Rect::Rect(), rectangularToGeographical(), QGpCompatibility::CompatMultiModalCurves::refinesToReport(), QGpCoreTools::Curve< pointType >::resample(), Process::run(), QGpCoreTools::Point::scalarProduct(), SciFigs::ImageLayer::scaling(), QGpCoreTools::Segment::set(), QGpCoreTools::Segment2D::set(), QGpCompatibility::CompatFunction::set(), MonoStation::AbstractSummary::setBubbleValues(), ToolSPAC::setCoArrayMap(), GeopsyCore::SEGYTraceHeader::setCoordinateFactor(), ArrayCore::FKStationSignals::setCurrentShift(), ArrayCore::FKStationSignals::setPhaseShift(), SciFigs::GraphicObject::setPrintAnchor(), GeopsyCore::SEGYTraceHeader::setReceiver(), QGpGuiTools::GeographicalReference::setReference(), GeopsyCore::SEGYTraceHeader::setSource(), ArrayCore::StationCouple::setStations(), TapePositioningSystem::Triangulator::solutions(), GeopsyCore::SignalHeaderObject::sourceY(), MonoStation::StatisticResults::studentTest(), CoordReader::terminate(), QGpCoreTools::Point::toString(), utmToGeographical(), utmZoneIndex(), QGpCoreWave::TheoreticalFK::value(), ArrayCore::FKHorizontal::value(), PhaseShifter::value(), QGpCoreTools::Point::vectorialProduct(), GeopsyCore::Signal::writeSac(), GeopsyCore::Signal::writeSeg2(), SciFigs::ParallelBand::xml_setProperty(), SciFigs::ImageLayer::ReferencePoint::xml_setProperty(), SciFigs::CircleViewer::xml_writeProperties(), SciFigs::GraphContentOptions::yr2s(), and QGpCoreTools::Plane::z().
{return _y;}
| double QGpCoreTools::Point2D::y | ( | const CurvePointOptions * | ) | const [inline] |
{return _y;}