00001 /* 00002 ------------------------------------------------------------------- 00003 00004 Copyright (C) 2006, 2007, 2008, Andrew W. Steiner 00005 00006 This file is part of O2scl. 00007 00008 O2scl is free software; you can redistribute it and/or modify 00009 it under the terms of the GNU General Public License as published by 00010 the Free Software Foundation; either version 3 of the License, or 00011 (at your option) any later version. 00012 00013 O2scl is distributed in the hope that it will be useful, 00014 but WITHOUT ANY WARRANTY; without even the implied warranty of 00015 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the 00016 GNU General Public License for more details. 00017 00018 You should have received a copy of the GNU General Public License 00019 along with O2scl. If not, see <http://www.gnu.org/licenses/>. 00020 00021 ------------------------------------------------------------------- 00022 */ 00023 #ifndef O2SCL_CONTOUR_H 00024 #define O2SCL_CONTOUR_H 00025 00026 #include <gsl/gsl_math.h> 00027 #include <o2scl/smart_interp.h> 00028 #include <o2scl/omatrix_tlate.h> 00029 #include <o2scl/pinside.h> 00030 00031 #ifndef DOXYGENP 00032 namespace o2scl { 00033 #endif 00034 00035 /** 00036 \brief Calculate contour lines from a two-dimensional data set 00037 00038 \b Basic \b Usage 00039 00040 - Specify the data as a two-dimensional square grid of "heights" 00041 with set_data(). 00042 - Specify the contour levels with set_levels(). 00043 - Compute the contours with calc_contours() which returns 00044 the number of contours (which is often larger than the number 00045 of contour levels, since one level can result in multiple contours). 00046 - Retrieve the contours individually using calls to 00047 get_contour(). 00048 00049 The data should be stored so that the y-index is first, i.e. 00050 data[iy][ix]. One can always switch \c x_fun and \c y_fun if 00051 this is not the case. The data is copied by set_data(), so 00052 changing the data will not change the contours unless set_data() 00053 is called again. The functions set_levels() and calc() can be 00054 called several times for the same data without calling 00055 set_data() again. 00056 00057 Linear interpolation is used to decide whether or not a line 00058 segment and a contour cross. This choice is intentional, since 00059 (in addition to making the algorithm much simpler) 00060 it is the user (and not the class) which is likely best able to 00061 refine the data. In case a simple refinement scheme is desired, 00062 the method regrid_data() is provided which uses cubic spline 00063 interpolation to refine the data and thus make the curves more 00064 continuous. 00065 00066 Since linear interpolation is used, the contour calculation 00067 implicitly assumes that there is not more than one intersection 00068 of any contour level with any line segment, so if this is the 00069 case, then either a more refined data set should be specified or 00070 regrid_data() should be used. For contours which do not close 00071 inside the region of interest, the results will always end at 00072 either the minimum or maximum values of \c x_fun or \c y_fun (no 00073 extrapolation is ever done). 00074 00075 As an example, for the function 00076 \f[ 00077 15 e^{-(x-20)^2/400-(y-5)^2/25}+ 00078 40 e^{-(x-70)^2/4900-(y-2)^2/4} 00079 \f] 00080 a 10x10 grid gives the contours: 00081 \image html contourg.png "contourg.png" 00082 \image latex contourg.eps "contourg.eps" width=9cm 00083 00084 While after a call to regrid_data(3,3), the contours 00085 are a little smoother: 00086 \image html contourg2.png "contourg2.png" 00087 \image latex contourg2.eps "contourg2.eps" width=9cm 00088 00089 Mathematica gives a similar result: 00090 \image html contourg3.png "contourg3.png" 00091 \image latex contourg3.eps "contourg3.eps" width=9cm 00092 00093 \b The \b Algorithm: 00094 00095 This works by viewing the data as defining a square 00096 two-dimensional grid. The function calc() exhaustively 00097 enumerates every line segment in the grid which involves a level 00098 crossing and then organizing the points defined by the 00099 intersection of a line segment with a level curve into a full 00100 contour. 00101 00102 \b Representing \b Contours \b by \b Fill \b Regions 00103 00104 If the user wants to "shade" the contours to provide a 00105 shaded or colored contour plot, then it is useful 00106 to provide a set of closed contours to be shaded instead 00107 of the (possibly) open contours returned by calc_contours(). 00108 After a call to calc_contours(), the function regions() can 00109 be used to organize the closed contours into regions to 00110 be shaded. Open contours are closed by adding points defining 00111 lines along the edges of the data and closed contours are 00112 inverted (if necessary) by adding an external line emanating 00113 from the closed contour and properly including the boundary 00114 edges. 00115 00116 \todo 00117 - Some contours which should be closed are not properly 00118 closed yet. See the tests for examples which fail. 00119 - Use twod_intp for regrid_data 00120 - Work on how memory is allocated 00121 - Include the functionality to provide regions to be shaded 00122 instead of just lines. This can be done by providing a method, 00123 i.e. regions() which converts the curves into regions. 00124 00125 */ 00126 class contour { 00127 00128 public: 00129 00130 contour(); 00131 ~contour(); 00132 00133 /// \name Basic usage 00134 //@{ 00135 00136 /** 00137 \brief Set the data 00138 00139 The types \c vec_t and \c mat_t can be any types which have \c 00140 operator[] and \c operator[][] for array and matrix indexing. 00141 00142 Note that this method copies all of the user-specified data to 00143 local storage so that changes in the data after calling this function 00144 will not be reflected in the contours that are generated. 00145 */ 00146 template<class vec_t, class mat_t> 00147 int set_data(size_t sizex, size_t sizey, const vec_t &x_fun, 00148 const vec_t &y_fun, const mat_t &udata) { 00149 int i; 00150 00151 if (sizex<2 || sizey<2) { 00152 set_err_ret("Not enough data (must be at least 2x2) in set_data().", 00153 gsl_einval); 00154 } 00155 00156 free_memory(); 00157 00158 nx=sizex; 00159 ny=sizey; 00160 xfun=new ovector(nx); 00161 yfun=new ovector(ny); 00162 for(int is=0;is<nx;is++) (*xfun)[is]=x_fun[is]; 00163 for(int is=0;is<ny;is++) (*yfun)[is]=y_fun[is]; 00164 00165 data=new ovector *[ny]; 00166 for(i=0;i<ny;i++) { 00167 data[i]=new ovector(nx); 00168 for(int j=0;j<nx;j++) { 00169 (*(data[i]))[j]=udata[i][j]; 00170 } 00171 } 00172 00173 return 0; 00174 } 00175 00176 /** 00177 \brief Set the contour levels 00178 00179 This is separate from the function calc_contours() so that 00180 the user can compute the contours for different data sets using 00181 the same levels 00182 */ 00183 template<class vec_t> int set_levels(size_t nlevels, vec_t &ulevels) { 00184 nlev=nlevels; 00185 levels.allocate(nlevels); 00186 for(size_t i=0;i<nlevels;i++) { 00187 levels[i]=ulevels[i]; 00188 } 00189 levels_set=true; 00190 return 0; 00191 } 00192 00193 /** 00194 \brief Calculate the contours 00195 00196 The function calc_contours() returns the total number of 00197 contours found. Since there may be more than one disconnected 00198 contours for the same contour level, or no contours for a 00199 given level, the total number of contours may be less than or 00200 greater than the number of levels given by set_levels(). 00201 00202 If an error occurs, zero is returned. 00203 */ 00204 template<class vec_t> 00205 int calc_contours(vec_t &new_levels, bool debug=false) { 00206 int i,j,k,ileft,iright,m; 00207 00208 for(i=0;i<nlev;i++) new_levels[i]=levels[i]; 00209 00210 if (nx==0) { 00211 set_err("Data not set in calc_contours().",gsl_einval); 00212 return 0; 00213 } 00214 if (levels_set==false) { 00215 set_err("Contour levels not set in calc_contours().",gsl_einval); 00216 return 0; 00217 } 00218 00219 // Free previously stored contours 00220 if (ncurves!=0) { 00221 csizes.clear(); 00222 vals.clear(); 00223 xc.clear(); 00224 yc.clear(); 00225 ncurves=0; 00226 } 00227 00228 // Free space for previously computed edges 00229 for(j=0;j<((int)(re.size()));j++) { 00230 delete re[j]; 00231 delete be[j]; 00232 delete redges[j]; 00233 delete bedges[j]; 00234 } 00235 re.clear(); 00236 be.clear(); 00237 redges.clear(); 00238 bedges.clear(); 00239 00240 // Temporary storage for the interpolation arrays 00241 int nint; 00242 double *xi, *yi; 00243 00244 // Some useful variables for looping over the lines 00245 bool foundline, linedone; 00246 int fp; 00247 00248 // Temporary storage for the curves 00249 std::vector<double> x, y; 00250 00251 // The interpolation object 00252 00253 linear_interp<ovector_view> it1; 00254 linear_interp<ovector_const_reverse> it2; 00255 linear_interp<ovector_const_subvector> it3; 00256 linear_interp<ovector_const_subvector_reverse> it4; 00257 sm_interp *si=new sm_interp(it1,it2,it3,it4); 00258 00259 // Start with zero curves 00260 ncurves=0; 00261 00262 // For each level 00263 for(i=0;i<nlev;i++) { 00264 00265 // Make space for the edges and add to the list 00266 rep=new omatrix_int(ny-1,nx); 00267 bep=new omatrix_int(ny,nx-1); 00268 redgesp=new omatrix(ny-1,nx); 00269 bedgesp=new omatrix(ny,nx-1); 00270 re.push_back(rep); 00271 be.push_back(bep); 00272 redges.push_back(redgesp); 00273 bedges.push_back(bedgesp); 00274 00275 if (verbose>0) { 00276 std::cout << "\nLooking for edges for level: " 00277 << new_levels[i] << std::endl; 00278 } 00279 00280 // Examine the each of the rows for an intersection 00281 find_intersections(i,new_levels[i]); 00282 00283 if (verbose>0) { 00284 std::cout << "\nInterpolating edge intersections for level: " 00285 << new_levels[i] << std::endl; 00286 } 00287 00288 // Process the right edges 00289 right_edges(new_levels[i],si); 00290 00291 // Process the bottom edges 00292 bottom_edges(new_levels[i],si); 00293 00294 if (verbose>0) { 00295 std::cout << "\nPiecing together contour lines for level: " 00296 << new_levels[i] << std::endl; 00297 } 00298 00299 // Now go through and one side of the line 00300 foundline=true; 00301 while(foundline==true) { 00302 for(j=0;j<nx;j++) { 00303 for(k=0;k<ny;k++) { 00304 foundline=false; 00305 00306 // A line beginning with a right edge 00307 if (k<ny-1 && (*rep)[k][j]==edge) { 00308 if (verbose>0) { 00309 std::cout << "Starting contour line: " << std::endl; 00310 std::cout << "(" << (*xfun)[j] << ", " << (*redgesp)[k][j] 00311 << ")" << std::endl; 00312 } 00313 x.push_back((*xfun)[j]); 00314 y.push_back((*redgesp)[k][j]); 00315 (*rep)[k][j]++; 00316 00317 // Go through both sides 00318 process_line(j,k,dright,x,y,true); 00319 if (verbose>0) { 00320 std::cout << "Computing other side of line." << std::endl; 00321 } 00322 process_line(j,k,dright,x,y,false); 00323 foundline=true; 00324 } 00325 00326 // A line beginning with a bottom edge 00327 if (j<nx-1 && foundline==false && (*bep)[k][j]==edge) { 00328 if (verbose>0) { 00329 std::cout << "Starting contour line: " << std::endl; 00330 std::cout << "(" << (*bedgesp)[k][j] << ", " << (*yfun)[k] 00331 << ")" << std::endl; 00332 } 00333 x.push_back((*bedgesp)[k][j]); 00334 y.push_back((*yfun)[k]); 00335 (*bep)[k][j]++; 00336 00337 // Go through both sides 00338 process_line(j,k,dbottom,x,y,true); 00339 if (verbose>0) { 00340 std::cout << "Computing other side of line." << std::endl; 00341 } 00342 process_line(j,k,dbottom,x,y,false); 00343 foundline=true; 00344 } 00345 00346 // Now repackage the line information in x and y 00347 // to xc and yc and update csizes, vals and ncurves: 00348 if (foundline==true) { 00349 csizes.push_back(x.size()); 00350 vals.push_back(new_levels[i]); 00351 ncurves++; 00352 std::vector<double> xv, yv; 00353 for(m=0;m<((int)x.size());m++) { 00354 xv.push_back(x[m]); 00355 yv.push_back(y[m]); 00356 } 00357 xc.push_back(xv); 00358 yc.push_back(yv); 00359 x.clear(); 00360 y.clear(); 00361 } 00362 } 00363 } 00364 } 00365 } 00366 00367 return ncurves; 00368 } 00369 00370 /// Return the number of points in the specified contour 00371 int get_contour_size(int index); 00372 00373 /** 00374 \brief Get a contour 00375 00376 Given the \c index, which is between 0 and the number of contours 00377 returned by calc_contours() minus 1 (inclusive), this 00378 returns the level for this contour with the x and y-values 00379 in \c x and \c y which are both of length \c csize. 00380 The vectors \c x and \c y must have been previously allocated. 00381 The user can obtain the necessary size for \c x and \c y by 00382 calling get_contour_size(). 00383 */ 00384 template<class vec_t> int get_contour(int index, double &val, int &csize, 00385 vec_t &x, vec_t &y) { 00386 int i; 00387 if (index>=ncurves) { 00388 set_err_ret("No curve with specified index in get_contour().", 00389 gsl_efailed); 00390 } 00391 00392 csize=csizes[index]; 00393 val=vals[index]; 00394 for(i=0;i<csize;i++) { 00395 x[i]=xc[index][i]; 00396 y[i]=yc[index][i]; 00397 } 00398 return 0; 00399 } 00400 00401 /** 00402 \brief Compute closed contour regions (unfinished) 00403 */ 00404 int regions(bool larger); 00405 00406 //@} 00407 00408 /// \name Regrid function 00409 //@{ 00410 /** 00411 \brief Regrid the data 00412 00413 The uses cubic spline interpolation to refine the data set, 00414 ideally used before attempting to calculate the contour 00415 levels. If the original number of data points is 00416 \f$(\mathrm{nx},\mathrm{ny})\f$, then the new number of data 00417 points is 00418 \f[ 00419 (\mathrm{xfact}~(\mathrm{nx}-1)+1, 00420 \mathrm{yfact}~(\mathrm{ny}-1)+1) 00421 \f] 00422 */ 00423 int regrid_data(size_t xfact, size_t yfact); 00424 //@} 00425 00426 /// \name Obtain internal data 00427 //@{ 00428 /** 00429 \brief Get the data 00430 00431 This is useful to see how the data has changed after 00432 a call to regrid_data(). 00433 */ 00434 int get_data(size_t &sizex, size_t &sizey, const ovector *&x_fun, 00435 const ovector *&y_fun, const ovector **&udata) { 00436 if (nx==0) { 00437 set_err_ret("Data not set in calc().",gsl_einval); 00438 } 00439 sizex=nx; 00440 sizey=ny; 00441 x_fun=xfun; 00442 y_fun=yfun; 00443 udata=(const ovector **)data; 00444 return 0; 00445 } 00446 00447 /** \brief Return the edges used to compute the contours 00448 */ 00449 int get_edges(const std::vector<omatrix_int *> *rints, 00450 const std::vector<omatrix_int *> *bints, 00451 const std::vector<omatrix *> *rpoints, 00452 const std::vector<omatrix *> *bpoints); 00453 00454 /** 00455 \brief Return the edges used to compute the contours 00456 00457 This function allocates memory for \c rints, \c bints, 00458 \c rpoints, and \c bpoints and fills them with a copy of 00459 the data that the class is using. As such, there is no 00460 need for them to be const. 00461 */ 00462 int get_edges_for_level(size_t nl, omatrix_int &rints, omatrix_int &bints, 00463 omatrix &rpoints, omatrix &bpoints); 00464 00465 //@} 00466 00467 /// Verbosity parameter 00468 int verbose; 00469 00470 /// (default \f$ 10^{-8} \f$) 00471 double lev_adjust; 00472 00473 #ifndef DOXYGEN_INTERNAL 00474 00475 protected: 00476 00477 /// \name Edge direction 00478 //@{ 00479 static const int dright=0; 00480 static const int dbottom=1; 00481 //@} 00482 00483 /// \name Edge status 00484 //@{ 00485 static const int empty=0; 00486 static const int edge=1; 00487 static const int contourp=2; 00488 static const int endpoint=3; 00489 //@} 00490 00491 /// \name Edge found or not found 00492 //@{ 00493 static const int efound=1; 00494 static const int enot_found=0; 00495 //@} 00496 00497 /// \name User-specified data 00498 //@{ 00499 int nx, ny; 00500 ovector *xfun, *yfun, **data; 00501 //@} 00502 00503 int new_debug; 00504 00505 /// \name User-specified contour levels 00506 //@{ 00507 int nlev; 00508 ovector levels; 00509 bool levels_set; 00510 //@} 00511 00512 /// \name Generated curves 00513 //@{ 00514 int ncurves; 00515 std::vector<int> csizes; 00516 std::vector<double> vals; 00517 std::vector< std::vector<double> > xc, yc; 00518 //@} 00519 00520 /// \name Storage for edges 00521 //@{ 00522 std::vector<omatrix *> redges, bedges; 00523 std::vector<omatrix_int *> re, be; 00524 omatrix_int *rep; 00525 omatrix_int *bep; 00526 omatrix *redgesp; 00527 omatrix *bedgesp; 00528 //@} 00529 00530 /// Object to find if a point is inside a polygon 00531 pinside pi; 00532 00533 /// Find next point starting from a point on a right edge 00534 int find_next_point_right(int j, int k, int &jnext, int &knext, 00535 int &dir_next, int nsw=1); 00536 00537 /// Find next point starting from a point on a bottom edge 00538 int find_next_point_bottom(int j, int k, int &jnext, int &knext, 00539 int &dir_next, int nsw=1); 00540 00541 /// Find all of the intersections of the edges with the contour level 00542 int find_intersections(size_t ilev, double &level); 00543 00544 /// Interpolate all right edge crossings 00545 int right_edges(double level, o2scl::sm_interp *si); 00546 00547 /// Interpolate all bottom edge crossings 00548 int bottom_edges(double level, o2scl::sm_interp *si); 00549 00550 /// Create a contour line from a starting edge 00551 int process_line(int j, int k, int dir, std::vector<double> &x, 00552 std::vector<double> &y, bool first=true); 00553 00554 /** 00555 \brief Test if a point is inside a closed contour (unfinished) 00556 00557 This returns true if the point (x1,y1) is "inside" the contour 00558 (i.e. a collection of line segments) given in \c x and \c 00559 y. The arrays \c x and \c y must be "ordered" so that adjacent 00560 points are placed at adjacent entries. The result is 00561 undefined if this is not the case or if the contours are not 00562 properly closed. The first and last points in \c x and \c y 00563 should be the same to indicate a closed contour. The values 00564 xscale and yscale should be an approximate scale for the 00565 contours \c x and \c y. 00566 00567 \note This function is deprecated and has been replaced by 00568 \ref pinside 00569 */ 00570 bool is_point_inside_old(double x1, double y1, const ovector_view &x, 00571 const ovector_view &y, double xscale=0.01, 00572 double yscale=0.01); 00573 00574 /// Free memory 00575 int free_memory(); 00576 00577 /** 00578 \brief Check if lines cross 00579 00580 Returns true if the line segment defined by (x1,y1) and (x2,y2) 00581 and the line segment defined by (x3,y3) and (x4,y4) have 00582 an intersection point that is between the endpoints of both 00583 segments. The function handles vertical and horizontal lines 00584 appropriately. This function will fail if x1==y1 and x2==y2 00585 or if x3==y3 and x4==y4, i.e. if the points do not really 00586 define a line. If the function fails, it returns false. 00587 00588 \note This function is deprecated and has been replaced by 00589 \ref pinside 00590 */ 00591 bool lines_cross_old(double x1, double y1, double x2, double y2, 00592 double x3, double y3, double x4, double y4); 00593 00594 /// Check to ensure the x- and y-arrays are monotonic 00595 int check_data(); 00596 00597 /** \brief Smooth the contours by adding internal points using cubic 00598 interpolation (this doesn't work) 00599 00600 This makes the contours smoother by adding internal points 00601 between each contour line segment determined by cubic spline 00602 interpolation. 00603 00604 For more accurate contours, it is better to provide the 00605 original data on a finer grid, or use regrid_data(). 00606 */ 00607 int smooth_contours(size_t nfact); 00608 00609 #endif 00610 00611 }; 00612 00613 #ifndef DOXYGENP 00614 } 00615 #endif 00616 00617 #endif 00618 00619
Documentation generated with Doxygen and provided under the GNU Free Documentation License. See License Information for details.
Project hosting provided by
,
O2scl Sourceforge Project Page