Geant4  10.02.p02
ptwXY_misc.cc
Go to the documentation of this file.
1 /*
2 # <<BEGIN-copyright>>
3 # <<END-copyright>>
4 */
5 
6 #include <stdio.h>
7 #include <stdlib.h>
8 #include <cmath>
9 #include <float.h>
10 
11 #include "ptwXY.h"
12 
13 #if defined __cplusplus
14 #include <cmath>
15 #include "G4Log.hh"
16 namespace GIDI {
17 using namespace GIDI;
18 #endif
19 
20 static nfu_status ptwXY_createFromFunctionBisect( ptwXYPoints *ptwXY, double x1, double y1, double x2, double y2, ptwXY_createFromFunction_callback func,
21  void *argList, int level, int checkForRoots, double eps );
22 static nfu_status ptwXY_createFromFunctionZeroCrossing( ptwXYPoints *ptwXY, double x1, double y1, double x2, double y2,
23  ptwXY_createFromFunction_callback func, void *argList, double eps );
24 static nfu_status ptwXY_applyFunction2( ptwXYPoints *ptwXY1, double y1, double y2, ptwXYPoint *p1, ptwXYPoint *p2, ptwXY_applyFunction_callback func,
25  void *argList, int level, int checkForRoots );
26 static nfu_status ptwXY_applyFunctionZeroCrossing( ptwXYPoints *ptwXY1, double y1, double y2, ptwXYPoint *p1, ptwXYPoint *p2,
27  ptwXY_applyFunction_callback func, void *argList );
28 /*
29 ************************************************************
30 */
31 void ptwXY_update_biSectionMax( ptwXYPoints *ptwXY1, double oldLength ) {
32 
33  ptwXY1->biSectionMax = ptwXY1->biSectionMax - 1.442695 * G4Log( ptwXY1->length / oldLength ); /* 1.442695 = 1 / std::log( 2. ) */
34  if( ptwXY1->biSectionMax < 0 ) ptwXY1->biSectionMax = 0;
35  if( ptwXY1->biSectionMax > ptwXY_maxBiSectionMax ) ptwXY1->biSectionMax = ptwXY_maxBiSectionMax;
36 }
37 /*
38 ************************************************************
39 */
40 ptwXYPoints *ptwXY_createFromFunction( int n, double *xs, ptwXY_createFromFunction_callback func, void *argList, double accuracy, int checkForRoots,
41  int biSectionMax, nfu_status *status ) {
42 
43  int64_t i;
44  double x1, y1, x2, y2, eps = ClosestAllowXFactor * DBL_EPSILON;
45  ptwXYPoints *ptwXY;
46  ptwXYPoint *p1, *p2;
47 
48  *status = nfu_Okay;
49  if( n < 2 ) { *status = nfu_tooFewPoints; return( NULL ); }
50  for( i = 1; i < n; i++ ) {
51  if( xs[i-1] >= xs[i] ) *status = nfu_XNotAscending;
52  }
53  if( *status == nfu_XNotAscending ) return( NULL );
54 
55  x1 = xs[0];
56  if( ( *status = func( x1, &y1, argList ) ) != nfu_Okay ) return( NULL );
57  if( ( ptwXY = ptwXY_new( ptwXY_interpolationLinLin, NULL, biSectionMax, accuracy, 500, 50, status, 0 ) ) == NULL ) return( NULL );
58  for( i = 1; i < n; i++ ) {
59  if( ( *status = ptwXY_setValueAtX_overrideIfClose( ptwXY, x1, y1, eps, 0 ) ) != nfu_Okay ) goto err;
60  x2 = xs[i];
61  if( ( *status = func( x2, &y2, argList ) ) != nfu_Okay ) goto err;
62  if( ( *status = ptwXY_createFromFunctionBisect( ptwXY, x1, y1, x2, y2, func, argList, 0, checkForRoots, eps ) ) != nfu_Okay ) goto err;
63  x1 = x2;
64  y1 = y2;
65  }
66  if( ( *status = ptwXY_setValueAtX_overrideIfClose( ptwXY, x2, y2, eps, 1 ) ) != nfu_Okay ) goto err;
67 
68  if( checkForRoots ) {
69  if( ( *status = ptwXY_simpleCoalescePoints( ptwXY ) ) != nfu_Okay ) goto err;
70  for( i = ptwXY->length - 1, p2 = NULL; i >= 0; i--, p2 = p1 ) { /* Work backward so lower points are still valid if a new point is added. */
71  p1 = &(ptwXY->points[i]);
72  if( p2 != NULL ) {
73  if( ( p1->y * p2->y ) < 0. ) {
74  if( ( *status = ptwXY_createFromFunctionZeroCrossing( ptwXY, p1->x, p1->y, p2->x, p2->y, func, argList, eps ) ) != nfu_Okay ) goto err;
75  }
76  }
77  }
78  }
79 
80  return( ptwXY );
81 
82 err:
83  ptwXY_free( ptwXY );
84  return( NULL );
85 }
86 /*
87 ************************************************************
88 */
89 ptwXYPoints *ptwXY_createFromFunction2( ptwXPoints *xs, ptwXY_createFromFunction_callback func, void *argList, double accuracy, int checkForRoots,
90  int biSectionMax, nfu_status *status ) {
91 
92  return( ptwXY_createFromFunction( (int) xs->length, xs->points, func, argList, accuracy, checkForRoots, biSectionMax, status ) );
93 }
94 /*
95 ************************************************************
96 */
97 static nfu_status ptwXY_createFromFunctionBisect( ptwXYPoints *ptwXY, double x1, double y1, double x2, double y2, ptwXY_createFromFunction_callback func,
98  void *argList, int level, int checkForRoots, double eps ) {
99 
100  nfu_status status;
101  double x, y, f;
102 
103  if( ( x2 - x1 ) < ClosestAllowXFactor * DBL_EPSILON * ( std::fabs( x1 ) + std::fabs( x2 ) ) ) return( nfu_Okay );
104  if( level >= ptwXY->biSectionMax ) return( nfu_Okay );
105  x = 0.5 * ( x1 + x2 );
106  if( ( status = ptwXY_interpolatePoint( ptwXY->interpolation, x, &y, x1, y1, x2, y2 ) ) != nfu_Okay ) return( status );
107  if( ( status = func( x, &f, argList ) ) != nfu_Okay ) return( status );
108  if( std::fabs( f - y ) <= 0.8 * std::fabs( f * ptwXY->accuracy ) ) return( nfu_Okay );
109  if( ( status = ptwXY_createFromFunctionBisect( ptwXY, x1, y1, x, f, func, argList, level + 1, checkForRoots, eps ) ) ) return( status );
110  if( ( status = ptwXY_setValueAtX_overrideIfClose( ptwXY, x, f, eps, 0 ) ) != nfu_Okay ) return( status );
111  return( ptwXY_createFromFunctionBisect( ptwXY, x, f, x2, y2, func, argList, level + 1, checkForRoots, eps ) );
112 }
113 /*
114 ************************************************************
115 */
116 static nfu_status ptwXY_createFromFunctionZeroCrossing( ptwXYPoints *ptwXY, double x1, double y1, double x2, double y2,
117  ptwXY_createFromFunction_callback func, void *argList, double eps ) {
118 
119  int i;
120  double x, y;
121  nfu_status status;
122 
123  for( i = 0; i < 16; i++ ) {
124  if( y2 == y1 ) break;
125  x = ( y2 * x1 - y1 * x2 ) / ( y2 - y1 );
126  if( x <= x1 ) x = x1 + 0.1 * ( x2 - x1 );
127  if( x >= x2 ) x = x2 - 0.1 * ( x2 - x1 );
128  if( ( status = func( x, &y, argList ) ) != nfu_Okay ) return( status );
129  if( y == 0 ) break;
130  if( y1 * y < 0 ) {
131  x2 = x;
132  y2 = y; }
133  else {
134  x1 = x;
135  y1 = y;
136  }
137  }
138  return( ptwXY_setValueAtX_overrideIfClose( ptwXY, x, 0., eps, 1 ) );
139 }
140 /*
141 ************************************************************
142 */
143 nfu_status ptwXY_applyFunction( ptwXYPoints *ptwXY1, ptwXY_applyFunction_callback func, void *argList, int checkForRoots ) {
144 
145  int64_t i, originalLength = ptwXY1->length, notFirstPass = 0;
146  double y1, y2 = 0;
147  nfu_status status;
148  ptwXYPoint p1, p2;
149 
150  checkForRoots = checkForRoots && ptwXY1->biSectionMax;
151  if( ptwXY1->status != nfu_Okay ) return( ptwXY1->status );
152  if( ptwXY1->interpolation == ptwXY_interpolationOther ) return( nfu_otherInterpolation );
153  if( ptwXY1->interpolation == ptwXY_interpolationFlat ) return( nfu_invalidInterpolation );
154  if( ( status = ptwXY_simpleCoalescePoints( ptwXY1 ) ) != nfu_Okay ) return( status );
155  for( i = originalLength - 1; i >= 0; i-- ) {
156  y1 = ptwXY1->points[i].y;
157  if( ( status = func( &(ptwXY1->points[i]), argList ) ) != nfu_Okay ) return( status );
158  p1 = ptwXY1->points[i];
159  if( notFirstPass ) {
160  if( ( status = ptwXY_applyFunction2( ptwXY1, y1, y2, &p1, &p2, func, argList, 0, checkForRoots ) ) != nfu_Okay ) return( status );
161  }
162  notFirstPass = 1;
163  p2 = p1;
164  y2 = y1;
165  }
166  ptwXY_update_biSectionMax( ptwXY1, (double) originalLength );
167  return( status );
168 }
169 /*
170 ************************************************************
171 */
172 static nfu_status ptwXY_applyFunction2( ptwXYPoints *ptwXY1, double y1, double y2, ptwXYPoint *p1, ptwXYPoint *p2, ptwXY_applyFunction_callback func,
173  void *argList, int level, int checkForRoots ) {
174 
175  double y;
176  ptwXYPoint p;
177  nfu_status status;
178 
179  if( ( p2->x - p1->x ) < ClosestAllowXFactor * DBL_EPSILON * ( std::fabs( p1->x ) + std::fabs( p2->x ) ) ) return( nfu_Okay );
180  if( level >= ptwXY1->biSectionMax ) goto checkForZeroCrossing;
181  p.x = 0.5 * ( p1->x + p2->x );
182  if( ( status = ptwXY_interpolatePoint( ptwXY1->interpolation, p.x, &y, p1->x, y1, p2->x, y2 ) ) != nfu_Okay ) return( status );
183  p.y = y;
184  if( ( status = func( &p, argList ) ) != nfu_Okay ) return( status );
185  if( std::fabs( ( p.x - p1->x ) * ( p2->y - p1->y ) + ( p2->x - p1->x ) * ( p1->y - p.y ) ) <= 0.8 * std::fabs( ( p2->x - p1->x ) * p.y * ptwXY1->accuracy ) )
186  goto checkForZeroCrossing;
187  if( ( status = ptwXY_setValueAtX( ptwXY1, p.x, p.y ) ) != nfu_Okay ) return( status );
188  if( ( status = ptwXY_applyFunction2( ptwXY1, y1, y, p1, &p, func, argList, level + 1, checkForRoots ) ) ) return( status );
189  return( ptwXY_applyFunction2( ptwXY1, y, y2, &p, p2, func, argList, level + 1, checkForRoots ) );
190 
191 checkForZeroCrossing:
192  if( checkForRoots && ( ( p1->y * p2->y ) < 0. ) ) return( ptwXY_applyFunctionZeroCrossing( ptwXY1, y1, y2, p1, p2, func, argList ) );
193  return( nfu_Okay );
194 }
195 /*
196 ************************************************************
197 */
198 static nfu_status ptwXY_applyFunctionZeroCrossing( ptwXYPoints *ptwXY1, double y1, double y2, ptwXYPoint *p1, ptwXYPoint *p2,
199  ptwXY_applyFunction_callback func, void *argList ) {
200 
201  int i;
202  double y, x1 = p1->x, x2 = p2->x, nY1 = p1->y, nY2 = p2->y, refY = 0.5 * ( std::fabs( p1->y ) + std::fabs( p2->y ) );
203  ptwXYPoint p;
204  nfu_status status;
205 
206  for( i = 0; i < 6; i++ ) {
207  if( nY2 == nY1 ) break;
208  p.x = ( nY2 * x1 - nY1 * x2 ) / ( nY2 - nY1 );
209  if( p.x <= x1 ) p.x = 0.5 * ( x1 + x2 );
210  if( p.x >= x2 ) p.x = 0.5 * ( x1 + x2 );
211  if( ( status = ptwXY_interpolatePoint( ptwXY1->interpolation, p.x, &y, p1->x, y1, p2->x, y2 ) ) != nfu_Okay ) return( status );
212  p.y = y;
213  if( ( status = func( &p, argList ) ) != nfu_Okay ) return( status );
214  if( p.y == 0 ) break;
215  if( 0.5 * refY < std::fabs( p.y ) ) break;
216  refY = std::fabs( p.y );
217  if( p1->y * p.y < 0 ) {
218  x2 = p.x;
219  nY2 = p.y; }
220  else {
221  x1 = p.x;
222  nY1 = p.y;
223  }
224  }
225  return( ptwXY_setValueAtX( ptwXY1, p.x, 0. ) );
226 }
227 /*
228 ************************************************************
229 */
230 ptwXYPoints *ptwXY_fromString( char const *str, ptwXY_interpolation interpolation, ptwXY_interpolationOtherInfo const *interpolationOtherInfo,
231  double biSectionMax, double accuracy, char **endCharacter, nfu_status *status ) {
232 
233  int64_t numberConverted;
234  double *doublePtr;
235  ptwXYPoints *ptwXY = NULL;
236 
237  if( ( *status = nfu_stringToListOfDoubles( str, &numberConverted, &doublePtr, endCharacter ) ) != nfu_Okay ) return( NULL );
238  *status = nfu_oddNumberOfValues;
239  if( ( numberConverted % 2 ) == 0 )
240  ptwXY = ptwXY_create( interpolation, interpolationOtherInfo, biSectionMax, accuracy, numberConverted, 10, numberConverted / 2, doublePtr, status, 0 );
241  nfu_free( doublePtr );
242  return( ptwXY );
243 }
244 /*
245 ************************************************************
246 */
247 void ptwXY_showInteralStructure( ptwXYPoints *ptwXY, FILE *f, int printPointersAsNull ) {
248 
249  int64_t i, n = ptwXY_getNonOverflowLength( ptwXY );
250  ptwXYPoint *point = ptwXY->points;
251  ptwXYOverflowPoint *overflowPoint;
252 
253  fprintf( f, "status = %d interpolation = %d length = %d allocatedSize = %d\n",
254  (int) ptwXY->status, (int) ptwXY->interpolation, (int) ptwXY->length, (int) ptwXY->allocatedSize );
255  fprintf( f, "userFlag = %d biSectionMax = %.8e accuracy = %.2e minFractional_dx = %.6e\n",
256  ptwXY->userFlag, ptwXY->biSectionMax, ptwXY->accuracy, ptwXY->minFractional_dx );
257  fprintf( f, "interpolationString = %s\n", ptwXY->interpolationOtherInfo.interpolationString );
258  fprintf( f, "getValueFunc is NULL = %d. argList is NULL = %d.\n",
259  ptwXY->interpolationOtherInfo.getValueFunc == NULL, ptwXY->interpolationOtherInfo.argList == NULL );
260  fprintf( f, " overflowLength = %d overflowAllocatedSize = %d mallocFailedSize = %d\n",
261  (int) ptwXY->overflowLength, (int) ptwXY->overflowAllocatedSize, (int) ptwXY->mallocFailedSize );
262  fprintf( f, " Points data, points = %20p\n", ( printPointersAsNull ? NULL : (void*)ptwXY->points ) );
263  for( i = 0; i < n; i++, point++ ) fprintf( f, " %14.7e %14.7e\n", point->x, point->y );
264  fprintf( f, " Overflow points data; %20p\n", ( printPointersAsNull ? NULL : (void*)&(ptwXY->overflowHeader) ) );
265  for( overflowPoint = ptwXY->overflowHeader.next; overflowPoint != &(ptwXY->overflowHeader); overflowPoint = overflowPoint->next ) {
266  fprintf( f, " %14.7e %14.7e %8d %20p %20p %20p\n", overflowPoint->point.x, overflowPoint->point.y, (int) overflowPoint->index,
267  (void*) ( printPointersAsNull ? NULL : overflowPoint ), (void*) ( printPointersAsNull ? NULL : overflowPoint->prior ),
268  (void*) ( printPointersAsNull ? NULL : overflowPoint->next ) );
269  }
270  fprintf( f, " Points in order\n" );
271  for( i = 0; i < ptwXY->length; i++ ) {
272  point = ptwXY_getPointAtIndex( ptwXY, i );
273  fprintf( f, " %14.7e %14.7e\n", point->x, point->y );
274  }
275 }
276 /*
277 ************************************************************
278 */
279 void ptwXY_simpleWrite( ptwXYPoints *ptwXY, FILE *f, char *format ) {
280 
281  int64_t i;
282  ptwXYPoint *point;
283 
284  for( i = 0; i < ptwXY->length; i++ ) {
285  point = ptwXY_getPointAtIndex( ptwXY, i );
286  fprintf( f, format, point->x, point->y );
287  }
288 }
289 /*
290 ************************************************************
291 */
292 void ptwXY_simplePrint( ptwXYPoints *ptwXY, char *format ) {
293 
294  ptwXY_simpleWrite( ptwXY, stdout, format );
295 }
296 
297 #if defined __cplusplus
298 }
299 #endif
void ptwXY_simpleWrite(ptwXYPoints *ptwXY, FILE *f, char *format)
Definition: ptwXY_misc.cc:279
static nfu_status ptwXY_applyFunction2(ptwXYPoints *ptwXY1, double y1, double y2, ptwXYPoint *p1, ptwXYPoint *p2, ptwXY_applyFunction_callback func, void *argList, int level, int checkForRoots)
Definition: ptwXY_misc.cc:172
ptwXYPoints * ptwXY_createFromFunction2(ptwXPoints *xs, ptwXY_createFromFunction_callback func, void *argList, double accuracy, int checkForRoots, int biSectionMax, nfu_status *status)
Definition: ptwXY_misc.cc:89
static nfu_status ptwXY_createFromFunctionZeroCrossing(ptwXYPoints *ptwXY, double x1, double y1, double x2, double y2, ptwXY_createFromFunction_callback func, void *argList, double eps)
Definition: ptwXY_misc.cc:116
static nfu_status ptwXY_createFromFunctionBisect(ptwXYPoints *ptwXY, double x1, double y1, double x2, double y2, ptwXY_createFromFunction_callback func, void *argList, int level, int checkForRoots, double eps)
Definition: ptwXY_misc.cc:97
ptwXYPoints * ptwXY_create(ptwXY_interpolation interpolation, ptwXY_interpolationOtherInfo const *interpolationOtherInfo, double biSectionMax, double accuracy, int64_t primarySize, int64_t secondarySize, int64_t length, double const *xy, nfu_status *status, int userFlag)
Definition: ptwXY_core.cc:108
int64_t ptwXY_getNonOverflowLength(ptwXYPoints const *ptwXY)
Definition: ptwXY_core.cc:590
nfu_status ptwXY_simpleCoalescePoints(ptwXYPoints *ptwXY)
Definition: ptwXY_core.cc:529
ptwXYPoints * ptwXY_free(ptwXYPoints *ptwXY)
Definition: ptwXY_core.cc:574
static const G4double eps
nfu_status ptwXY_setValueAtX(ptwXYPoints *ptwXY, double x, double y)
Definition: ptwXY_core.cc:876
nfu_status ptwXY_applyFunction(ptwXYPoints *ptwXY1, ptwXY_applyFunction_callback func, void *argList, int checkForRoots)
Definition: ptwXY_misc.cc:143
nfu_status nfu_stringToListOfDoubles(char const *str, int64_t *numberConverted, double **doublePtr, char **endCharacter)
static nfu_status ptwXY_applyFunctionZeroCrossing(ptwXYPoints *ptwXY1, double y1, double y2, ptwXYPoint *p1, ptwXYPoint *p2, ptwXY_applyFunction_callback func, void *argList)
Definition: ptwXY_misc.cc:198
nfu_status ptwXY_setValueAtX_overrideIfClose(ptwXYPoints *ptwXY, double x, double y, double eps, int override)
Definition: ptwXY_core.cc:883
#define DBL_EPSILON
Definition: templates.hh:87
const G4int n
G4double G4Log(G4double x)
Definition: G4Log.hh:230
void ptwXY_update_biSectionMax(ptwXYPoints *ptwXY1, double oldLength)
Definition: ptwXY_misc.cc:31
const G4double x[NPOINTSGL]
void ptwXY_showInteralStructure(ptwXYPoints *ptwXY, FILE *f, int printPointersAsNull)
Definition: ptwXY_misc.cc:247
ptwXYPoints * ptwXY_createFromFunction(int n, double *xs, ptwXY_createFromFunction_callback func, void *argList, double accuracy, int checkForRoots, int biSectionMax, nfu_status *status)
Definition: ptwXY_misc.cc:40
ptwXYPoint * ptwXY_getPointAtIndex(ptwXYPoints *ptwXY, int64_t index)
Definition: ptwXY_core.cc:675
void ptwXY_simplePrint(ptwXYPoints *ptwXY, char *format)
Definition: ptwXY_misc.cc:292
void * nfu_free(void *p)
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)
Definition: ptwXY_core.cc:29
nfu_status ptwXY_interpolatePoint(ptwXY_interpolation interpolation, double x, double *y, double x1, double y1, double x2, double y2)
ptwXYPoints * ptwXY_fromString(char const *str, ptwXY_interpolation interpolation, ptwXY_interpolationOtherInfo const *interpolationOtherInfo, double biSectionMax, double accuracy, char **endCharacter, nfu_status *status)
Definition: ptwXY_misc.cc:230