11#if defined __cplusplus
21static double ptwXY_flatInterpolationToLinear_eps(
double px,
double eps );
23static nfu_status ptwXY_LogLogToLinLin(
ptwXYPoints *desc,
double x1,
double y1,
double x2,
double y2,
int depth );
24static nfu_status ptwXY_LinLogToLinLin(
ptwXYPoints *desc,
double x1,
double y1,
double x2,
double y2,
int depth );
25static nfu_status ptwXY_LogLinToLinLin(
ptwXYPoints *desc,
double x1,
double y1,
double x2,
double y2,
int depth );
26static nfu_status ptwXY_otherToLinLin(
ptwXYPoints *desc,
double x1,
double y1,
double x2,
double y2,
int depth );
39 *y = 0.5 * ( y1 + y2 ); }
45 switch( interpolation ) {
47 *y = ( y1 * ( x2 - x ) + y2 * ( x - x1 ) ) / ( x2 - x1 );
51 *y = ( y1 *
G4Log( x2 / x ) + y2 *
G4Log( x / x1 ) ) /
G4Log( x2 / x1 );
55 *y =
G4Exp( (
G4Log( y1 ) * ( x2 - x ) +
G4Log( y2 ) * ( x - x1 ) ) / ( x2 - x1 ) );
87 if( ( lowerEps < 0 ) || ( upperEps < 0 ) || ( ( lowerEps == 0 ) && ( upperEps == 0 ) ) )
return( NULL );
88 if( ( lowerEps != 0 ) && ( lowerEps <
minEps ) ) lowerEps =
minEps;
89 if( ( upperEps != 0 ) && ( upperEps <
minEps ) ) upperEps =
minEps;
91 length = ptwXY->
length * ( 1 + ( lowerEps == 0 ? 0 : 1 ) + ( lowerEps == 0 ? 0 : 1 ) );
96 for( i = 0; i < ptwXY->
length; i++, p3++ ) {
99 x = ptwXY_flatInterpolationToLinear_eps( p2->x, -lowerEps );
107 x = ptwXY_flatInterpolationToLinear_eps( p2->x, upperEps );
117 if( ( lowerEps != 0 ) && ( p1->
y != p2->y ) ) {
118 x = ptwXY_flatInterpolationToLinear_eps( p2->x, -lowerEps );
137static double ptwXY_flatInterpolationToLinear_eps(
double px,
double eps ) {
142 x = ( 1 - eps ) * px; }
144 x = ( 1 + eps ) * px; }
168 func = ptwXY_LogLogToLinLin;
break;
170 func = ptwXY_LinLogToLinLin;
break;
172 func = ptwXY_LogLinToLinLin;
break;
183 if( func == NULL )
return( NULL );
187 if( accuracy < ptwXY->accuracy ) accuracy = ptwXY->
accuracy;
192 *status = ptwXY_toOtherInterpolation2( n1, ptwXY, func );
205 double x1, y1, x2, y2;
211 for( i = 1; i < src->
length; i++ ) {
214 if( ( x1 != x2 ) && ( y1 != y2 ) ) {
215 if( ( status = func( desc, x1, y1, x2, y2, 0 ) ) !=
nfu_Okay )
break;
225static nfu_status ptwXY_LogLogToLinLin(
ptwXYPoints *desc,
double x1,
double y1,
double x2,
double y2,
int depth ) {
228 double x, y, u, u2 = x2 / x1, v2 = y2 / y1, logYXs, logXs =
G4Log( u2 ), logYs =
G4Log( v2 ), vLin, vLog, w;
230 logYXs = logYs / logXs;
232 if( depth > 16 )
return(
nfu_Okay );
233 if( std::fabs( logYXs - 1 ) < 1e-5 ) {
234 u = 0.5 * ( 1 + u2 );
235 w = ( logYXs - 1 ) * logXs;
236 vLog = u * ( 1. + w * ( 1 + 0.5 * w ) ); }
238 u = logYXs * ( u2 - v2 ) / ( ( 1 - logYXs ) * ( v2 - 1 ) );
241 vLin = ( u2 - u + v2 * ( u - 1 ) ) / ( u2 - 1 );
242 if( std::fabs( vLog - vLin ) <= ( vLog * desc->
accuracy ) )
return( status );
246 if( ( status = ptwXY_LogLogToLinLin( desc, x1, y1, x, y, depth + 1 ) ) !=
nfu_Okay )
return( status );
247 return( ptwXY_LogLogToLinLin( desc, x, y, x2, y2, depth + 1 ) );
252static nfu_status ptwXY_LinLogToLinLin(
ptwXYPoints *desc,
double x1,
double y1,
double x2,
double y2,
int depth ) {
255 double x, y, logYs =
G4Log( y2 / y1 ), yLinLin;
257 if( depth > 16 )
return(
nfu_Okay );
258 x = ( x2 - x1 ) / ( y2 - y1 ) * ( ( y2 - y1 ) / logYs - y1 ) + x1;
259 y = y1 *
G4Exp( logYs / ( x2 - x1 ) * ( x - x1 ) );
260 yLinLin = ( y1 * ( x2 - x ) + y2 * ( x - x1 ) ) / ( x2 - x1 );
261 if( std::fabs( y - yLinLin ) <= ( y * desc->
accuracy ) )
return( status );
263 if( ( status = ptwXY_LinLogToLinLin( desc, x1, y1, x, y, depth + 1 ) ) !=
nfu_Okay )
return( status );
264 return( ptwXY_LinLogToLinLin( desc, x, y, x2, y2, depth + 1 ) );
269static nfu_status ptwXY_LogLinToLinLin(
ptwXYPoints *desc,
double x1,
double y1,
double x2,
double y2,
int depth ) {
272 double x = std::sqrt( x2 * x1 ), y, logXs =
G4Log( x2 / x1 ), yLinLin;
274 if( depth > 16 )
return(
nfu_Okay );
276 x = ( y1 * x2 - y2 * x1 ) / ( y1 * logXs + ( y2 - y1 ) * ( std::log( x / x1 ) - 1 ) );
278 y = ( y2 - y1 ) *
G4Log( x / x1 ) / logXs + y1;
279 yLinLin = ( y1 * ( x2 - x ) + y2 * ( x - x1 ) ) / ( x2 - x1 );
280 if( std::fabs( y - yLinLin ) <= ( y * desc->
accuracy ) )
return( status );
282 if( ( status = ptwXY_LogLinToLinLin( desc, x1, y1, x, y, depth + 1 ) ) !=
nfu_Okay )
return( status );
283 return( ptwXY_LogLinToLinLin( desc, x, y, x2, y2, depth + 1 ) );
288static nfu_status ptwXY_otherToLinLin(
ptwXYPoints *desc,
double x1,
double y1,
double x2,
double y2,
int depth ) {
291 double x = 0.5 * ( x1 + x2 ), y, yLinLin;
295 if( depth > 16 )
return(
nfu_Okay );
296 if( ( status = getValueFunc( argList, x, &y, x1, y1, x2, y2 ) ) !=
nfu_Okay )
return( status );
297 yLinLin = ( y1 * ( x2 - x ) + y2 * ( x - x1 ) ) / ( x2 - x1 );
298 if( std::fabs( y - yLinLin ) <= ( y * desc->
accuracy ) )
return( status );
300 if( ( status = ptwXY_otherToLinLin( desc, x1, y1, x, y, depth + 1 ) ) !=
nfu_Okay )
return( status );
301 return( ptwXY_otherToLinLin( desc, x, y, x2, y2, depth + 1 ) );
311 double xMin, xMax, dx, inverseDx;
314 if( ptwXY->
length < 2 )
return( NULL );
315 if( ( n =
ptwXY_clone( ptwXY, status ) ) == NULL )
return( NULL );
317 xMin = n->points[0].x;
318 xMax = n->points[n->length-1].x;
321 for( i = 0, p = n->points; i < n->length; i++, p++ ) {
322 p->
x = ( p->
x - xMin ) * inverseDx;
325 n->points[n->length-1].x = 1.;
336 double dx, inverseDx, xLast = 0.;
339 if( ptwXY->
length < 2 )
return( NULL );
340 if( ( n =
ptwXY_clone( ptwXY, status ) ) == NULL )
return( NULL );
345 for( i = 0, p2 = p = n->points; i < length; ++i, ++p ) {
346 p2->
x = p->
x * dx + xMin;
348 if( std::fabs( p2->
x - xLast ) <= 10. *
DBL_EPSILON * ( std::fabs( p2->
x ) + std::fabs( xLast ) ) ) {
353 p2->
y = p->
y * inverseDx;
357 n->points[n->length-1].x = xMax;
368 ptwXYPoints *n1 = NULL, *n2 = NULL, *a = NULL, *r = NULL;
370 double f, g, xMin, xMax;
374 if( w < w1 )
return( NULL );
378 if( w > w2 )
return( NULL );
383 f = ( w - w1 ) / ( w2 - w1 );
385 for( i = 0, p = n1->
points; i < n1->length; i++, p++ ) p->
y *= g;
386 for( i = 0, p = n2->points; i < n2->length; i++, p++ ) p->
y *= f;
404#if defined __cplusplus
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
G4double G4Log(G4double x)
static G4Pow * GetInstance()
G4double powA(G4double A, G4double y) const
@ nfu_invalidInterpolation
@ nfu_unsupportedInterpolationConversion
enum nfu_status_e nfu_status
nfu_status ptwXY_setValueAtX(ptwXYPoints *ptwXY, double x, double y)
ptwXYPoints * ptwXY_toUnitbase(ptwXYPoints *ptwXY, nfu_status *status)
ptwXYPoints * ptwXY_unitbaseInterpolate(double w, double w1, ptwXYPoints *ptwXY1, double w2, ptwXYPoints *ptwXY2, nfu_status *status)
nfu_status ptwXY_simpleCoalescePoints(ptwXYPoints *ptwXY)
ptwXYPoints * ptwXY_new(ptwXY_interpolation interpolation, ptwXY_interpolationOtherInfo const *interpolationOtherInfo, double biSectionMax, double accuracy, int64_t primarySize, int64_t secondarySize, nfu_status *status, int userFlag)
enum ptwXY_interpolation_e ptwXY_interpolation
@ ptwXY_interpolationFlat
@ ptwXY_interpolationLinLog
@ ptwXY_interpolationLogLog
@ ptwXY_interpolationLinLin
@ ptwXY_interpolationOther
@ ptwXY_interpolationLogLin
nfu_status(* ptwXY_getValue_callback)(void *argList, double x, double *y, double x1, double y1, double x2, double y2)
ptwXYPoints * ptwXY_fromUnitbase(ptwXYPoints *ptwXY, double xMin, double xMax, nfu_status *status)
ptwXYPoints * ptwXY_cloneToInterpolation(ptwXYPoints *ptwXY, ptwXY_interpolation interpolationTo, nfu_status *status)
ptwXYPoints * ptwXY_free(ptwXYPoints *ptwXY)
ptwXYPoints * ptwXY_flatInterpolationToLinear(ptwXYPoints *ptwXY, double lowerEps, double upperEps, nfu_status *status)
ptwXYPoints * ptwXY_add_ptwXY(ptwXYPoints *ptwXY1, ptwXYPoints *ptwXY2, nfu_status *status)
ptwXYPoints * ptwXY_clone(ptwXYPoints *ptwXY, nfu_status *status)
nfu_status ptwXY_interpolatePoint(ptwXY_interpolation interpolation, double x, double *y, double x1, double y1, double x2, double y2)
ptwXYPoints * ptwXY_toOtherInterpolation(ptwXYPoints *ptwXY, ptwXY_interpolation interpolation, double accuracy, nfu_status *status)
nfu_status(* interpolation_func)(ptwXYPoints *desc, double x1, double y1, double x2, double y2, int depth)
ptwXY_interpolation interpolation
ptwXY_interpolationOtherInfo interpolationOtherInfo
ptwXY_getValue_callback getValueFunc