# include # include # include # include # include # include using namespace std; int main ( int argc, char *argv[] ); char ch_cap ( char c ); bool ch_eqi ( char c1, char c2 ); int ch_to_digit ( char c ); int diaedg ( double x0, double y0, double x1, double y1, double x2, double y2, double x3, double y3 ); int dtris2 ( int point_num, int base, double point_xy[], int *tri_num, int tri_vert[], int tri_nabe[] ); int file_column_count ( string input_file_name ); int file_row_count ( string input_file_name ); int i4_max ( int i1, int i2 ); int i4_min ( int i1, int i2 ); int i4_modp ( int i, int j ); int i4_sign ( int i ); int i4_wrap ( int ival, int ilo, int ihi ); void i4mat_transpose_print_some ( int m, int n, int a[], int ilo, int jlo, int ihi, int jhi, string title ); void i4mat_write ( string output_filename, int m, int n, int table[] ); int *i4vec_indicator_new ( int n ); int i4vec_min ( int n, int a[] ); int lrline ( double xu, double yu, double xv1, double yv1, double xv2, double yv2, double dv ); bool perm_check ( int n, int p[], int base ); void perm_inverse ( int n, int p[] ); double r8_abs ( double x ); double r8_epsilon ( ); double r8_huge ( ); double r8_max ( double x, double y ); double r8_min ( double x, double y ); void r82vec_permute ( int n, int p[], int base, double a[] ); int *r82vec_sort_heap_index_a ( int n, int base, double a[] ); double *r8mat_data_read ( string input_filename, int m, int n ); void r8mat_header_read ( string input_filename, int *m, int *n ); void r8mat_transpose_print_some ( int m, int n, double a[], int ilo, int jlo, int ihi, int jhi, string title ); int s_len_trim ( string s ); double s_to_r8 ( string s, int *lchar, bool *error ); bool s_to_r8vec ( string s, int n, double rvec[] ); int s_word_count ( string s ); int swapec ( int i, int *top, int *btri, int *bedg, int point_num, double point_xy[], int tri_num, int tri_vert[], int tri_nabe[], int stack[] ); void timestamp ( ); void vbedg ( double x, double y, int point_num, double point_xy[], int tri_num, int tri_vert[], int tri_nabe[], int *ltri, int *ledg, int *rtri, int *redg ); //****************************************************************************80 int main ( int argc, char *argv[] ) //****************************************************************************80 // // Purpose: // // MAIN is the main program for TABLE_DELAUNAY. // // Discussion: // // TABLE_DELAUNAY computes the Delaunay triangulation of a TABLE dataset. // // The dataset is simply a set of points in the plane. // // Thus, given a set of points V1, V2, ..., VN, we apply a standard // Delaunay triangulation. The Delaunay triangulation is an organization // of the data into triples, forming a triangulation of the data, with // the property that the circumcircle of each triangle never contains // another data point. // // Usage: // // table_delaunay prefix // // where: // // 'prefix' is the common prefix for the node and triangle files: // // * prefix_nodes.txt, the node coordinates (input). // * prefix_elements.txt, the nodes that make up each triangle (output). // // Licensing: // // This code is distributed under the GNU LGPL license. // // Modified: // // 22 June 2006 // // Author: // // John Burkardt // { int base = 1; int node_dim; string node_filename; int node_num; double *node_xy; string prefix; string triangle_filename; int *triangle_neighbor; int *triangle_node; int triangle_num; int triangle_order; cout << "\n"; timestamp ( ); cout << "\n"; cout << "TABLE_DELAUNAY\n"; cout << " C++ version:\n"; cout << "\n"; cout << " Compiled on " << __DATE__ << " at " << __TIME__ << ".\n"; cout << "\n"; cout << " Read a TABLE dataset of N points in 2 dimensions,\n"; cout << " Compute the Delaunay triangulation.\n"; cout << " Write an integer TABLE dataset of the triangulation.\n"; // // First argument is the file prefix. // if ( argc <= 1 ) { cout << "\n"; cout << "TABLE_DELAUNAY:\n"; cout << " Please enter the file prefix.\n"; cin >> prefix; } else { prefix = argv[1]; } // // Create the filenames. // node_filename = prefix + "_nodes.txt"; triangle_filename = prefix + "_elements.txt"; // // Read the point coordinates. // r8mat_header_read ( node_filename, &node_dim, &node_num ); cout << "\n"; cout << " Read the header of \"" << node_filename << "\"\n"; cout << "\n"; cout << " Node dimension NODE_DIM = " << node_dim << "\n"; cout << " Node number NODE_NUM = " << node_num << "\n"; if ( node_dim != 2 ) { cout << "\n"; cout << "TABLE_DELAUNAY - Fatal error!\n"; cout << " The node dimension is not 2.\n"; exit ( 1 ); } node_xy = r8mat_data_read ( node_filename, node_dim, node_num ); cout << "\n"; cout << " Read the data of \"" << node_filename << "\"\n"; r8mat_transpose_print_some ( node_dim, node_num, node_xy, 1, 1, node_dim, 5, " Initial portion of node data:" ); // // Determine the Delaunay triangulation. // triangle_order = 3; triangle_node = new int[triangle_order*3*node_num]; triangle_neighbor = new int[triangle_order*3*node_num]; dtris2 ( node_num, base, node_xy, &triangle_num, triangle_node, triangle_neighbor ); // // Print a portion of the triangulation. // cout << "\n"; cout << " Computed the triangulation.\n"; cout << " Number of triangles is " << triangle_num << "\n"; i4mat_transpose_print_some ( triangle_order, triangle_num, triangle_node, 1, 1, 3, 5, " Initial portion of triangulation data:" ); // // Write the triangulation to a file. // i4mat_write ( triangle_filename, triangle_order, triangle_num, triangle_node ); cout << "\n"; cout << " Wrote the triangulation data to \"" << triangle_filename << "\".\n"; // // Terminate execution. // cout << "\n"; cout << "TABLE_DELAUNAY:\n"; cout << " Normal end of execution.\n"; cout << "\n"; timestamp ( ); delete [] node_xy; delete [] triangle_node; delete [] triangle_neighbor; return 0; } //****************************************************************************80 char ch_cap ( char c ) //****************************************************************************80 // // Purpose: // // CH_CAP capitalizes a single character. // // Discussion: // // This routine should be equivalent to the library "toupper" function. // // Licensing: // // This code is distributed under the GNU LGPL license. // // Modified: // // 19 July 1998 // // Author: // // John Burkardt // // Parameters: // // Input, char C, the character to capitalize. // // Output, char CH_CAP, the capitalized character. // { if ( 97 <= c && c <= 122 ) { c = c - 32; } return c; } //****************************************************************************80 bool ch_eqi ( char c1, char c2 ) //****************************************************************************80 // // Purpose: // // CH_EQI is true if two characters are equal, disregarding case. // // Licensing: // // This code is distributed under the GNU LGPL license. // // Modified: // // 13 June 2003 // // Author: // // John Burkardt // // Parameters: // // Input, char C1, char C2, the characters to compare. // // Output, bool CH_EQI, is true if the two characters are equal, // disregarding case. // { if ( 97 <= c1 && c1 <= 122 ) { c1 = c1 - 32; } if ( 97 <= c2 && c2 <= 122 ) { c2 = c2 - 32; } return ( c1 == c2 ); } //****************************************************************************80 int ch_to_digit ( char c ) //****************************************************************************80 // // Purpose: // // CH_TO_DIGIT returns the integer value of a base 10 digit. // // Example: // // C DIGIT // --- ----- // '0' 0 // '1' 1 // ... ... // '9' 9 // ' ' 0 // 'X' -1 // // Licensing: // // This code is distributed under the GNU LGPL license. // // Modified: // // 13 June 2003 // // Author: // // John Burkardt // // Parameters: // // Input, char C, the decimal digit, '0' through '9' or blank are legal. // // Output, int CH_TO_DIGIT, the corresponding integer value. If C was // 'illegal', then DIGIT is -1. // { int digit; if ( '0' <= c && c <= '9' ) { digit = c - '0'; } else if ( c == ' ' ) { digit = 0; } else { digit = -1; } return digit; } //****************************************************************************80 int diaedg ( double x0, double y0, double x1, double y1, double x2, double y2, double x3, double y3 ) //****************************************************************************80 // // Purpose: // // DIAEDG chooses a diagonal edge. // // Discussion: // // The routine determines whether 0--2 or 1--3 is the diagonal edge // that should be chosen, based on the circumcircle criterion, where // (X0,Y0), (X1,Y1), (X2,Y2), (X3,Y3) are the vertices of a simple // quadrilateral in counterclockwise order. // // Licensing: // // This code is distributed under the GNU LGPL license. // // Modified: // // 28 August 2003 // // Author: // // Original FORTRAN77 version by Barry Joe. // C++ version by John Burkardt. // // Reference: // // Barry Joe, // GEOMPACK - a software package for the generation of meshes // using geometric algorithms, // Advances in Engineering Software, // Volume 13, pages 325-331, 1991. // // Parameters: // // Input, double X0, Y0, X1, Y1, X2, Y2, X3, Y3, the coordinates of the // vertices of a quadrilateral, given in counter clockwise order. // // Output, int DIAEDG, chooses a diagonal: // +1, if diagonal edge 02 is chosen; // -1, if diagonal edge 13 is chosen; // 0, if the four vertices are cocircular. // { double ca; double cb; double dx10; double dx12; double dx30; double dx32; double dy10; double dy12; double dy30; double dy32; double s; double tol; double tola; double tolb; int value; tol = 100.0 * r8_epsilon ( ); dx10 = x1 - x0; dy10 = y1 - y0; dx12 = x1 - x2; dy12 = y1 - y2; dx30 = x3 - x0; dy30 = y3 - y0; dx32 = x3 - x2; dy32 = y3 - y2; tola = tol * r8_max ( fabs ( dx10 ), r8_max ( fabs ( dy10 ), r8_max ( fabs ( dx30 ), fabs ( dy30 ) ) ) ); tolb = tol * r8_max ( fabs ( dx12 ), r8_max ( fabs ( dy12 ), r8_max ( fabs ( dx32 ), fabs ( dy32 ) ) ) ); ca = dx10 * dx30 + dy10 * dy30; cb = dx12 * dx32 + dy12 * dy32; if ( tola < ca && tolb < cb ) { value = -1; } else if ( ca < -tola && cb < -tolb ) { value = 1; } else { tola = r8_max ( tola, tolb ); s = ( dx10 * dy30 - dx30 * dy10 ) * cb + ( dx32 * dy12 - dx12 * dy32 ) * ca; if ( tola < s ) { value = -1; } else if ( s < -tola ) { value = 1; } else { value = 0; } } return value; } //****************************************************************************80 int dtris2 ( int point_num, int base, double point_xy[], int *tri_num, int tri_vert[], int tri_nabe[] ) //****************************************************************************80 // // Purpose: // // DTRIS2 constructs a Delaunay triangulation of 2D vertices. // // Discussion: // // The routine constructs the Delaunay triangulation of a set of 2D vertices // using an incremental approach and diagonal edge swaps. Vertices are // first sorted in lexicographically increasing (X,Y) order, and // then are inserted one at a time from outside the convex hull. // // Licensing: // // This code is distributed under the GNU LGPL license. // // Modified: // // 03 June 2009 // // Author: // // Original FORTRAN77 version by Barry Joe. // C++ version by John Burkardt. // // Reference: // // Barry Joe, // GEOMPACK - a software package for the generation of meshes // using geometric algorithms, // Advances in Engineering Software, // Volume 13, pages 325-331, 1991. // // Parameters: // // Input, int POINT_NUM, the number of vertices. // // Input, int BASE, the base for the indexing of TRI_VERT. // 0, use 0-based indexing. // 1, use 1-based indexing. // // Input/output, double POINT_XY[POINT_NUM*2], the coordinates of the vertices. // On output, the vertices have been sorted into dictionary order. // // Output, int *TRI_NUM, the number of triangles in the triangulation; // TRI_NUM is equal to 2*POINT_NUM - NB - 2, where NB is the number // of boundary vertices. // // Output, int TRI_VERT[TRI_NUM*3], the nodes that make up each triangle. // The elements are indices of POINT_XY. The vertices of the triangles are // in counter clockwise order. // // Output, int TRI_NABE[TRI_NUM*3], the triangle neighbor list. // Positive elements are indices of TIL; negative elements are used for links // of a counter clockwise linked list of boundary edges; LINK = -(3*I + J-1) // where I, J = triangle, edge index; TRI_NABE[I,J] refers to // the neighbor along edge from vertex J to J+1 (mod 3). // // Output, int DTRIS2, is 0 for no error. { double cmax; int e; int error; int i; int *indx; int j; int k; int l; int ledg; int lr; int ltri; int m; int m1; int m2; int n; int redg; int rtri; int *stack; int t; double tol; int top; stack = new int[point_num]; tol = 100.0 * r8_epsilon ( ); // // Sort the vertices by increasing (x,y). // indx = r82vec_sort_heap_index_a ( point_num, base, point_xy ); r82vec_permute ( point_num, indx, base, point_xy ); // // Make sure that the data points are "reasonably" distinct. // m1 = 1; for ( i = 2; i <= point_num; i++ ) { m = m1; m1 = i; k = -1; for ( j = 0; j <= 1; j++ ) { cmax = r8_max ( fabs ( point_xy[2*(m-1)+j] ), fabs ( point_xy[2*(m1-1)+j] ) ); if ( tol * ( cmax + 1.0 ) < fabs ( point_xy[2*(m-1)+j] - point_xy[2*(m1-1)+j] ) ) { k = j; break; } } if ( k == -1 ) { cout << "\n"; cout << "DTRIS2 - Fatal error!\n"; cout << " Fails for point number I = " << i << "\n"; cout << " M = " << m << "\n"; cout << " M1 = " << m1 << "\n"; cout << " X,Y(M) = " << point_xy[2*(m-1)+0] << " " << point_xy[2*(m-1)+1] << "\n"; cout << " X,Y(M1) = " << point_xy[2*(m1-1)+0] << " " << point_xy[2*(m1-1)+1] << "\n"; delete [] stack; return 224; } } // // Starting from points M1 and M2, search for a third point M that // makes a "healthy" triangle (M1,M2,M) // m1 = 1; m2 = 2; j = 3; for ( ; ; ) { if ( point_num < j ) { cout << "\n"; cout << "DTRIS2 - Fatal error!\n"; delete [] stack; return 225; } m = j; lr = lrline ( point_xy[2*(m-1)+0], point_xy[2*(m-1)+1], point_xy[2*(m1-1)+0], point_xy[2*(m1-1)+1], point_xy[2*(m2-1)+0], point_xy[2*(m2-1)+1], 0.0 ); if ( lr != 0 ) { break; } j = j + 1; } // // Set up the triangle information for (M1,M2,M), and for any other // triangles you created because points were collinear with M1, M2. // *tri_num = j - 2; if ( lr == -1 ) { tri_vert[3*0+0] = m1; tri_vert[3*0+1] = m2; tri_vert[3*0+2] = m; tri_nabe[3*0+2] = -3; for ( i = 2; i <= *tri_num; i++ ) { m1 = m2; m2 = i+1; tri_vert[3*(i-1)+0] = m1; tri_vert[3*(i-1)+1] = m2; tri_vert[3*(i-1)+2] = m; tri_nabe[3*(i-1)+0] = -3 * i; tri_nabe[3*(i-1)+1] = i; tri_nabe[3*(i-1)+2] = i - 1; } tri_nabe[3*(*tri_num-1)+0] = -3 * (*tri_num) - 1; tri_nabe[3*(*tri_num-1)+1] = -5; ledg = 2; ltri = *tri_num; } else { tri_vert[3*0+0] = m2; tri_vert[3*0+1] = m1; tri_vert[3*0+2] = m; tri_nabe[3*0+0] = -4; for ( i = 2; i <= *tri_num; i++ ) { m1 = m2; m2 = i+1; tri_vert[3*(i-1)+0] = m2; tri_vert[3*(i-1)+1] = m1; tri_vert[3*(i-1)+2] = m; tri_nabe[3*(i-2)+2] = i; tri_nabe[3*(i-1)+0] = -3 * i - 3; tri_nabe[3*(i-1)+1] = i - 1; } tri_nabe[3*(*tri_num-1)+2] = -3 * (*tri_num); tri_nabe[3*0+1] = -3 * (*tri_num) - 2; ledg = 2; ltri = 1; } // // Insert the vertices one at a time from outside the convex hull, // determine visible boundary edges, and apply diagonal edge swaps until // Delaunay triangulation of vertices (so far) is obtained. // top = 0; for ( i = j+1; i <= point_num; i++ ) { m = i; m1 = tri_vert[3*(ltri-1)+ledg-1]; if ( ledg <= 2 ) { m2 = tri_vert[3*(ltri-1)+ledg]; } else { m2 = tri_vert[3*(ltri-1)+0]; } lr = lrline ( point_xy[2*(m-1)+0], point_xy[2*(m-1)+1], point_xy[2*(m1-1)+0], point_xy[2*(m1-1)+1], point_xy[2*(m2-1)+0], point_xy[2*(m2-1)+1], 0.0 ); if ( 0 < lr ) { rtri = ltri; redg = ledg; ltri = 0; } else { l = -tri_nabe[3*(ltri-1)+ledg-1]; rtri = l / 3; redg = (l % 3) + 1; } vbedg ( point_xy[2*(m-1)+0], point_xy[2*(m-1)+1], point_num, point_xy, *tri_num, tri_vert, tri_nabe, <ri, &ledg, &rtri, &redg ); n = *tri_num + 1; l = -tri_nabe[3*(ltri-1)+ledg-1]; for ( ; ; ) { t = l / 3; e = ( l % 3 ) + 1; l = -tri_nabe[3*(t-1)+e-1]; m2 = tri_vert[3*(t-1)+e-1]; if ( e <= 2 ) { m1 = tri_vert[3*(t-1)+e]; } else { m1 = tri_vert[3*(t-1)+0]; } *tri_num = *tri_num + 1; tri_nabe[3*(t-1)+e-1] = *tri_num; tri_vert[3*(*tri_num-1)+0] = m1; tri_vert[3*(*tri_num-1)+1] = m2; tri_vert[3*(*tri_num-1)+2] = m; tri_nabe[3*(*tri_num-1)+0] = t; tri_nabe[3*(*tri_num-1)+1] = *tri_num - 1; tri_nabe[3*(*tri_num-1)+2] = *tri_num + 1; top = top + 1; if ( point_num < top ) { cout << "\n"; cout << "DTRIS2 - Fatal error!\n"; cout << " Stack overflow.\n"; delete [] stack; return 8; } stack[top-1] = *tri_num; if ( t == rtri && e == redg ) { break; } } tri_nabe[3*(ltri-1)+ledg-1] = -3 * n - 1; tri_nabe[3*(n-1)+1] = -3 * (*tri_num) - 2; tri_nabe[3*(*tri_num-1)+2] = -l; ltri = n; ledg = 2; error = swapec ( m, &top, <ri, &ledg, point_num, point_xy, *tri_num, tri_vert, tri_nabe, stack ); if ( error != 0 ) { cout << "\n"; cout << "DTRIS2 - Fatal error!\n"; cout << " Error return from SWAPEC.\n"; delete [] stack; return error; } } // // Now account for the sorting that we did. // for ( i = 0; i < 3; i++ ) { for ( j = 0; j < *tri_num; j++ ) { tri_vert[i+j*3] = indx [ tri_vert[i+j*3] - 1 ]; } } perm_inverse ( point_num, indx ); r82vec_permute ( point_num, indx, base, point_xy ); delete [] indx; delete [] stack; return 0; } //****************************************************************************80 int file_column_count ( string input_filename ) //****************************************************************************80 // // Purpose: // // FILE_COLUMN_COUNT counts the number of columns in the first line of a file. // // Discussion: // // The file is assumed to be a simple text file. // // Most lines of the file is presumed to consist of COLUMN_NUM words, separated // by spaces. There may also be some blank lines, and some comment lines, // which have a "#" in column 1. // // The routine tries to find the first non-comment non-blank line and // counts the number of words in that line. // // If all lines are blanks or comments, it goes back and tries to analyze // a comment line. // // Licensing: // // This code is distributed under the GNU LGPL license. // // Modified: // // 23 February 2009 // // Author: // // John Burkardt // // Parameters: // // Input, string INPUT_FILENAME, the name of the file. // // Output, int FILE_COLUMN_COUNT, the number of columns assumed // to be in the file. // { int column_num; ifstream input; bool got_one; char line[255]; // // Open the file. // input.open ( input_filename.c_str ( ) ); if ( !input ) { column_num = -1; cerr << "\n"; cerr << "FILE_COLUMN_COUNT - Fatal error!\n"; cerr << " Could not open the file:\n"; cerr << " \"" << input_filename << "\"\n"; return column_num; } // // Read one line, but skip blank lines and comment lines. // got_one = false; for ( ; ; ) { input.getline ( line, sizeof ( line ) ); if ( input.eof ( ) ) { break; } if ( s_len_trim ( line ) == 0 ) { continue; } if ( line[0] == '#' ) { continue; } got_one = true; break; } if ( !got_one ) { input.close ( ); input.open ( input_filename.c_str ( ) ); for ( ; ; ) { input.getline ( line, sizeof ( line ) ); if ( input.eof ( ) ) { break; } if ( s_len_trim ( line ) == 0 ) { continue; } got_one = true; break; } } input.close ( ); if ( !got_one ) { cerr << "\n"; cerr << "FILE_COLUMN_COUNT - Warning!\n"; cerr << " The file does not seem to contain any data.\n"; return -1; } column_num = s_word_count ( line ); return column_num; } //****************************************************************************80 int file_row_count ( string input_filename ) //****************************************************************************80 // // Purpose: // // FILE_ROW_COUNT counts the number of row records in a file. // // Discussion: // // It does not count lines that are blank, or that begin with a // comment symbol '#'. // // Licensing: // // This code is distributed under the GNU LGPL license. // // Modified: // // 23 February 2009 // // Author: // // John Burkardt // // Parameters: // // Input, string INPUT_FILENAME, the name of the input file. // // Output, int FILE_ROW_COUNT, the number of rows found. // { int comment_num; ifstream input; char line[255]; int record_num; int row_num; row_num = 0; comment_num = 0; record_num = 0; input.open ( input_filename.c_str ( ) ); if ( !input ) { cerr << "\n"; cerr << "FILE_ROW_COUNT - Fatal error!\n"; cerr << " Could not open the input file: \"" << input_filename << "\"\n"; return (-1); } for ( ; ; ) { input.getline ( line, sizeof ( line ) ); if ( input.eof ( ) ) { break; } record_num = record_num + 1; if ( line[0] == '#' ) { comment_num = comment_num + 1; continue; } if ( s_len_trim ( line ) == 0 ) { comment_num = comment_num + 1; continue; } row_num = row_num + 1; } input.close ( ); return row_num; } //****************************************************************************80 int i4_max ( int i1, int i2 ) //****************************************************************************80 // // Purpose: // // I4_MAX returns the maximum of two I4's. // // Licensing: // // This code is distributed under the GNU LGPL license. // // Modified: // // 13 October 1998 // // Author: // // John Burkardt // // Parameters: // // Input, int I1, I2, are two integers to be compared. // // Output, int I4_MAX, the larger of I1 and I2. // { int value; if ( i2 < i1 ) { value = i1; } else { value = i2; } return value; } //****************************************************************************80 int i4_min ( int i1, int i2 ) //****************************************************************************80 // // Purpose: // // I4_MIN returns the minimum of two I4's. // // Licensing: // // This code is distributed under the GNU LGPL license. // // Modified: // // 13 October 1998 // // Author: // // John Burkardt // // Parameters: // // Input, int I1, I2, two integers to be compared. // // Output, int I4_MIN, the smaller of I1 and I2. // { int value; if ( i1 < i2 ) { value = i1; } else { value = i2; } return value; } //****************************************************************************80 int i4_modp ( int i, int j ) //****************************************************************************80 // // Purpose: // // I4_MODP returns the nonnegative remainder of I4 division. // // Discussion: // // If // NREM = I4_MODP ( I, J ) // NMULT = ( I - NREM ) / J // then // I = J * NMULT + NREM // where NREM is always nonnegative. // // The MOD function computes a result with the same sign as the // quantity being divided. Thus, suppose you had an angle A, // and you wanted to ensure that it was between 0 and 360. // Then mod(A,360) would do, if A was positive, but if A // was negative, your result would be between -360 and 0. // // On the other hand, I4_MODP(A,360) is between 0 and 360, always. // // I J MOD I4_MODP I4_MODP Factorization // // 107 50 7 7 107 = 2 * 50 + 7 // 107 -50 7 7 107 = -2 * -50 + 7 // -107 50 -7 43 -107 = -3 * 50 + 43 // -107 -50 -7 43 -107 = 3 * -50 + 43 // // Licensing: // // This code is distributed under the GNU LGPL license. // // Modified: // // 26 May 1999 // // Author: // // John Burkardt // // Parameters: // // Input, int I, the number to be divided. // // Input, int J, the number that divides I. // // Output, int I4_MODP, the nonnegative remainder when I is // divided by J. // { int value; if ( j == 0 ) { cout << "\n"; cout << "I4_MODP - Fatal error!\n"; cout << " I4_MODP ( I, J ) called with J = " << j << "\n"; exit ( 1 ); } value = i % j; if ( value < 0 ) { value = value + abs ( j ); } return value; } //****************************************************************************80 int i4_sign ( int i ) //****************************************************************************80 // // Purpose: // // I4_SIGN returns the sign of an I4. // // Licensing: // // This code is distributed under the GNU LGPL license. // // Modified: // // 27 March 2004 // // Author: // // John Burkardt // // Parameters: // // Input, int I, the integer whose sign is desired. // // Output, int I4_SIGN, the sign of I. { int value; if ( i < 0 ) { value = -1; } else { value = 1; } return value; } //****************************************************************************80 int i4_wrap ( int ival, int ilo, int ihi ) //****************************************************************************80 // // Purpose: // // I4_WRAP forces an I4 to lie between given limits by wrapping. // // Example: // // ILO = 4, IHI = 8 // // I Value // // -2 8 // -1 4 // 0 5 // 1 6 // 2 7 // 3 8 // 4 4 // 5 5 // 6 6 // 7 7 // 8 8 // 9 4 // 10 5 // 11 6 // 12 7 // 13 8 // 14 4 // // Licensing: // // This code is distributed under the GNU LGPL license. // // Modified: // // 19 August 2003 // // Author: // // John Burkardt // // Parameters: // // Input, int IVAL, an integer value. // // Input, int ILO, IHI, the desired bounds for the integer value. // // Output, int I4_WRAP, a "wrapped" version of IVAL. // { int jhi; int jlo; int value; int wide; jlo = i4_min ( ilo, ihi ); jhi = i4_max ( ilo, ihi ); wide = jhi + 1 - jlo; if ( wide == 1 ) { value = jlo; } else { value = jlo + i4_modp ( ival - jlo, wide ); } return value; } //****************************************************************************80 void i4mat_transpose_print_some ( int m, int n, int a[], int ilo, int jlo, int ihi, int jhi, string title ) //****************************************************************************80 // // Purpose: // // I4MAT_TRANSPOSE_PRINT_SOME prints some of an I4MAT, transposed. // // Discussion: // // An I4MAT is an MxN array of I4's, stored by (I,J) -> [I+J*M]. // // Licensing: // // This code is distributed under the GNU LGPL license. // // Modified: // // 14 June 2005 // // Author: // // John Burkardt // // Parameters: // // Input, int M, the number of rows of the matrix. // M must be positive. // // Input, int N, the number of columns of the matrix. // N must be positive. // // Input, int A[M*N], the matrix. // // Input, int ILO, JLO, IHI, JHI, designate the first row and // column, and the last row and column to be printed. // // Input, string TITLE, a title. // { # define INCX 10 int i; int i2hi; int i2lo; int j; int j2hi; int j2lo; cout << "\n"; cout << title << "\n"; // // Print the columns of the matrix, in strips of INCX. // for ( i2lo = ilo; i2lo <= ihi; i2lo = i2lo + INCX ) { i2hi = i2lo + INCX - 1; i2hi = i4_min ( i2hi, m ); i2hi = i4_min ( i2hi, ihi ); cout << "\n"; // // For each row I in the current range... // // Write the header. // cout << " Row: "; for ( i = i2lo; i <= i2hi; i++ ) { cout << setw(6) << i << " "; } cout << "\n"; cout << " Col\n"; cout << "\n"; // // Determine the range of the rows in this strip. // j2lo = i4_max ( jlo, 1 ); j2hi = i4_min ( jhi, n ); for ( j = j2lo; j <= j2hi; j++ ) { // // Print out (up to INCX) entries in column J, that lie in the current strip. // cout << setw(5) << j << " "; for ( i = i2lo; i <= i2hi; i++ ) { cout << setw(6) << a[i-1+(j-1)*m] << " "; } cout << "\n"; } } return; # undef INCX } //****************************************************************************80 void i4mat_write ( string output_filename, int m, int n, int table[] ) //****************************************************************************80 // // Purpose: // // I4MAT_WRITE writes an I4MAT file with no header. // // Licensing: // // This code is distributed under the GNU LGPL license. // // Modified: // // 01 June 2009 // // Author: // // John Burkardt // // Parameters: // // Input, string OUTPUT_FILENAME, the output filename. // // Input, int M, the spatial dimension. // // Input, int N, the number of points. // // Input, int TABLE[M*N], the table data. // { int i; int j; ofstream output; // // Open the file. // output.open ( output_filename.c_str ( ) ); if ( !output ) { cerr << "\n"; cerr << "I4MAT_WRITE - Fatal error!\n"; cerr << " Could not open the output file.\n"; return; } // // Write the data. // for ( j = 0; j < n; j++ ) { for ( i = 0; i < m; i++ ) { output << " " << setw(10) << table[i+j*m]; } output << "\n"; } // // Close the file. // output.close ( ); return; } //****************************************************************************80 int *i4vec_indicator_new ( int n ) //****************************************************************************80 // // Purpose: // // I4VEC_INDICATOR_NEW sets an I4VEC to the indicator vector. // // Discussion: // // An I4VEC is a vector of I4's. // // Licensing: // // This code is distributed under the GNU LGPL license. // // Modified: // // 03 June 2009 // // Author: // // John Burkardt // // Parameters: // // Input, int N, the number of elements of A. // // Output, int I4VEC_INDICATOR_NEW[N], the array. // { int *a; int i; a = new int[n]; for ( i = 0; i < n; i++ ) { a[i] = i + 1; } return a; } //****************************************************************************80 int i4vec_min ( int n, int a[] ) //****************************************************************************80 // // Purpose: // // I4VEC_MIN returns the value of the minimum element in an I4VEC. // // Discussion: // // An I4VEC is a vector of I4's. // // Licensing: // // This code is distributed under the GNU LGPL license. // // Modified: // // 17 May 2003 // // Author: // // John Burkardt // // Parameters: // // Input, int N, the number of entries in the array. // // Input, int A[N], the array to be checked. // // Output, int I4VEC_MIN, the value of the minimum element. This // is set to 0 if N <= 0. // { int i; int value; if ( n <= 0 ) { return 0; } value = a[0]; for ( i = 1; i < n; i++ ) { if ( a[i] < value ) { value = a[i]; } } return value; } //****************************************************************************80 int lrline ( double xu, double yu, double xv1, double yv1, double xv2, double yv2, double dv ) //****************************************************************************80 // // Purpose: // // LRLINE determines where a point lies in relation to a directed line. // // Discussion: // // LRLINE determines whether a point is to the left of, right of, // or on a directed line parallel to a line through given points. // // Licensing: // // This code is distributed under the GNU LGPL license. // // Modified: // // 22 June 2009 // // Author: // // Original FORTRAN77 version by Barry Joe. // C++ version by John Burkardt. // // Reference: // // Barry Joe, // GEOMPACK - a software package for the generation of meshes // using geometric algorithms, // Advances in Engineering Software, // Volume 13, pages 325-331, 1991. // // Parameters: // // Input, double XU, YU, XV1, YV1, XV2, YV2, are vertex coordinates; the // directed line is parallel to and at signed distance DV to the left of // the directed line from (XV1,YV1) to (XV2,YV2); (XU,YU) is the vertex for // which the position relative to the directed line is to be determined. // // Input, double DV, the signed distance, positive for left. // // Output, int LRLINE, is +1, 0, or -1 depending on whether (XU,YU) is // to the right of, on, or left of the directed line. LRLINE is 0 if // the line degenerates to a point. // { double dx; double dxu; double dy; double dyu; double t; double tol; double tolabs; int value; tol = 100.0 * r8_epsilon ( ); dx = xv2 - xv1; dy = yv2 - yv1; dxu = xu - xv1; dyu = yu - yv1; tolabs = tol * r8_max ( fabs ( dx ), r8_max ( fabs ( dy ), r8_max ( fabs ( dxu ), r8_max ( fabs ( dyu ), fabs ( dv ) ) ) ) ); t = dy * dxu - dx * dyu + dv * sqrt ( dx * dx + dy * dy ); if ( tolabs < t ) { value = 1; } else if ( -tolabs <= t ) { value = 0; } else if ( t < -tolabs ) { value = -1; } return value; } //****************************************************************************80 bool perm_check ( int n, int p[], int base ) //****************************************************************************80 // // Purpose: // // PERM_CHECK checks that a vector represents a permutation. // // Discussion: // // The routine verifies that each of the integers from BASE to // to BASE+N-1 occurs among the N entries of the permutation. // // Set the input quantity BASE to 0, if P is a 0-based permutation, // or to 1 if P is a 1-based permutation. // // Licensing: // // This code is distributed under the GNU LGPL license. // // Modified: // // 30 October 2008 // // Author: // // John Burkardt // // Parameters: // // Input, int N, the number of entries. // // Input, int P[N], the array to check. // // Input, int BASE, the index base. // // Output, bool PERM_CHECK, is TRUE if the permutation is OK. // { bool found; int i; int seek; for ( seek = base; seek < base + n; seek++ ) { found = false; for ( i = 0; i < n; i++ ) { if ( p[i] == seek ) { found = true; break; } } if ( !found ) { return false; } } return true; } //****************************************************************************80 void perm_inverse ( int n, int p[] ) //****************************************************************************80 // // Purpose: // // PERM_INVERSE inverts a permutation "in place". // // Discussion: // // This algorithm assumes that the entries in the permutation vector are // strictly positive. In particular, the value 0 must no occur. // // When necessary, this function shifts the data temporarily so that // this requirement is satisfied. // // Licensing: // // This code is distributed under the GNU LGPL license. // // Modified: // // 13 January 2004 // // Author: // // John Burkardt // // Parameters: // // Input, int N, the number of objects being permuted. // // Input/output, int P[N], the permutation, in standard index form. // On output, P describes the inverse permutation // { int base; int i; int i0; int i1; int i2; int is; int p_min; if ( n <= 0 ) { cout << "\n"; cout << "PERM_INVERSE - Fatal error!\n"; cout << " Input value of N = " << n << "\n"; exit ( 1 ); } // // Find the least value, and shift data so it begins at 1. // p_min = i4vec_min ( n, p ); base = 1; for ( i = 0; i < n; i++ ) { p[i] = p[i] - p_min + base; } // // Now we can safely check the permutation. // if ( !perm_check ( n, p, base ) ) { cerr << "\n"; cerr << "PERM_INVERSE - Fatal error!\n"; cerr << " PERM_CHECK rejects this permutation.\n"; exit ( 1 ); } // // Now we can invert the permutation. // is = 1; for ( i = 1; i <= n; i++ ) { i1 = p[i-1]; while ( i < i1 ) { i2 = p[i1-1]; p[i1-1] = -i2; i1 = i2; } is = - i4_sign ( p[i-1] ); p[i-1] = i4_sign ( is ) * abs ( p[i-1] ); } for ( i = 1; i <= n; i++ ) { i1 = - p[i-1]; if ( 0 <= i1 ) { i0 = i; for ( ; ; ) { i2 = p[i1-1]; p[i1-1] = i0; if ( i2 < 0 ) { break; } i0 = i1; i1 = i2; } } } // // Now we can restore the permutation. // for ( i = 0; i < n; i++ ) { p[i] = p[i] + p_min - base; } return; } //****************************************************************************80 double r8_abs ( double x ) //****************************************************************************80 // // Purpose: // // R8_ABS returns the absolute value of an R8. // // Licensing: // // This code is distributed under the GNU LGPL license. // // Modified: // // 14 November 2006 // // Author: // // John Burkardt // // Parameters: // // Input, double X, the quantity whose absolute value is desired. // // Output, double R8_ABS, the absolute value of X. // { double value; if ( 0.0 <= x ) { value = x; } else { value = - x; } return value; } //****************************************************************************80 double r8_epsilon ( ) //****************************************************************************80 // // Purpose: // // R8_EPSILON returns the R8 roundoff unit. // // Discussion: // // The roundoff unit is a number R which is a power of 2 with the // property that, to the precision of the computer's arithmetic, // 1 < 1 + R // but // 1 = ( 1 + R / 2 ) // // Licensing: // // This code is distributed under the GNU LGPL license. // // Modified: // // 01 September 2012 // // Author: // // John Burkardt // // Parameters: // // Output, double R8_EPSILON, the R8 round-off unit. // { const double value = 2.220446049250313E-016; return value; } //****************************************************************************80 double r8_huge ( ) //****************************************************************************80 // // Purpose: // // R8_HUGE returns a "huge" R8. // // Discussion: // // The value returned by this function is NOT required to be the // maximum representable R8. This value varies from machine to machine, // from compiler to compiler, and may cause problems when being printed. // We simply want a "very large" but non-infinite number. // // Licensing: // // This code is distributed under the GNU LGPL license. // // Modified: // // 06 October 2007 // // Author: // // John Burkardt // // Parameters: // // Output, double R8_HUGE, a "huge" R8 value. // { double value; value = 1.0E+30; return value; } //****************************************************************************80 double r8_max ( double x, double y ) //****************************************************************************80 // // Purpose: // // R8_MAX returns the maximum of two R8's. // // Licensing: // // This code is distributed under the GNU LGPL license. // // Modified: // // 18 August 2004 // // Author: // // John Burkardt // // Parameters: // // Input, double X, Y, the quantities to compare. // // Output, double R8_MAX, the maximum of X and Y. // { double value; if ( y < x ) { value = x; } else { value = y; } return value; } //****************************************************************************80 double r8_min ( double x, double y ) //****************************************************************************80 // // Purpose: // // R8_MIN returns the minimum of two R8's. // // Licensing: // // This code is distributed under the GNU LGPL license. // // Modified: // // 31 August 2004 // // Author: // // John Burkardt // // Parameters: // // Input, double X, Y, the quantities to compare. // // Output, double R8_MIN, the minimum of X and Y. // { double value; if ( y < x ) { value = y; } else { value = x; } return value; } //****************************************************************************80 int r8_nint ( double x ) //****************************************************************************80 // // Purpose: // // R8_NINT returns the nearest integer to an R8. // // Example: // // X Value // // 1.3 1 // 1.4 1 // 1.5 1 or 2 // 1.6 2 // 0.0 0 // -0.7 -1 // -1.1 -1 // -1.6 -2 // // Licensing: // // This code is distributed under the GNU LGPL license. // // Modified: // // 26 August 2004 // // Author: // // John Burkardt // // Parameters: // // Input, double X, the value. // // Output, int R8_NINT, the nearest integer to X. // { int value; if ( x < 0.0 ) { value = - ( int ) ( r8_abs ( x ) + 0.5 ); } else { value = ( int ) ( r8_abs ( x ) + 0.5 ); } return value; } //****************************************************************************80 void r82vec_permute ( int n, int p[], int base, double a[] ) //****************************************************************************80 // // Purpose: // // R82VEC_PERMUTE permutes an R82VEC in place. // // Discussion: // // An R82VEC is a vector whose entries are R82's. // An R82 is a vector of type double precision with two entries. // An R82VEC may be stored as a 2 by N array. // // This routine permutes an array of real "objects", but the same // logic can be used to permute an array of objects of any arithmetic // type, or an array of objects of any complexity. The only temporary // storage required is enough to store a single object. The number // of data movements made is N + the number of cycles of order 2 or more, // which is never more than N + N/2. // // Example: // // Input: // // N = 5 // P = ( 2, 4, 5, 1, 3 ) // A = ( 1.0, 2.0, 3.0, 4.0, 5.0 ) // (11.0, 22.0, 33.0, 44.0, 55.0 ) // // Output: // // A = ( 2.0, 4.0, 5.0, 1.0, 3.0 ) // ( 22.0, 44.0, 55.0, 11.0, 33.0 ). // // Licensing: // // This code is distributed under the GNU LGPL license. // // Modified: // // 30 October 2008 // // Author: // // John Burkardt // // Parameters: // // Input, int N, the number of objects. // // Input, int P[N], the permutation. P(I) = J means // that the I-th element of the output array should be the J-th // element of the input array. // // Input, int BASE, is 0 for a 0-based permutation and 1 for a 1-based permutation. // // Input/output, double A[2*N], the array to be permuted. // { double a_temp[2]; int i; int iget; int iput; int istart; if ( !perm_check ( n, p, base ) ) { cerr << "\n"; cerr << "R82VEC_PERMUTE - Fatal error!\n"; cerr << " PERM_CHECK rejects this permutation.\n"; exit ( 1 ); } // // In order for the sign negation trick to work, we need to assume that the // entries of P are strictly positive. Presumably, the lowest number is BASE. // So temporarily add 1-BASE to each entry to force positivity. // for ( i = 0; i < n; i++ ) { p[i] = p[i] + 1 - base; } // // Search for the next element of the permutation that has not been used. // for ( istart = 1; istart <= n; istart++ ) { if ( p[istart-1] < 0 ) { continue; } else if ( p[istart-1] == istart ) { p[istart-1] = - p[istart-1]; continue; } else { a_temp[0] = a[0+(istart-1)*2]; a_temp[1] = a[1+(istart-1)*2]; iget = istart; // // Copy the new value into the vacated entry. // for ( ; ; ) { iput = iget; iget = p[iget-1]; p[iput-1] = - p[iput-1]; if ( iget < 1 || n < iget ) { cout << "\n"; cout << "R82VEC_PERMUTE - Fatal error!\n"; cout << " Entry IPUT = " << iput << " of the permutation has\n"; cout << " an illegal value IGET = " << iget << ".\n"; exit ( 1 ); } if ( iget == istart ) { a[0+(iput-1)*2] = a_temp[0]; a[1+(iput-1)*2] = a_temp[1]; break; } a[0+(iput-1)*2] = a[0+(iget-1)*2]; a[1+(iput-1)*2] = a[1+(iget-1)*2]; } } } // // Restore the signs of the entries. // for ( i = 0; i < n; i++ ) { p[i] = - p[i]; } // // Restore the base of the entries. // for ( i = 0; i < n; i++ ) { p[i] = p[i] - 1 + base; } return; } //****************************************************************************80 int *r82vec_sort_heap_index_a ( int n, int base, double a[] ) //****************************************************************************80 // // Purpose: // // R82VEC_SORT_HEAP_INDEX_A does an indexed heap ascending sort of an R82VEC. // // Discussion: // // An R82VEC is a vector whose entries are R82's. // An R82 is a vector of type double precision with two entries. // An R82VEC may be stored as a 2 by N array. // // The sorting is not actually carried out. Rather an index array is // created which defines the sorting. This array may be used to sort // or index the array, or to sort or index related arrays keyed on the // original array. // // Once the index array is computed, the sorting can be carried out // "implicitly: // // a(*,indx(*)) // // or explicitly, by the call // // r82vec_permute ( n, indx, base, a ) // // after which a(*,*) is sorted. // // Licensing: // // This code is distributed under the GNU LGPL license. // // Modified: // // 03 June 2009 // // Author: // // John Burkardt // // Parameters: // // Input, int N, the number of entries in the array. // // Input, int BASE, the desired indexing for the sort index: // 0 for 0-based indexing, // 1 for 1-based indexing. // // Input, double A[2*N], an array to be index-sorted. // // Output, int R82VEC_SORT_HEAP_INDEX_A[N], the sort index. The // I-th element of the sorted array is A(0:1,R8VEC_SORT_HEAP_INDEX_A(I)). // { double aval[2]; int i; int *indx; int indxt; int ir; int j; int l; if ( n < 1 ) { return NULL; } indx = new int[n]; for ( i = 0; i < n; i++ ) { indx[i] = i; } if ( n == 1 ) { indx[0] = indx[0] + base; return indx; } l = n / 2 + 1; ir = n; for ( ; ; ) { if ( 1 < l ) { l = l - 1; indxt = indx[l-1]; aval[0] = a[0+indxt*2]; aval[1] = a[1+indxt*2]; } else { indxt = indx[ir-1]; aval[0] = a[0+indxt*2]; aval[1] = a[1+indxt*2]; indx[ir-1] = indx[0]; ir = ir - 1; if ( ir == 1 ) { indx[0] = indxt; break; } } i = l; j = l + l; while ( j <= ir ) { if ( j < ir ) { if ( a[0+indx[j-1]*2] < a[0+indx[j]*2] || ( a[0+indx[j-1]*2] == a[0+indx[j]*2] && a[1+indx[j-1]*2] < a[1+indx[j]*2] ) ) { j = j + 1; } } if ( aval[0] < a[0+indx[j-1]*2] || ( aval[0] == a[0+indx[j-1]*2] && aval[1] < a[1+indx[j-1]*2] ) ) { indx[i-1] = indx[j-1]; i = j; j = j + j; } else { j = ir + 1; } } indx[i-1] = indxt; } // // Take care of the base. // for ( i = 0; i < n; i++ ) { indx[i] = indx[i] + base; } return indx; } //****************************************************************************80 double *r8mat_data_read ( string input_filename, int m, int n ) //****************************************************************************80 // // Purpose: // // R8MAT_DATA_READ reads the data from an R8MAT file. // // Discussion: // // The file is assumed to contain one record per line. // // Records beginning with '#' are comments, and are ignored. // Blank lines are also ignored. // // Each line that is not ignored is assumed to contain exactly (or at least) // M real numbers, representing the coordinates of a point. // // There are assumed to be exactly (or at least) N such records. // // Licensing: // // This code is distributed under the GNU LGPL license. // // Modified: // // 23 February 2009 // // Author: // // John Burkardt // // Parameters: // // Input, string INPUT_FILENAME, the name of the input file. // // Input, int M, the number of spatial dimensions. // // Input, int N, the number of points. The program // will stop reading data once N values have been read. // // Output, double R8MAT_DATA_READ[M*N], the table data. // { bool error; ifstream input; int i; int j; string line; double *table; double *x; input.open ( input_filename.c_str ( ) ); if ( !input ) { cerr << "\n"; cerr << "R8MAT_DATA_READ - Fatal error!\n"; cerr << " Could not open the input file: \"" << input_filename << "\"\n"; return NULL; } table = new double[m*n]; x = new double[m]; j = 0; while ( j < n ) { getline ( input, line ); if ( input.eof ( ) ) { break; } if ( line[0] == '#' || s_len_trim ( line ) == 0 ) { continue; } error = s_to_r8vec ( line, m, x ); if ( error ) { continue; } for ( i = 0; i < m; i++ ) { table[i+j*m] = x[i]; } j = j + 1; } input.close ( ); delete [] x; return table; } //****************************************************************************80 void r8mat_header_read ( string input_filename, int *m, int *n ) //****************************************************************************80 // // Purpose: // // R8MAT_HEADER_READ reads the header from an R8MAT file. // // Licensing: // // This code is distributed under the GNU LGPL license. // // Modified: // // 23 February 2009 // // Author: // // John Burkardt // // Parameters: // // Input, string INPUT_FILENAME, the name of the input file. // // Output, int *M, the number of spatial dimensions. // // Output, int *N, the number of points. // { *m = file_column_count ( input_filename ); if ( *m <= 0 ) { cerr << "\n"; cerr << "R8MAT_HEADER_READ - Fatal error!\n"; cerr << " FILE_COLUMN_COUNT failed.\n"; *n = -1; return; } *n = file_row_count ( input_filename ); if ( *n <= 0 ) { cerr << "\n"; cerr << "R8MAT_HEADER_READ - Fatal error!\n"; cerr << " FILE_ROW_COUNT failed.\n"; return; } return; } //****************************************************************************80 void r8mat_transpose_print_some ( int m, int n, double a[], int ilo, int jlo, int ihi, int jhi, string title ) //****************************************************************************80 // // Purpose: // // R8MAT_TRANSPOSE_PRINT_SOME prints some of an R8MAT, transposed. // // Discussion: // // An R8MAT is a doubly dimensioned array of R8 values, stored as a vector // in column-major order. // // Licensing: // // This code is distributed under the GNU LGPL license. // // Modified: // // 10 September 2009 // // Author: // // John Burkardt // // Parameters: // // Input, int M, N, the number of rows and columns. // // Input, double A[M*N], an M by N matrix to be printed. // // Input, int ILO, JLO, the first row and column to print. // // Input, int IHI, JHI, the last row and column to print. // // Input, string TITLE, a title. // { # define INCX 5 int i; int i2; int i2hi; int i2lo; int inc; int j; int j2hi; int j2lo; cout << "\n"; cout << title << "\n"; for ( i2lo = i4_max ( ilo, 1 ); i2lo <= i4_min ( ihi, m ); i2lo = i2lo + INCX ) { i2hi = i2lo + INCX - 1; i2hi = i4_min ( i2hi, m ); i2hi = i4_min ( i2hi, ihi ); inc = i2hi + 1 - i2lo; cout << "\n"; cout << " Row: "; for ( i = i2lo; i <= i2hi; i++ ) { cout << setw(7) << i << " "; } cout << "\n"; cout << " Col\n"; cout << "\n"; j2lo = i4_max ( jlo, 1 ); j2hi = i4_min ( jhi, n ); for ( j = j2lo; j <= j2hi; j++ ) { cout << setw(5) << j << " "; for ( i2 = 1; i2 <= inc; i2++ ) { i = i2lo - 1 + i2; cout << setw(14) << a[(i-1)+(j-1)*m]; } cout << "\n"; } } return; # undef INCX } //****************************************************************************80 int s_len_trim ( string s ) //****************************************************************************80 // // Purpose: // // S_LEN_TRIM returns the length of a string to the last nonblank. // // Licensing: // // This code is distributed under the GNU LGPL license. // // Modified: // // 05 July 2009 // // Author: // // John Burkardt // // Parameters: // // Input, string S, a string. // // Output, int S_LEN_TRIM, the length of the string to the last nonblank. // If S_LEN_TRIM is 0, then the string is entirely blank. // { int n; n = s.length ( ); while ( 0 < n ) { if ( s[n-1] != ' ' ) { return n; } n = n - 1; } return n; } //****************************************************************************80 double s_to_r8 ( string s, int *lchar, bool *error ) //****************************************************************************80 // // Purpose: // // S_TO_R8 reads an R8 from a string. // // Discussion: // // This routine will read as many characters as possible until it reaches // the end of the string, or encounters a character which cannot be // part of the real number. // // Legal input is: // // 1 blanks, // 2 '+' or '-' sign, // 2.5 spaces // 3 integer part, // 4 decimal point, // 5 fraction part, // 6 'E' or 'e' or 'D' or 'd', exponent marker, // 7 exponent sign, // 8 exponent integer part, // 9 exponent decimal point, // 10 exponent fraction part, // 11 blanks, // 12 final comma or semicolon. // // with most quantities optional. // // Example: // // S R // // '1' 1.0 // ' 1 ' 1.0 // '1A' 1.0 // '12,34,56' 12.0 // ' 34 7' 34.0 // '-1E2ABCD' -100.0 // '-1X2ABCD' -1.0 // ' 2E-1' 0.2 // '23.45' 23.45 // '-4.2E+2' -420.0 // '17d2' 1700.0 // '-14e-2' -0.14 // 'e2' 100.0 // '-12.73e-9.23' -12.73 * 10.0**(-9.23) // // Licensing: // // This code is distributed under the GNU LGPL license. // // Modified: // // 05 July 2009 // // Author: // // John Burkardt // // Parameters: // // Input, string S, the string containing the // data to be read. Reading will begin at position 1 and // terminate at the end of the string, or when no more // characters can be read to form a legal real. Blanks, // commas, or other nonnumeric data will, in particular, // cause the conversion to halt. // // Output, int *LCHAR, the number of characters read from // the string to form the number, including any terminating // characters such as a trailing comma or blanks. // // Output, bool *ERROR, is true if an error occurred. // // Output, double S_TO_R8, the real value that was read from the string. // { char c; int ihave; int isgn; int iterm; int jbot; int jsgn; int jtop; int nchar; int ndig; double r; double rbot; double rexp; double rtop; char TAB = 9; nchar = s_len_trim ( s ); *error = false; r = 0.0; *lchar = -1; isgn = 1; rtop = 0.0; rbot = 1.0; jsgn = 1; jtop = 0; jbot = 1; ihave = 1; iterm = 0; for ( ; ; ) { c = s[*lchar+1]; *lchar = *lchar + 1; // // Blank or TAB character. // if ( c == ' ' || c == TAB ) { if ( ihave == 2 ) { } else if ( ihave == 6 || ihave == 7 ) { iterm = 1; } else if ( 1 < ihave ) { ihave = 11; } } // // Comma. // else if ( c == ',' || c == ';' ) { if ( ihave != 1 ) { iterm = 1; ihave = 12; *lchar = *lchar + 1; } } // // Minus sign. // else if ( c == '-' ) { if ( ihave == 1 ) { ihave = 2; isgn = -1; } else if ( ihave == 6 ) { ihave = 7; jsgn = -1; } else { iterm = 1; } } // // Plus sign. // else if ( c == '+' ) { if ( ihave == 1 ) { ihave = 2; } else if ( ihave == 6 ) { ihave = 7; } else { iterm = 1; } } // // Decimal point. // else if ( c == '.' ) { if ( ihave < 4 ) { ihave = 4; } else if ( 6 <= ihave && ihave <= 8 ) { ihave = 9; } else { iterm = 1; } } // // Exponent marker. // else if ( ch_eqi ( c, 'E' ) || ch_eqi ( c, 'D' ) ) { if ( ihave < 6 ) { ihave = 6; } else { iterm = 1; } } // // Digit. // else if ( ihave < 11 && '0' <= c && c <= '9' ) { if ( ihave <= 2 ) { ihave = 3; } else if ( ihave == 4 ) { ihave = 5; } else if ( ihave == 6 || ihave == 7 ) { ihave = 8; } else if ( ihave == 9 ) { ihave = 10; } ndig = ch_to_digit ( c ); if ( ihave == 3 ) { rtop = 10.0 * rtop + ( double ) ndig; } else if ( ihave == 5 ) { rtop = 10.0 * rtop + ( double ) ndig; rbot = 10.0 * rbot; } else if ( ihave == 8 ) { jtop = 10 * jtop + ndig; } else if ( ihave == 10 ) { jtop = 10 * jtop + ndig; jbot = 10 * jbot; } } // // Anything else is regarded as a terminator. // else { iterm = 1; } // // If we haven't seen a terminator, and we haven't examined the // entire string, go get the next character. // if ( iterm == 1 || nchar <= *lchar + 1 ) { break; } } // // If we haven't seen a terminator, and we have examined the // entire string, then we're done, and LCHAR is equal to NCHAR. // if ( iterm != 1 && (*lchar) + 1 == nchar ) { *lchar = nchar; } // // Number seems to have terminated. Have we got a legal number? // Not if we terminated in states 1, 2, 6 or 7! // if ( ihave == 1 || ihave == 2 || ihave == 6 || ihave == 7 ) { *error = true; return r; } // // Number seems OK. Form it. // if ( jtop == 0 ) { rexp = 1.0; } else { if ( jbot == 1 ) { rexp = pow ( 10.0, jsgn * jtop ); } else { rexp = jsgn * jtop; rexp = rexp / jbot; rexp = pow ( 10.0, rexp ); } } r = isgn * rexp * rtop / rbot; return r; } //****************************************************************************80 bool s_to_r8vec ( string s, int n, double rvec[] ) //****************************************************************************80 // // Purpose: // // S_TO_R8VEC reads an R8VEC from a string. // // Licensing: // // This code is distributed under the GNU LGPL license. // // Modified: // // 05 July 2009 // // Author: // // John Burkardt // // Parameters: // // Input, string S, the string to be read. // // Input, int N, the number of values expected. // // Output, double RVEC[N], the values read from the string. // // Output, bool S_TO_R8VEC, is true if an error occurred. // { int begin; bool error; int i; int lchar; int length; begin = 0; length = s.length ( ); error = 0; for ( i = 0; i < n; i++ ) { rvec[i] = s_to_r8 ( s.substr(begin,length), &lchar, &error ); if ( error ) { return error; } begin = begin + lchar; length = length - lchar; } return error; } //****************************************************************************80 int s_word_count ( string s ) //****************************************************************************80 // // Purpose: // // S_WORD_COUNT counts the number of "words" in a string. // // Licensing: // // This code is distributed under the GNU LGPL license. // // Modified: // // 05 July 2009 // // Author: // // John Burkardt // // Parameters: // // Input, string S, the string to be examined. // // Output, int S_WORD_COUNT, the number of "words" in the string. // Words are presumed to be separated by one or more blanks. // { bool blank; int char_count; int i; int word_count; word_count = 0; blank = true; char_count = s.length ( ); for ( i = 0; i < char_count; i++ ) { if ( isspace ( s[i] ) ) { blank = true; } else if ( blank ) { word_count = word_count + 1; blank = false; } } return word_count; } //****************************************************************************80 int swapec ( int i, int *top, int *btri, int *bedg, int point_num, double point_xy[], int tri_num, int tri_vert[], int tri_nabe[], int stack[] ) //****************************************************************************80 // // Purpose: // // SWAPEC swaps diagonal edges until all triangles are Delaunay. // // Discussion: // // The routine swaps diagonal edges in a 2D triangulation, based on // the empty circumcircle criterion, until all triangles are Delaunay, // given that I is the index of the new vertex added to the triangulation. // // Licensing: // // This code is distributed under the GNU LGPL license. // // Modified: // // 03 September 2003 // // Author: // // Original FORTRAN77 version by Barry Joe. // C++ version by John Burkardt. // // Reference: // // Barry Joe, // GEOMPACK - a software package for the generation of meshes // using geometric algorithms, // Advances in Engineering Software, // Volume 13, pages 325-331, 1991. // // Parameters: // // Input, int I, the index of the new vertex. // // Input/output, int *TOP, the index of the top of the stack. // On output, TOP is zero. // // Input/output, int *BTRI, *BEDG; on input, if positive, are the // triangle and edge indices of a boundary edge whose updated indices // must be recorded. On output, these may be updated because of swaps. // // Input, int POINT_NUM, the number of points. // // Input, double POINT_XY[POINT_NUM*2], the coordinates of the points. // // Input, int TRI_NUM, the number of triangles. // // Input/output, int TRI_VERT[TRI_NUM*3], the triangle incidence list. // May be updated on output because of swaps. // // Input/output, int TRI_NABE[TRI_NUM*3], the triangle neighbor list; // negative values are used for links of the counter-clockwise linked // list of boundary edges; May be updated on output because of swaps. // // LINK = -(3*I + J-1) where I, J = triangle, edge index. // // Workspace, int STACK[MAXST]; on input, entries 1 through TOP // contain the indices of initial triangles (involving vertex I) // put in stack; the edges opposite I should be in interior; entries // TOP+1 through MAXST are used as a stack. // // Output, int SWAPEC, is set to 8 for abnormal return. // { int a; int b; int c; int e; int ee; int em1; int ep1; int f; int fm1; int fp1; int l; int r; int s; int swap; int t; int tt; int u; double x; double y; // // Determine whether triangles in stack are Delaunay, and swap // diagonal edge of convex quadrilateral if not. // x = point_xy[2*(i-1)+0]; y = point_xy[2*(i-1)+1]; for ( ; ; ) { if ( *top <= 0 ) { break; } t = stack[(*top)-1]; *top = *top - 1; if ( tri_vert[3*(t-1)+0] == i ) { e = 2; b = tri_vert[3*(t-1)+2]; } else if ( tri_vert[3*(t-1)+1] == i ) { e = 3; b = tri_vert[3*(t-1)+0]; } else { e = 1; b = tri_vert[3*(t-1)+1]; } a = tri_vert[3*(t-1)+e-1]; u = tri_nabe[3*(t-1)+e-1]; if ( tri_nabe[3*(u-1)+0] == t ) { f = 1; c = tri_vert[3*(u-1)+2]; } else if ( tri_nabe[3*(u-1)+1] == t ) { f = 2; c = tri_vert[3*(u-1)+0]; } else { f = 3; c = tri_vert[3*(u-1)+1]; } swap = diaedg ( x, y, point_xy[2*(a-1)+0], point_xy[2*(a-1)+1], point_xy[2*(c-1)+0], point_xy[2*(c-1)+1], point_xy[2*(b-1)+0], point_xy[2*(b-1)+1] ); if ( swap == 1 ) { em1 = i4_wrap ( e - 1, 1, 3 ); ep1 = i4_wrap ( e + 1, 1, 3 ); fm1 = i4_wrap ( f - 1, 1, 3 ); fp1 = i4_wrap ( f + 1, 1, 3 ); tri_vert[3*(t-1)+ep1-1] = c; tri_vert[3*(u-1)+fp1-1] = i; r = tri_nabe[3*(t-1)+ep1-1]; s = tri_nabe[3*(u-1)+fp1-1]; tri_nabe[3*(t-1)+ep1-1] = u; tri_nabe[3*(u-1)+fp1-1] = t; tri_nabe[3*(t-1)+e-1] = s; tri_nabe[3*(u-1)+f-1] = r; if ( 0 < tri_nabe[3*(u-1)+fm1-1] ) { *top = *top + 1; stack[(*top)-1] = u; } if ( 0 < s ) { if ( tri_nabe[3*(s-1)+0] == u ) { tri_nabe[3*(s-1)+0] = t; } else if ( tri_nabe[3*(s-1)+1] == u ) { tri_nabe[3*(s-1)+1] = t; } else { tri_nabe[3*(s-1)+2] = t; } *top = *top + 1; if ( point_num < *top ) { return 8; } stack[(*top)-1] = t; } else { if ( u == *btri && fp1 == *bedg ) { *btri = t; *bedg = e; } l = - ( 3 * t + e - 1 ); tt = t; ee = em1; while ( 0 < tri_nabe[3*(tt-1)+ee-1] ) { tt = tri_nabe[3*(tt-1)+ee-1]; if ( tri_vert[3*(tt-1)+0] == a ) { ee = 3; } else if ( tri_vert[3*(tt-1)+1] == a ) { ee = 1; } else { ee = 2; } } tri_nabe[3*(tt-1)+ee-1] = l; } if ( 0 < r ) { if ( tri_nabe[3*(r-1)+0] == t ) { tri_nabe[3*(r-1)+0] = u; } else if ( tri_nabe[3*(r-1)+1] == t ) { tri_nabe[3*(r-1)+1] = u; } else { tri_nabe[3*(r-1)+2] = u; } } else { if ( t == *btri && ep1 == *bedg ) { *btri = u; *bedg = f; } l = - ( 3 * u + f - 1 ); tt = u; ee = fm1; while ( 0 < tri_nabe[3*(tt-1)+ee-1] ) { tt = tri_nabe[3*(tt-1)+ee-1]; if ( tri_vert[3*(tt-1)+0] == b ) { ee = 3; } else if ( tri_vert[3*(tt-1)+1] == b ) { ee = 1; } else { ee = 2; } } tri_nabe[3*(tt-1)+ee-1] = l; } } } return 0; } //****************************************************************************80 void timestamp ( ) //****************************************************************************80 // // Purpose: // // TIMESTAMP prints the current YMDHMS date as a time stamp. // // Example: // // May 31 2001 09:45:54 AM // // Licensing: // // This code is distributed under the GNU LGPL license. // // Modified: // // 21 August 2002 // // Author: // // John Burkardt // // Parameters: // // None // { # define TIME_SIZE 29 static char time_buffer[TIME_SIZE]; const struct tm *tm; time_t now; now = time ( NULL ); tm = localtime ( &now ); strftime ( time_buffer, TIME_SIZE, "%d %B %Y %I:%M:%S %p", tm ); cout << time_buffer << "\n"; return; # undef TIME_SIZE } //****************************************************************************80 void vbedg ( double x, double y, int point_num, double point_xy[], int tri_num, int tri_vert[], int tri_nabe[], int *ltri, int *ledg, int *rtri, int *redg ) //****************************************************************************80 // // Purpose: // // VBEDG determines which boundary edges are visible to a point. // // Discussion: // // The point (X,Y) is assumed to be outside the convex hull of the // region covered by the 2D triangulation. // // Licensing: // // This code is distributed under the GNU LGPL license. // // Modified: // // 30 December 2008 // // Author: // // Original FORTRAN77 version by Barry Joe. // C++ version by John Burkardt. // // Reference: // // Barry Joe, // GEOMPACK - a software package for the generation of meshes // using geometric algorithms, // Advances in Engineering Software, // Volume 13, pages 325-331, 1991. // // Modified: // // 02 September 2003 // // Parameters: // // Input, double X, Y, the coordinates of a point outside the convex hull // of the current triangulation. // // Input, int POINT_NUM, the number of points. // // Input, double POINT_XY[POINT_NUM*2], the coordinates of the vertices. // // Input, int TRI_NUM, the number of triangles. // // Input, int TRI_VERT[TRI_NUM*3], the triangle incidence list. // // Input, int TRI_NABE[TRI_NUM*3], the triangle neighbor list; negative // values are used for links of a counter clockwise linked list of boundary // edges; // LINK = -(3*I + J-1) where I, J = triangle, edge index. // // Input/output, int *LTRI, *LEDG. If LTRI != 0 then these values are // assumed to be already computed and are not changed, else they are updated. // On output, LTRI is the index of boundary triangle to the left of the // leftmost boundary triangle visible from (X,Y), and LEDG is the boundary // edge of triangle LTRI to the left of the leftmost boundary edge visible // from (X,Y). 1 <= LEDG <= 3. // // Input/output, int *RTRI. On input, the index of the boundary triangle // to begin the search at. On output, the index of the rightmost boundary // triangle visible from (X,Y). // // Input/output, int *REDG, the edge of triangle RTRI that is visible // from (X,Y). 1 <= REDG <= 3. // { int a; double ax; double ay; int b; double bx; double by; bool done; int e; int l; int lr; int t; // // Find the rightmost visible boundary edge using links, then possibly // leftmost visible boundary edge using triangle neighbor information. // if ( *ltri == 0 ) { done = false; *ltri = *rtri; *ledg = *redg; } else { done = true; } for ( ; ; ) { l = -tri_nabe[3*((*rtri)-1)+(*redg)-1]; t = l / 3; e = 1 + l % 3; a = tri_vert[3*(t-1)+e-1]; if ( e <= 2 ) { b = tri_vert[3*(t-1)+e]; } else { b = tri_vert[3*(t-1)+0]; } ax = point_xy[2*(a-1)+0]; ay = point_xy[2*(a-1)+1]; bx = point_xy[2*(b-1)+0]; by = point_xy[2*(b-1)+1]; lr = lrline ( x, y, ax, ay, bx, by, 0.0 ); if ( lr <= 0 ) { break; } *rtri = t; *redg = e; } if ( done ) { return; } t = *ltri; e = *ledg; for ( ; ; ) { b = tri_vert[3*(t-1)+e-1]; e = i4_wrap ( e-1, 1, 3 ); while ( 0 < tri_nabe[3*(t-1)+e-1] ) { t = tri_nabe[3*(t-1)+e-1]; if ( tri_vert[3*(t-1)+0] == b ) { e = 3; } else if ( tri_vert[3*(t-1)+1] == b ) { e = 1; } else { e = 2; } } a = tri_vert[3*(t-1)+e-1]; ax = point_xy[2*(a-1)+0]; ay = point_xy[2*(a-1)+1]; bx = point_xy[2*(b-1)+0]; by = point_xy[2*(b-1)+1]; lr = lrline ( x, y, ax, ay, bx, by, 0.0 ); if ( lr <= 0 ) { break; } } *ltri = t; *ledg = e; return; }