# include # include # include # include # include # include using namespace std; # include "fem_sample.H" //****************************************************************************80 void basis_mn_t3 ( double t[2*3], int n, double p[], double phi[], double dphidx[], double dphidy[] ) //****************************************************************************80 // // Purpose: // // BASIS_MN_T3: all bases at N points for a T3 element. // // Discussion: // // The routine is given the coordinates of the vertices of a triangle. // It works directly with these coordinates, and does not refer to a // reference element. // // The sides of the triangle DO NOT have to lie along a coordinate // axis. // // The routine evaluates the basis functions associated with each vertex, // and their derivatives with respect to X and Y. // // Physical Element T3: // // 3 // / \ // / \ // / \ // / \ // 1---------2 // // Modified: // // 09 February 2006 // // Author: // // John Burkardt // // Parameters: // // Input, double T[2*3], the coordinates of the vertices // of the triangle. It is common to list these points in counter clockwise // order. // // Input, int N, the number of evaluation points. // // Input, double P[2*N], the points where the basis functions // are to be evaluated. // // Output, double PHI[3*N], the value of the basis functions // at the evaluation points. // // Output, double DPHIDX[3*N], DPHIDY[3*N], the value of the // derivatives at the evaluation points. // // Local parameters: // // Local, double AREA, is (twice) the area of the triangle. // { double area; int i; int j; double temp; area = t[0+0*2] * ( t[1+1*2] - t[1+2*2] ) + t[0+1*2] * ( t[1+2*2] - t[1+0*2] ) + t[0+2*2] * ( t[1+0*2] - t[1+1*2] ); if ( area == 0.0 ) { cout << "\n"; cout << "BASIS_MN_T3 - Fatal error!\n"; cout << " Element has zero area.\n"; exit ( 1 ); } for ( j = 0; j < n; j++ ) { phi[0+j*3] = ( ( t[0+2*2] - t[0+1*2] ) * ( p[1+j*2] - t[1+1*2] ) - ( t[1+2*2] - t[1+1*2] ) * ( p[0+j*2] - t[0+1*2] ) ); dphidx[0+j*3] = - ( t[1+2*2] - t[1+1*2] ); dphidy[0+j*3] = ( t[0+2*2] - t[0+1*2] ); phi[1+j*3] = ( ( t[0+0*2] - t[0+2*2] ) * ( p[1+j*2] - t[1+2*2] ) - ( t[1+0*2] - t[1+2*2] ) * ( p[0+j*2] - t[0+2*2] ) ); dphidx[1+j*3] = - ( t[1+0*2] - t[1+2*2] ); dphidy[1+j*3] = ( t[0+0*2] - t[0+2*2] ); phi[2+j*3] = ( ( t[0+1*2] - t[0+0*2] ) * ( p[1+j*2] - t[1+0*2] ) - ( t[1+1*2] - t[1+0*2] ) * ( p[0+j*2] - t[0+0*2] ) ); dphidx[2+j*3] = - ( t[1+1*2] - t[1+0*2] ); dphidy[2+j*3] = ( t[0+1*2] - t[0+0*2] ); } // // Normalize. // for ( j = 0; j < n; j++ ) { for ( i = 0; i < 3; i++ ) { phi[i+j*3] = phi[i+j*3] / area; dphidx[i+j*3] = dphidx[i+j*3] / area; dphidy[i+j*3] = dphidy[i+j*3] / area; } } return; } //****************************************************************************80 void basis_mn_t6 ( double t[2*6], int n, double p[], double phi[], double dphidx[], double dphidy[] ) //****************************************************************************80 // // Purpose: // // BASIS_MN_T6: all bases at N points for a T6 element. // // Discussion: // // The routine is given the coordinates of the vertices and midside // nodes of a triangle. It works directly with these coordinates, and does // not refer to a reference element. // // This routine requires that the midside nodes be "in line" // with the vertices, that is, that the sides of the triangle be // straight. However, the midside nodes do not actually have to // be halfway along the side of the triangle. // // The physical element T6: // // This picture indicates the assumed ordering of the six nodes // of the triangle. // // | // | // | 3 // | / \ // | / \ // Y 6 5 // | / \ // | / \ // | 1-----4-----2 // | // +--------X--------> // // Modified: // // 12 February 2006 // // Author: // // John Burkardt // // Parameters: // // Input, double T[2*6], the nodal oordinates of the element. // It is common to list these points in counter clockwise order. // // Input, int N, the number of evaluation points. // // Input, double P[2*N], the coordinates of the point where // the basis functions are to be evaluated. // // Output, double PHI[6*N], the value of the basis functions at P. // // Output, double DPHIDX[6*N], DPHIDY[6*N], the value of the X // and Y derivatives of the basis functions at P. // { double gn; double gx; double hn; double hx; int j; for ( j = 0; j < n; j++ ) { // // Basis function 1: PHI(X,Y) = G(3,2) * H(6,4) / normalization. // gx = ( p[0+j*2] - t[0+1*2] ) * ( t[1+2*2] - t[1+1*2] ) - ( t[0+2*2] - t[0+1*2] ) * ( p[1+j*2] - t[1+1*2] ); gn = ( t[0+0*2] - t[0+1*2] ) * ( t[1+2*2] - t[1+1*2] ) - ( t[0+2*2] - t[0+1*2] ) * ( t[1+0*2] - t[1+1*2] ); hx = ( p[0+j*2] - t[0+3*2] ) * ( t[1+5*2] - t[1+3*2] ) - ( t[0+5*2] - t[0+3*2] ) * ( p[1+j*2] - t[1+3*2] ); hn = ( t[0+0*2] - t[0+3*2] ) * ( t[1+5*2] - t[1+3*2] ) - ( t[0+5*2] - t[0+3*2] ) * ( t[1+0*2] - t[1+3*2] ); phi[0+j*6] = ( gx * hx ) / ( gn * hn ); dphidx[0+j*6] = ( ( t[1+2*2] - t[1+1*2] ) * hx + gx * ( t[1+5*2] - t[1+3*2] ) ) / ( gn * hn ); dphidy[0+j*6] = -( ( t[0+2*2] - t[0+1*2] ) * hx + gx * ( t[0+5*2] - t[0+3*2] ) ) / ( gn * hn ); // // Basis function 2: PHI(X,Y) = G(3,1) * H(4,5) / normalization. // gx = ( p[0+j*2] - t[0+0*2] ) * ( t[1+2*2] - t[1+0*2] ) - ( t[0+2*2] - t[0+0*2] ) * ( p[1+j*2] - t[1+0*2] ); gn = ( t[0+1*2] - t[0+0*2] ) * ( t[1+2*2] - t[1+0*2] ) - ( t[0+2*2] - t[0+0*2] ) * ( t[1+1*2] - t[1+0*2] ); hx = ( p[0+j*2] - t[0+4*2] ) * ( t[1+3*2] - t[1+4*2] ) - ( t[0+3*2] - t[0+4*2] ) * ( p[1+j*2] - t[1+4*2] ); hn = ( t[0+1*2] - t[0+4*2] ) * ( t[1+3*2] - t[1+4*2] ) - ( t[0+3*2] - t[0+4*2] ) * ( t[1+1*2] - t[1+4*2] ); phi[1+j*6] = ( gx * hx ) / ( gn * hn ); dphidx[1+j*6] = ( ( t[1+2*2] - t[1+0*2] ) * hx + gx * ( t[1+3*2] - t[1+4*2] ) ) / ( gn * hn ); dphidy[1+j*6] = -( ( t[0+2*2] - t[0+0*2] ) * hx + gx * ( t[0+3*2] - t[0+4*2] ) ) / ( gn * hn ); // // Basis function 3: PHI(X,Y) = G(1,2) * H(5,6) / normalization. // gx = ( p[0+j*2] - t[0+1*2] ) * ( t[1+0*2] - t[1+1*2] ) - ( t[0+0*2] - t[0+1*2] ) * ( p[1+j*2] - t[1+1*2] ); gn = ( t[0+2*2] - t[0+1*2] ) * ( t[1+0*2] - t[1+1*2] ) - ( t[0+0*2] - t[0+1*2] ) * ( t[1+2*2] - t[1+1*2] ); hx = ( p[0+j*2] - t[0+5*2] ) * ( t[1+4*2] - t[1+5*2] ) - ( t[0+4*2] - t[0+5*2] ) * ( p[1+j*2] - t[1+5*2] ); hn = ( t[0+2*2] - t[0+5*2] ) * ( t[1+4*2] - t[1+5*2] ) - ( t[0+4*2] - t[0+5*2] ) * ( t[1+2*2] - t[1+5*2] ); phi[2+j*6] = ( gx * hx ) / ( gn * hn ); dphidx[2+j*6] = ( ( t[1+0*2] - t[1+1*2] ) * hx + gx * ( t[1+4*2] - t[1+5*2] ) ) / ( gn * hn ); dphidy[2+j*6] = -( ( t[0+0*2] - t[0+1*2] ) * hx + gx * ( t[0+4*2] - t[0+5*2] ) ) / ( gn * hn ); // // Basis function 4: PHI(X,Y) = G(1,3) * H(2,3) / normalization. // gx = ( p[0+j*2] - t[0+2*2] ) * ( t[1+0*2] - t[1+2*2] ) - ( t[0+0*2] - t[0+2*2] ) * ( p[1+j*2] - t[1+2*2] ); gn = ( t[0+3*2] - t[0+2*2] ) * ( t[1+0*2] - t[1+2*2] ) - ( t[0+0*2] - t[0+2*2] ) * ( t[1+3*2] - t[1+2*2] ); hx = ( p[0+j*2] - t[0+2*2] ) * ( t[1+1*2] - t[1+2*2] ) - ( t[0+1*2] - t[0+2*2] ) * ( p[1+j*2] - t[1+2*2] ); hn = ( t[0+3*2] - t[0+2*2] ) * ( t[1+1*2] - t[1+2*2] ) - ( t[0+1*2] - t[0+2*2] ) * ( t[1+3*2] - t[1+2*2] ); phi[3+j*6] = ( gx * hx ) / ( gn * hn ); dphidx[3+j*6] = ( ( t[1+0*2] - t[1+2*2] ) * hx + gx * ( t[1+1*2] - t[1+2*2] ) ) / ( gn * hn ); dphidy[3+j*6] = -( ( t[0+0*2] - t[0+2*2] ) * hx + gx * ( t[0+1*2] - t[0+2*2] ) ) / ( gn * hn ); // // Basis function 5: PHI(X,Y) = G(2,1) * H(3,1) / normalization. // gx = ( p[0+j*2] - t[0+0*2] ) * ( t[1+1*2] - t[1+0*2] ) - ( t[0+1*2] - t[0+0*2] ) * ( p[1+j*2] - t[1+0*2] ); gn = ( t[0+4*2] - t[0+0*2] ) * ( t[1+1*2] - t[1+0*2] ) - ( t[0+1*2] - t[0+0*2] ) * ( t[1+4*2] - t[1+0*2] ); hx = ( p[0+j*2] - t[0+0*2] ) * ( t[1+2*2] - t[1+0*2] ) - ( t[0+2*2] - t[0+0*2] ) * ( p[1+j*2] - t[1+0*2] ); hn = ( t[0+4*2] - t[0+0*2] ) * ( t[1+2*2] - t[1+0*2] ) - ( t[0+2*2] - t[0+0*2] ) * ( t[1+4*2] - t[1+0*2] ); phi[4+j*6] = ( gx * hx ) / ( gn * hn ); dphidx[4+j*6] = ( ( t[1+1*2] - t[1+0*2] ) * hx + gx * ( t[1+2*2] - t[1+0*2] ) ) / ( gn * hn ); dphidy[4+j*6] = -( ( t[0+1*2] - t[0+0*2] ) * hx + gx * ( t[0+2*2] - t[0+0*2] ) ) / ( gn * hn ); // // Basis function 6: PHI(X,Y) = G(1,2) * H(3,2) / normalization. // gx = ( p[0+j*2] - t[0+1*2] ) * ( t[1+0*2] - t[1+1*2] ) - ( t[0+0*2] - t[0+1*2] ) * ( p[1+j*2] - t[1+1*2] ); gn = ( t[0+5*2] - t[0+1*2] ) * ( t[1+0*2] - t[1+1*2] ) - ( t[0+0*2] - t[0+1*2] ) * ( t[1+5*2] - t[1+1*2] ); hx = ( p[0+j*2] - t[0+1*2] ) * ( t[1+2*2] - t[1+1*2] ) - ( t[0+2*2] - t[0+1*2] ) * ( p[1+j*2] - t[1+1*2] ); hn = ( t[0+5*2] - t[0+1*2] ) * ( t[1+2*2] - t[1+1*2] ) - ( t[0+2*2] - t[0+1*2] ) * ( t[1+5*2] - t[1+1*2] ); phi[5+j*6] = ( gx * hx ) / ( gn * hn ); dphidx[5+j*6] = ( ( t[1+0*2] - t[1+1*2] ) * hx + gx * ( t[1+2*2] - t[1+1*2] ) ) / ( gn * hn ); dphidy[5+j*6] = -( ( t[0+0*2] - t[0+1*2] ) * hx + gx * ( t[0+2*2] - t[0+1*2] ) ) / ( gn * hn ); } return; } //****************************************************************************80* char ch_cap ( char ch ) //****************************************************************************80* // // Purpose: // // CH_CAP capitalizes a single character. // // Discussion: // // This routine should be equivalent to the library "toupper" function. // // Modified: // // 19 July 1998 // // Author: // // John Burkardt // // Parameters: // // Input, char CH, the character to capitalize. // // Output, char CH_CAP, the capitalized character. // { if ( 97 <= ch && ch <= 122 ) { ch = ch - 32; } return ch; } //****************************************************************************80 bool ch_eqi ( char ch1, char ch2 ) //****************************************************************************80 // // Purpose: // // CH_EQI is true if two characters are equal, disregarding case. // // Modified: // // 13 June 2003 // // Author: // // John Burkardt // // Parameters: // // Input, char CH1, CH2, the characters to compare. // // Output, bool CH_EQI, is true if the two characters are equal, // disregarding case. // { if ( 97 <= ch1 && ch1 <= 122 ) { ch1 = ch1 - 32; } if ( 97 <= ch2 && ch2 <= 122 ) { ch2 = ch2 - 32; } return ( ch1 == ch2 ); } //****************************************************************************80 int ch_to_digit ( char ch ) //****************************************************************************80 // // Purpose: // // CH_TO_DIGIT returns the integer value of a base 10 digit. // // Example: // // CH DIGIT // --- ----- // '0' 0 // '1' 1 // ... ... // '9' 9 // ' ' 0 // 'X' -1 // // Modified: // // 13 June 2003 // // Author: // // John Burkardt // // Parameters: // // Input, char CH, the decimal digit, '0' through '9' or blank are legal. // // Output, int CH_TO_DIGIT, the corresponding integer value. If the character was // 'illegal', then DIGIT is -1. // { int digit; if ( '0' <= ch && ch <= '9' ) { digit = ch - '0'; } else if ( ch == ' ' ) { digit = 0; } else { digit = -1; } return digit; } //****************************************************************************80 double *dtable_data_read ( char *input_filename, int m, int n ) //****************************************************************************80 // // Purpose: // // DTABLE_DATA_READ reads the data from a DTABLE file. // // Discussion: // // The file is assumed to contain one record per line. // // Records beginning with the '#' character 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. // // Modified: // // 27 January 2005 // // Author: // // John Burkardt // // Parameters: // // Input, char *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 DTABLE_DATA_READ[M*N], the table data. // { bool error; ifstream input; int i; int j; char line[255]; double *table; double *x; input.open ( input_filename ); if ( !input ) { cout << "\n"; cout << "DTABLE_DATA_READ - Fatal error!\n"; cout << " 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 ) { input.getline ( line, sizeof ( 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 dtable_header_read ( char *input_filename, int *m, int *n ) //****************************************************************************80 // // Purpose: // // DTABLE_HEADER_READ reads the header from a DTABLE file. // // Modified: // // 04 June 2004 // // Author: // // John Burkardt // // Parameters: // // Input, char *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 ) { cout << "\n"; cout << "DTABLE_HEADER_READ - Fatal error!\n"; cout << " FILE_COLUMN_COUNT failed.\n"; *n = -1; return; } *n = file_row_count ( input_filename ); if ( *n <= 0 ) { cout << "\n"; cout << "DTABLE_HEADER_READ - Fatal error!\n"; cout << " FILE_ROW_COUNT failed.\n"; return; } return; } //****************************************************************************80 int file_column_count ( char *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. // // Modified: // // 13 June 2003 // // Author: // // John Burkardt // // Parameters: // // Input, char *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[256]; // // Open the file. // input.open ( input_filename ); if ( !input ) { column_num = -1; cout << "\n"; cout << "FILE_COLUMN_COUNT - Fatal error!\n"; cout << " Could not open the file:\n"; cout << " \"" << 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 ); 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 ) { cout << "\n"; cout << "FILE_COLUMN_COUNT - Warning!\n"; cout << " 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 ( char *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 '#'. // // Modified: // // 13 June 2003 // // Author: // // John Burkardt // // Parameters: // // Input, char *INPUT_FILENAME, the name of the input file. // // Output, int FILE_ROW_COUNT, the number of rows found. // { int bad_num; int comment_num; ifstream input; int i; char line[100]; int record_num; int row_num; row_num = 0; comment_num = 0; record_num = 0; bad_num = 0; input.open ( input_filename ); if ( !input ) { cout << "\n"; cout << "FILE_ROW_COUNT - Fatal error!\n"; cout << " 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; } //******************************************************************** unsigned long get_seed ( void ) //******************************************************************** // // Purpose: // // GET_SEED returns a random seed for the random number generator. // // Modified: // // 17 September 2003 // // Author: // // John Burkardt // // Parameters: // // Output, unsigned long GET_SEED, a random seed value. // { # define UNSIGNED_LONG_MAX 4294967295UL time_t clock; int i; int hours; int minutes; int seconds; struct tm *lt; static unsigned long seed = 0; time_t tloc; // // If the internal seed is 0, generate a value based on the time. // if ( seed == 0 ) { clock = time ( &tloc ); lt = localtime ( &clock ); // // Extract HOURS. // hours = lt->tm_hour; // // In case of 24 hour clocks, shift so that HOURS is between 1 and 12. // if ( 12 < hours ) { hours = hours - 12; } // // Move HOURS to 0, 1, ..., 11 // hours = hours - 1; minutes = lt->tm_min; seconds = lt->tm_sec; seed = seconds + 60 * ( minutes + 60 * hours ); // // We want values in [1,43200], not [0,43199]. // seed = seed + 1; // // Remap SEED from [1,43200] to [1,UNSIGNED_LONG_MAX]. // seed = ( unsigned long ) ( ( ( double ) seed ) * ( ( double ) UNSIGNED_LONG_MAX ) / ( 60.0 * 60.0 * 12.0 ) ); } // // Never use a seed of 0. // if ( seed == 0 ) { seed = 1; } return seed; # undef UNSIGNED_LONG_MAX } //****************************************************************************80 void grid_sample ( int node_num, double node_xy[], int triangle_order, int triangle_num, int triangle_node[], int triangle_neighbor[], double coef_node[], int x_number, int y_number, double x_vec[], double y_vec[], double u_grid[] ) //****************************************************************************80 // // Purpose: // // GRID_SAMPLE samples a (scalar) FE function on a T3 or T6 triangulation. // // Discussion: // // Note that the values of U returned are true values of the underlying // finite element function. They are NOT produced by constructing some // other function that interpolates the data at the finite element nodes // (something which MATLAB's griddata function can easily do.) Instead, // each sampling node is located within one of the associated finite // element triangles, and the finite element function is developed and // evaluated there. // // MATLAB's scattered data interpolation is wonderful, but it cannot // be guaranteed to reproduce the finite element function corresponding // to nodal data. This routine can. // // So if you are using finite elements, then using THIS routine // (but not MATLAB's griddata function), what you see is what you have// // // Modified: // // 13 October 2006 // // Author: // // John Burkardt // // Parameters: // // Input, int NODE_NUM, the number of nodes. // // Input, double NODE_XY[2*NODE_NUM], the coordinates of the nodes. // // Input, int TRIANGLE_ORDER, the order of the triangles, either // 3 or 6. // // Input, int TRIANGLE_NUM, the number of triangles. // // Input, int TRIANGLE_NODE[TRIANGLE_ORDER*TRIANGLE_NUM], the // nodes that make up each triangle. // // Input, int TRIANGLE_NEIGHBOR[3*TRIANGLE_NUM], the index of the // neighboring triangle on each side, or -1 if no neighbor there. // // Input, double COEF_NODE[NODE_NUM], the finite element // coefficient value at each node. // // Input, int X_NUMBER, Y_NUMBER, the number of sampling points in the // X and Y directions. // // Output, double X_VEC[X_NUMBER], Y_VEC[Y_NUMBER], // U_GRID[X_NUMBER*Y_NUMBER], the coordinates of the X and Y sampling // points, and the value of the finite element function U at those // sampling points. // { double *b; double *dbdx; double *dbdy; int edge; int i; double inc; int j; int node; int order; double p_xy[2]; int t; int *t_node; double *t_xy; double x_frac; double x_hi; double x_lo; double y_frac; double y_hi; double y_lo; b = new double[triangle_order]; dbdx = new double[triangle_order]; dbdy = new double[triangle_order]; t_node = new int[triangle_order]; t_xy = new double[2*triangle_order]; // // Generate the X and Y coordinates of the grid. // x_lo = r8_huge ( ); x_hi = -r8_huge ( ); y_lo = r8_huge ( ); y_hi = -r8_huge ( ); for ( node = 0; node < node_num; node++ ) { x_lo = r8_min ( x_lo, node_xy[0+node*2] ); x_hi = r8_max ( x_hi, node_xy[0+node*2] ); y_lo = r8_min ( y_lo, node_xy[1+node*2] ); y_hi = r8_max ( y_hi, node_xy[1+node*2] ); } if ( x_number <= 1 && y_number <= 1 ) { x_vec[0] = 0.5 * ( x_lo + x_hi ); y_vec[0] = 0.5 * ( y_lo + y_hi ); } else if ( y_number < x_number ) { inc = ( x_hi - x_lo ) / ( double ) ( x_number - 1 ); y_frac = 0.5 * ( ( y_hi - y_lo ) - ( double ) ( y_number - 1 ) * inc ); for ( i = 0; i < x_number; i++ ) { x_vec[i] = ( ( double ) ( x_number - i - 1 ) * x_lo + ( double ) ( i ) * x_hi ) / ( double ) ( x_number - 1 ); } for ( j = 0; j < y_number; j++ ) { y_vec[j] = ( ( double ) ( y_number - j - 1 ) * ( y_lo + y_frac ) + ( double ) ( j ) * ( y_hi - y_frac ) ) / ( double ) ( y_number - 1 ); } } else { inc = ( y_hi - y_lo ) / ( double ) ( y_number - 1 ); x_frac = 0.5 * ( ( x_hi - x_lo ) - ( double ) ( x_number - 1 ) * inc ); for ( i = 0; i < x_number; i++ ) { x_vec[i] = ( ( double ) ( x_number - i - 1 ) * ( x_lo + x_frac ) + ( double ) ( i ) * ( x_hi - x_frac ) ) / ( double ) ( x_number - 1 ); } for ( j = 0; j < y_number; j++ ) { y_vec[j] = ( ( double ) ( y_number - j - 1 ) * y_lo + ( double ) ( j ) * y_hi ) / ( double ) ( y_number - 1 ); } } // // Generate points on the uniform grid, find them in the triangulation, // and evaluate. // for ( j = 0; j < y_number; j++ ) { p_xy[1] = y_vec[j]; for ( i = 0; i < x_number; i++ ) { p_xy[0] = x_vec[i]; triangulation_search ( node_num, node_xy, triangle_order, triangle_num, triangle_node, triangle_neighbor, p_xy, &t, &edge ); // // Evaluate U by getting element_node nodes, coefficients, // basis function values. // // NOTE: // If EDGE is negative, then point P_XY is actually outside the triangulation. // For now, we simply evaluate U as though the basis function in the triangle // extended on out to the point. But this is probably unwise in most cases. // The value should probably be set to a flag value, or truncated to the // value at the nearest point on the triangle to P... // // Also note that in the following loop, we assume that TRIANGLE_NODE // is 1-based, and that for convenience we can get away with making // T_NODE 0-based. // for ( order = 0; order < triangle_order; order++ ) { t_node[order] = triangle_node[order+(t-1)*triangle_order] - 1; } for ( order = 0; order < triangle_order; order++ ) { t_xy[0+order*2] = node_xy[0+t_node[order]*2]; t_xy[1+order*2] = node_xy[1+t_node[order]*2]; } if ( triangle_order == 3 ) { basis_mn_t3 ( t_xy, 1, p_xy, b, dbdx, dbdy ); } else if ( triangle_order == 6 ) { basis_mn_t6 ( t_xy, 1, p_xy, b, dbdx, dbdy ); } u_grid[i+j*x_number] = 0.0; for ( order = 0; order < triangle_order; order++ ) { u_grid[i+j*x_number] = u_grid[i+j*x_number] + coef_node[t_node[order]] * b[order]; } } } delete [] b; delete [] dbdx; delete [] dbdy; delete [] t_node; delete [] t_xy; return; } //****************************************************************************80 void grid_sample_size ( int node_num, double node_xy[], int number, int *x_number, int *y_number ) //****************************************************************************80 // // Purpose: // // GRID_SAMPLE_SIZE determines sizes for a uniform grid. // // Discussion: // // We are given a set of node coordinates. We are going to want // to sample points in the rectangle defined by the range of the // nodes. We want to have NUMBER equally spaced points in the // widest of the two coordinate directions. // // This routine determines X_NUMBER and Y_NUMBER, the number of // nodes to use in the X and Y directions, so that this requirement // is met. Aside from special cases, one of X_NUMBER and Y_NUMBER // will be equal to NUMBER. // // Modified: // // 12 October 2006 // // Author: // // John Burkardt // // Parameters: // // Input, int NODE_NUM, the number of nodes. // // Input, double NODE_XY[2*NODE_NUM], the coordinates of the nodes. // // Input, int NUMBER, the number of sampling points desired in // the widest direction. NUMBER should be at least 1. // // Output, int *NUMBER_X, *NUMBER_Y, the number of sampling points in the // X and Y directions. // { double inc; int node; double x_hi; double x_lo; double y_hi; double y_lo; // // A special case. // if ( number <= 1 ) { *x_number = 1; *y_number = 1; } // // Determine the coordinate ranges for the nodes. // else { x_lo = r8_huge ( ); x_hi = -r8_huge ( ); y_lo = r8_huge ( ); y_hi = -r8_huge ( ); for ( node = 0; node < node_num; node++ ) { x_lo = r8_min ( x_lo, node_xy[0+node*2] ); x_hi = r8_max ( x_hi, node_xy[0+node*2] ); y_lo = r8_min ( y_lo, node_xy[1+node*2] ); y_hi = r8_max ( y_hi, node_xy[1+node*2] ); } // // We want NUMBER nodes in the widest direction. // if ( y_hi - y_lo < x_hi - x_lo ) { *x_number = number; inc = ( x_hi - x_lo ) / ( double ) ( *x_number - 1 ); *y_number = 1 + r8_nint ( ( y_hi - y_lo ) / inc ); } else { *y_number = number; inc = ( y_hi - y_lo ) / ( double ) ( *y_number - 1 ); *x_number = 1 + r8_nint ( ( x_hi - x_lo ) / inc ); } } return; } //****************************************************************************80 int i4_max ( int i1, int i2 ) //****************************************************************************80 // // Purpose: // // I4_MAX returns the maximum of two integers. // // Modified: // // 05 May 2003 // // Author: // // John Burkardt // // Parameters: // // Input, int I1, I2, two integers to be compared. // // Output, int I4_MAX, the larger of I1 and I2. // { if ( i2 < i1 ) { return i1; } else { return i2; } } //****************************************************************************80 int i4_min ( int i1, int i2 ) //****************************************************************************80 // // Purpose: // // I4_MIN returns the smaller of two integers. // // Modified: // // 05 May 2003 // // Author: // // John Burkardt // // Parameters: // // Input, int I1, I2, two integers to be compared. // // Output, int I4_MIN, the smaller of I1 and I2. // { if ( i1 < i2 ) { return i1; } else { return i2; } } //****************************************************************************80 int i4_uniform ( int a, int b, int *seed ) //****************************************************************************80 // // Purpose: // // I4_UNIFORM returns a scaled pseudorandom I4. // // Discussion: // // The pseudorandom number should be uniformly distributed // between A and B. // // Modified: // // 12 November 2006 // // Author: // // John Burkardt // // Reference: // // Paul Bratley, Bennett Fox, Linus Schrage, // A Guide to Simulation, // Springer Verlag, pages 201-202, 1983. // // Pierre L'Ecuyer, // Random Number Generation, // in Handbook of Simulation, // edited by Jerry Banks, // Wiley Interscience, page 95, 1998. // // Bennett Fox, // Algorithm 647: // Implementation and Relative Efficiency of Quasirandom // Sequence Generators, // ACM Transactions on Mathematical Software, // Volume 12, Number 4, pages 362-376, 1986. // // Peter Lewis, Allen Goodman, James Miller // A Pseudo-Random Number Generator for the System/360, // IBM Systems Journal, // Volume 8, pages 136-143, 1969. // // Parameters: // // Input, int A, B, the limits of the interval. // // Input/output, int *SEED, the "seed" value, which should NOT be 0. // On output, SEED has been updated. // // Output, int I4_UNIFORM, a number between A and B. // { int k; float r; int value; if ( *seed == 0 ) { cerr << "\n"; cerr << "I4_UNIFORM - Fatal error!\n"; cerr << " Input value of SEED = 0.\n"; exit ( 1 ); } k = *seed / 127773; *seed = 16807 * ( *seed - k * 127773 ) - k * 2836; if ( *seed < 0 ) { *seed = *seed + 2147483647; } r = ( float ) ( *seed ) * 4.656612875E-10; // // Scale R to lie between A-0.5 and B+0.5. // r = ( 1.0 - r ) * ( ( float ) ( i4_min ( a, b ) ) - 0.5 ) + r * ( ( float ) ( i4_max ( a, b ) ) + 0.5 ); // // Use rounding to convert R to an integer between A and B. // value = r4_nint ( r ); value = i4_max ( value, i4_min ( a, b ) ); value = i4_min ( value, i4_max ( a, b ) ); return value; } //****************************************************************************80 int i4row_compare ( int m, int n, int a[], int i, int j ) //****************************************************************************80 // // Purpose: // // I4ROW_COMPARE compares two rows of a integer array. // // Discussion: // // The two dimensional information is stored in a one dimensional array, // by columns. The entry A(I,J) is stored in A[I+J*M]. // // The input arguments I and J are row indices. They DO NOT use the // C convention of starting at 0, but rather start at 1. // // Example: // // Input: // // M = 3, N = 4, I = 2, J = 3 // // A = ( // 1 2 3 4 // 5 6 7 8 // 9 10 11 12 ) // // Output: // // I4ROW_COMPARE = -1 // // Modified: // // 13 October 2006 // // Author: // // John Burkardt // // Parameters: // // Input, int M, N, the number of rows and columns. // // Input, int A[M*N], the array of data. // // Input, int I, J, the rows to be compared. // I and J must be between 1 and M. // // Output, int I4ROW_COMPARE, the results of the comparison: // -1, row I < row J, // 0, row I = row J, // +1, row J < row I. // { int k; // // Check that I and J are legal. // if ( i < 1 ) { cout << "\n"; cout << "I4ROW_COMPARE - Fatal error!\n"; cout << " Row index I is less than 1.\n"; cout << " I = " << i << "\n"; exit ( 1 ); } else if ( m < i ) { cout << "\n"; cout << "I4ROW_COMPARE - Fatal error!\n"; cout << " Row index I is out of bounds.\n"; cout << " I = " << i << "\n"; cout << " Maximum legal value is M = " << m << "\n"; exit ( 1 ); } if ( j < 1 ) { cout << "\n"; cout << "I4ROW_COMPARE - Fatal error!\n"; cout << " Row index J is less than 1.\n"; cout << " J = " << j << "\n"; exit ( 1 ); } else if ( m < j ) { cout << "\n"; cout << "I4ROW_COMPARE - Fatal error!\n"; cout << " Row index J is out of bounds.\n"; cout << " J = " << j << "\n"; cout << " Maximum legal value is M = " << m << "\n"; exit ( 1 ); } if ( i == j ) { return 0; } for ( k = 0; k < n; k++ ) { if ( a[(i-1)+k*m] < a[(j-1)+k*m] ) { return -1; } else if ( a[(j-1)+k*m] < a[(i-1)+k*m] ) { return +1; } } return 0; } //****************************************************************************80 void i4row_sort_a ( int m, int n, int a[] ) //****************************************************************************80 // // Purpose: // // I4ROW_SORT_A ascending sorts the rows of an I4ROW. // // Definition: // // In lexicographic order, the statement "X < Y", applied to two // vectors X and Y of length M, means that there is some index I, with // 1 <= I <= M, with the property that // // X(J) = Y(J) for J < I, // and // X(I) < Y(I). // // In other words, X is less than Y if, at the first index where they // differ, the X value is less than the Y value. // // Example: // // Input: // // M = 5, N = 3 // // A = // 3 2 1 // 2 4 3 // 3 1 8 // 2 4 2 // 1 9 9 // // Output: // // A = // 1 9 9 // 2 4 2 // 2 4 3 // 3 1 8 // 3 2 1 // // Modified: // // 17 September 2003 // // Author: // // John Burkardt // // Parameters: // // Input, int M, the number of rows of A. // // Input, int N, the number of columns of A. // // Input/output, int A[M*N]. // On input, the array of M rows of N-vectors. // On output, the rows of A have been sorted in ascending // lexicographic order. // { int i; int indx; int isgn; int j; // // Initialize. // i = 0; indx = 0; isgn = 0; j = 0; // // Call the external heap sorter. // for ( ; ; ) { sort_heap_external ( m, &indx, &i, &j, isgn ); // // Interchange the I and J objects. // if ( 0 < indx ) { i4row_swap ( m, n, a, i, j ); } // // Compare the I and J objects. // else if ( indx < 0 ) { isgn = i4row_compare ( m, n, a, i, j ); } else if ( indx == 0 ) { break; } } return; } //****************************************************************************80 void i4row_swap ( int m, int n, int a[], int irow1, int irow2 ) //****************************************************************************80 // // Purpose: // // I4ROW_SWAP swaps two rows of an I4ROW. // // Discussion: // // The two dimensional information is stored as a one dimensional // array, by columns. // // The row indices are 1 based, NOT 0 based! However, a preprocessor // variable, called OFFSET, can be reset from 1 to 0 if you wish to // use 0-based indices. // // Modified: // // 16 September 2003 // // Author: // // John Burkardt // // Parameters: // // Input, int M, N, the number of rows and columns. // // Input/output, int A[M*N], an array of data. // // Input, int IROW1, IROW2, the two rows to swap. // These indices should be between 1 and M. // { # define OFFSET 1 int j; int t; // // Check. // if ( irow1 < 0+OFFSET || m-1+OFFSET < irow1 ) { cout << "\n"; cout << "I4ROW_SWAP - Fatal error!\n"; cout << " IROW1 is out of range.\n"; exit ( 1 ); } if ( irow2 < 0+OFFSET || m-1+OFFSET < irow2 ) { cout << "\n"; cout << "I4ROW_SWAP - Fatal error!\n"; cout << " IROW2 is out of range.\n"; exit ( 1 ); } if ( irow1 == irow2 ) { return; } for ( j = 0; j < n; j++ ) { t = a[irow1-OFFSET+j*m]; a[irow1-OFFSET+j*m] = a[irow2-OFFSET+j*m]; a[irow2-OFFSET+j*m] = t; } return; # undef OFFSET } //****************************************************************************80 int *itable_data_read ( char *input_filename, int m, int n ) //****************************************************************************80 // // Purpose: // // ITABLE_DATA_READ reads data from an ITABLE file. // // Discussion: // // The file is assumed to contain one record per line. // // Records beginning with the '#' character 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. // // Modified: // // 11 October 2003 // // Author: // // John Burkardt // // Parameters: // // Input, char *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, int ITABLE_DATA_READ[M*N], the table data. // { bool error; ifstream input; int i; int j; char line[255]; int *table; int *x; input.open ( input_filename ); if ( !input ) { cout << "\n"; cout << "ITABLE_DATA_READ - Fatal error!\n"; cout << " Could not open the input file: \"" << input_filename << "\"\n"; return NULL; } table = new int[m*n]; x = new int[m]; j = 0; while ( j < n ) { input.getline ( line, sizeof ( line ) ); if ( input.eof ( ) ) { break; } if ( line[0] == '#' || s_len_trim ( line ) == 0 ) { continue; } error = s_to_i4vec ( 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 itable_header_read ( char *input_filename, int *m, int *n ) //****************************************************************************80 // // Purpose: // // ITABLE_HEADER_READ reads the header from an ITABLE file. // // Modified: // // 04 June 2004 // // Author: // // John Burkardt // // Parameters: // // Input, char *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 ) { cout << "\n"; cout << "ITABLE_HEADER_READ - Fatal error!\n"; cout << " FILE_COLUMN_COUNT failed.\n"; *n = -1; return; } *n = file_row_count ( input_filename ); if ( *n <= 0 ) { cout << "\n"; cout << "ITABLE_HEADER_READ - Fatal error!\n"; cout << " FILE_ROW_COUNT failed.\n"; return; } return; } //********************************************************************* float r4_abs ( float x ) //********************************************************************* // // Purpose: // // R4_ABS returns the absolute value of an R4. // // Modified: // // 01 December 2006 // // Author: // // John Burkardt // // Parameters: // // Input, float X, the quantity whose absolute value is desired. // // Output, float R4_ABS, the absolute value of X. // { float value; if ( 0.0 <= x ) { value = x; } else { value = -x; } return value; } //****************************************************************************80 int r4_nint ( float x ) //****************************************************************************80 // // Purpose: // // R4_NINT returns the nearest integer to an R4. // // Examples: // // X R4_NINT // // 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 // // Modified: // // 14 November 2006 // // Author: // // John Burkardt // // Parameters: // // Input, float X, the value. // // Output, int R4_NINT, the nearest integer to X. // { int value; if ( x < 0.0 ) { value = - ( int ) ( r4_abs ( x ) + 0.5 ); } else { value = ( int ) ( r4_abs ( x ) + 0.5 ); } return value; } //****************************************************************************80 float r4_uniform_01 ( int *seed ) //****************************************************************************80 // // Purpose: // // R4_UNIFORM_01 returns a unit pseudorandom R4. // // Discussion: // // This routine implements the recursion // // seed = 16807 * seed mod ( 2**31 - 1 ) // r4_uniform_01 = seed / ( 2**31 - 1 ) // // The integer arithmetic never requires more than 32 bits, // including a sign bit. // // If the initial seed is 12345, then the first three computations are // // Input Output R4_UNIFORM_01 // SEED SEED // // 12345 207482415 0.096616 // 207482415 1790989824 0.833995 // 1790989824 2035175616 0.947702 // // Modified: // // 16 November 2004 // // Author: // // John Burkardt // // Reference: // // Paul Bratley, Bennett Fox, Linus Schrage, // A Guide to Simulation, // Springer Verlag, pages 201-202, 1983. // // Pierre L'Ecuyer, // Random Number Generation, // in Handbook of Simulation // edited by Jerry Banks, // Wiley Interscience, page 95, 1998. // // Bennett Fox, // Algorithm 647: // Implementation and Relative Efficiency of Quasirandom // Sequence Generators, // ACM Transactions on Mathematical Software, // Volume 12, Number 4, pages 362-376, 1986. // // Peter Lewis, Allen Goodman, James Miller, // A Pseudo-Random Number Generator for the System/360, // IBM Systems Journal, // Volume 8, pages 136-143, 1969. // // Parameters: // // Input/output, int *SEED, the "seed" value. Normally, this // value should not be 0. On output, SEED has been updated. // // Output, float R4_UNIFORM_01, a new pseudorandom variate, strictly between // 0 and 1. // { int k; float r; if ( *seed == 0 ) { cerr << "\n"; cerr << "R4_UNIFORM_01 - Fatal error!\n"; cerr << " Input value of SEED = 0.\n"; exit ( 1 ); } k = *seed / 127773; *seed = 16807 * ( *seed - k * 127773 ) - k * 2836; if ( *seed < 0 ) { *seed = *seed + 2147483647; } // // Although SEED can be represented exactly as a 32 bit integer, // it generally cannot be represented exactly as a 32 bit real number! // r = ( float ) ( *seed ) * 4.656612875E-10; return r; } //****************************************************************************80 double r8_huge ( void ) //****************************************************************************80 // // Purpose: // // R8_HUGE returns a "huge" R8. // // Discussion: // // HUGE_VAL is the largest representable legal double precision number, // and is usually defined in math.h, or sometimes in stdlib.h. // // Modified: // // 31 August 2004 // // Author: // // John Burkardt // // Parameters: // // Output, double R8_HUGE, a "huge" R8 value. // { return HUGE_VAL; } //********************************************************************* double r8_max ( double x, double y ) //********************************************************************* // // Purpose: // // R8_MAX returns the maximum of two double precision values. // // 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. // { if ( y < x ) { return x; } else { return y; } } //********************************************************************* double r8_min ( double x, double y ) //********************************************************************* // // Purpose: // // R8_MIN returns the minimum of two double precision values. // // 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. // { if ( y < x ) { return y; } else { return x; } } //****************************************************************************80 int r8_nint ( double x ) //****************************************************************************80 // // Purpose: // // R8_NINT returns the nearest integer to an R8. // // Examples: // // 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 // // Modified: // // 26 August 2004 // // Author: // // John Burkardt // // Parameters: // // Input, double X, the value. // // Output, int R8_NINT, the nearest integer to X. // { int s; int value; if ( x < 0.0 ) { s = -1; } else { s = 1; } value = s * ( int ) ( fabs ( x ) + 0.5 ); return value; } //****************************************************************************80 int s_len_trim ( char *s ) //****************************************************************************80 // // Purpose: // // S_LEN_TRIM returns the length of a string to the last nonblank. // // Modified: // // 26 April 2003 // // Author: // // John Burkardt // // Parameters: // // Input, char *S, a pointer to 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; char *t; n = strlen ( s ); t = s + strlen ( s ) - 1; while ( 0 < n ) { if ( *t != ' ' ) { return n; } t--; n--; } return n; } //****************************************************************************80 int s_to_i4 ( char *s, int *last, bool *error ) //****************************************************************************80 // // Purpose: // // S_TO_I4 reads an integer value from a string. // // Modified: // // 13 June 2003 // // Author: // // John Burkardt // // Parameters: // // Input, char *S, a string to be examined. // // Output, int *LAST, the last character of S used to make IVAL. // // Output, bool *ERROR is TRUE if an error occurred. // // Output, int *S_TO_I4, the integer value read from the string. // If the string is blank, then IVAL will be returned 0. // { char c; int i; int isgn; int istate; int ival; *error = false; istate = 0; isgn = 1; i = 0; ival = 0; while ( *s ) { c = s[i]; i = i + 1; // // Haven't read anything. // if ( istate == 0 ) { if ( c == ' ' ) { } else if ( c == '-' ) { istate = 1; isgn = -1; } else if ( c == '+' ) { istate = 1; isgn = + 1; } else if ( '0' <= c && c <= '9' ) { istate = 2; ival = c - '0'; } else { *error = true; return ival; } } // // Have read the sign, expecting digits. // else if ( istate == 1 ) { if ( c == ' ' ) { } else if ( '0' <= c && c <= '9' ) { istate = 2; ival = c - '0'; } else { *error = true; return ival; } } // // Have read at least one digit, expecting more. // else if ( istate == 2 ) { if ( '0' <= c && c <= '9' ) { ival = 10 * (ival) + c - '0'; } else { ival = isgn * ival; *last = i - 1; return ival; } } } // // If we read all the characters in the string, see if we're OK. // if ( istate == 2 ) { ival = isgn * ival; *last = s_len_trim ( s ); } else { *error = true; *last = 0; } return ival; } //****************************************************************************80 bool s_to_i4vec ( char *s, int n, int ivec[] ) //****************************************************************************80 // // Purpose: // // S_TO_I4VEC reads an integer vector from a string. // // Modified: // // 08 September 2003 // // Author: // // John Burkardt // // Parameters: // // Input, char *S, the string to be read. // // Input, int N, the number of values expected. // // Output, int IVEC[N], the values read from the string. // // Output, bool S_TO_I4VEC, is TRUE if an error occurred. // { bool error; int i; int lchar; double x; error = false; for ( i = 0; i < n; i++ ) { ivec[i] = s_to_i4 ( s, &lchar, &error ); if ( error ) { cout << "\n"; cout << "S_TO_I4VEC - Fatal error!\n"; cout << " S_TO_I4 returned error while reading item " << i << "\n"; return error; } s = s + lchar; } return error; } //****************************************************************************80 double s_to_r8 ( char *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. // // Examples: // // 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) // // Modified: // // 07 August 2003 // // Author: // // John Burkardt // // Parameters: // // Input, char *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 ( char *s, int n, double rvec[] ) //****************************************************************************80 // // Purpose: // // S_TO_R8VEC reads an R8VEC from a string. // // Modified: // // 19 February 2001 // // Author: // // John Burkardt // // Parameters: // // Input, char *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. // { bool error; int i; int lchar; double x; for ( i = 0; i < n; i++ ) { rvec[i] = s_to_r8 ( s, &lchar, &error ); if ( error ) { return error; } s = s + lchar; } return error; } //****************************************************************************80 int s_word_count ( char *s ) //****************************************************************************80 // // Purpose: // // S_WORD_COUNT counts the number of "words" in a string. // // Modified: // // 08 February 2005 // // Author: // // John Burkardt // // Parameters: // // Input, char *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 i; int nword; nword = 0; blank = true; while ( *s ) { if ( *s == ' ' ) { blank = true; } else if ( blank ) { nword = nword + 1; blank = false; } *s++; } return nword; } //******************************************************************** void sort_heap_external ( int n, int *indx, int *i, int *j, int isgn ) //******************************************************************** // // Purpose: // // SORT_HEAP_EXTERNAL externally sorts a list of items into ascending order. // // Discussion: // // The actual list is not passed to the routine. Hence it may // consist of integers, reals, numbers, names, etc. The user, // after each return from the routine, will be asked to compare or // interchange two items. // // The current version of this code mimics the FORTRAN version, // so the values of I and J, in particular, are FORTRAN indices. // // Modified: // // 05 February 2004 // // Reference: // // Albert Nijenhuis, Herbert Wilf, // Combinatorial Algorithms, // Academic Press, 1978, second edition, // ISBN 0-12-519260-6. // // Parameters: // // Input, int N, the length of the input list. // // Input/output, int *INDX. // The user must set INDX to 0 before the first call. // On return, // if INDX is greater than 0, the user must interchange // items I and J and recall the routine. // If INDX is less than 0, the user is to compare items I // and J and return in ISGN a negative value if I is to // precede J, and a positive value otherwise. // If INDX is 0, the sorting is done. // // Output, int *I, *J. On return with INDX positive, // elements I and J of the user's list should be // interchanged. On return with INDX negative, elements I // and J are to be compared by the user. // // Input, int ISGN. On return with INDX negative, the // user should compare elements I and J of the list. If // item I is to precede item J, set ISGN negative, // otherwise set ISGN positive. // { static int i_save = 0; static int j_save = 0; static int k = 0; static int k1 = 0; static int n1 = 0; // // INDX = 0: This is the first call. // if ( *indx == 0 ) { i_save = 0; j_save = 0; k = n / 2; k1 = k; n1 = n; } // // INDX < 0: The user is returning the results of a comparison. // else if ( *indx < 0 ) { if ( *indx == -2 ) { if ( isgn < 0 ) { i_save = i_save + 1; } j_save = k1; k1 = i_save; *indx = -1; *i = i_save; *j = j_save; return; } if ( 0 < isgn ) { *indx = 2; *i = i_save; *j = j_save; return; } if ( k <= 1 ) { if ( n1 == 1 ) { i_save = 0; j_save = 0; *indx = 0; } else { i_save = n1; j_save = 1; n1 = n1 - 1; *indx = 1; } *i = i_save; *j = j_save; return; } k = k - 1; k1 = k; } // // 0 < INDX: the user was asked to make an interchange. // else if ( *indx == 1 ) { k1 = k; } for ( ; ; ) { i_save = 2 * k1; if ( i_save == n1 ) { j_save = k1; k1 = i_save; *indx = -1; *i = i_save; *j = j_save; return; } else if ( i_save <= n1 ) { j_save = i_save + 1; *indx = -2; *i = i_save; *j = j_save; return; } if ( k <= 1 ) { break; } k = k - 1; k1 = k; } if ( n1 == 1 ) { i_save = 0; j_save = 0; *indx = 0; *i = i_save; *j = j_save; } else { i_save = n1; j_save = 1; n1 = n1 - 1; *indx = 1; *i = i_save; *j = j_save; } return; } //********************************************************************** void timestamp ( void ) //********************************************************************** // // Purpose: // // TIMESTAMP prints the current YMDHMS date as a time stamp. // // Example: // // May 31 2001 09:45:54 AM // // Modified: // // 02 October 2003 // // Author: // // John Burkardt // // Parameters: // // None // { # define TIME_SIZE 40 static char time_buffer[TIME_SIZE]; const struct tm *tm; size_t len; time_t now; now = time ( NULL ); tm = localtime ( &now ); len = strftime ( time_buffer, TIME_SIZE, "%d %B %Y %I:%M:%S %p", tm ); cout << time_buffer << "\n"; return; # undef TIME_SIZE } //****************************************************************************80 void triangulation_order3_neighbor_triangles ( int triangle_num, int triangle_node[], int triangle_neighbor[] ) //****************************************************************************80 // // Purpose: // // TRIANGULATION_ORDER3_NEIGHBOR_TRIANGLES determines triangle neighbors. // // Discussion: // // A triangulation of a set of nodes can be completely described by // the coordinates of the nodes, and the list of nodes that make up // each triangle. However, in some cases, it is necessary to know // triangle adjacency information, that is, which triangle, if any, // is adjacent to a given triangle on a particular side. // // This routine creates a data structure recording this information. // // The primary amount of work occurs in sorting a list of 3 * TRIANGLE_NUM // data items. // // Example: // // The input information from TRIANGLE_NODE: // // Triangle Nodes // -------- --------------- // 1 3 4 1 // 2 3 1 2 // 3 3 2 8 // 4 2 1 5 // 5 8 2 13 // 6 8 13 9 // 7 3 8 9 // 8 13 2 5 // 9 9 13 7 // 10 7 13 5 // 11 6 7 5 // 12 9 7 6 // 13 10 9 6 // 14 6 5 12 // 15 11 6 12 // 16 10 6 11 // // The output information in TRIANGLE_NEIGHBOR: // // Triangle Neighboring Triangles // -------- --------------------- // // 1 -1 -1 2 // 2 1 4 3 // 3 2 5 7 // 4 2 -1 8 // 5 3 8 6 // 6 5 9 7 // 7 3 6 -1 // 8 5 4 10 // 9 6 10 12 // 10 9 8 11 // 11 12 10 14 // 12 9 11 13 // 13 -1 12 16 // 14 11 -1 15 // 15 16 14 -1 // 16 13 15 -1 // // Modified: // // 17 October 2003 // // Author: // // John Burkardt // // Parameters: // // Input, int TRIANGLE_NUM, the number of triangles. // // Input, int TRIANGLE_NODE[3*TRIANGLE_NUM], the nodes that make up each triangle. // // Output, int TRIANGLE_NEIGHBOR[3*TRIANGLE_NUM], the three triangles that are direct // neighbors of a given triangle. TRIANGLE_NEIGHBOR(1,I) is the index of the triangle // which touches side 1, defined by nodes 2 and 3, and so on. TRIANGLE_NEIGHBOR(1,I) // is negative if there is no neighbor on that side. In this case, that // side of the triangle lies on the boundary of the triangulation. // { int i; int irow; int j; int k; int *row; int side1; int side2; int tri; int tri1; int tri2; row = new int[3*triangle_num*4]; // // Step 1. // From the list of nodes for triangle T, of the form: (I,J,K) // construct the three neighbor relations: // // (I,J,1,T) or (J,I,1,T), // (J,K,2,T) or (K,J,2,T), // (K,I,3,T) or (I,K,3,T) // // where we choose (I,J,1,T) if I < J, or else (J,I,1,T) // for ( tri = 0; tri < triangle_num; tri++ ) { i = triangle_node[0+tri*3]; j = triangle_node[1+tri*3]; k = triangle_node[2+tri*3]; if ( i < j ) { row[3*tri+0*3*triangle_num] = i; row[3*tri+1*3*triangle_num] = j; row[3*tri+2*3*triangle_num] = 1; row[3*tri+3*3*triangle_num] = tri + 1; } else { row[3*tri+0*3*triangle_num] = j; row[3*tri+1*3*triangle_num] = i; row[3*tri+2*3*triangle_num] = 1; row[3*tri+3*3*triangle_num] = tri + 1; } if ( j < k ) { row[3*tri+1+0*3*triangle_num] = j; row[3*tri+1+1*3*triangle_num] = k; row[3*tri+1+2*3*triangle_num] = 2; row[3*tri+1+3*3*triangle_num] = tri + 1; } else { row[3*tri+1+0*3*triangle_num] = k; row[3*tri+1+1*3*triangle_num] = j; row[3*tri+1+2*3*triangle_num] = 2; row[3*tri+1+3*3*triangle_num] = tri + 1; } if ( k < i ) { row[3*tri+2+0*3*triangle_num] = k; row[3*tri+2+1*3*triangle_num] = i; row[3*tri+2+2*3*triangle_num] = 3; row[3*tri+2+3*3*triangle_num] = tri + 1; } else { row[3*tri+2+0*3*triangle_num] = i; row[3*tri+2+1*3*triangle_num] = k; row[3*tri+2+2*3*triangle_num] = 3; row[3*tri+2+3*3*triangle_num] = tri + 1; } } // // Step 2. Perform an ascending dictionary sort on the neighbor relations. // We only intend to sort on columns 1 and 2; the routine we call here // sorts on columns 1 through 4 but that won't hurt us. // // What we need is to find cases where two triangles share an edge. // Say they share an edge defined by the nodes I and J. Then there are // two rows of ROW that start out ( I, J, ?, ? ). By sorting ROW, // we make sure that these two rows occur consecutively. That will // make it easy to notice that the triangles are neighbors. // i4row_sort_a ( 3*triangle_num, 4, row ); // // Step 3. Neighboring triangles show up as consecutive rows with // identical first two entries. Whenever you spot this happening, // make the appropriate entries in TRIANGLE_NEIGHBOR. // for ( j = 0; j < triangle_num; j++ ) { for ( i = 0; i < 3; i++ ) { triangle_neighbor[i+j*3] = -1; } } irow = 1; for ( ; ; ) { if ( 3 * triangle_num <= irow ) { break; } if ( row[irow-1+0*3*triangle_num] != row[irow+0*3*triangle_num] || row[irow-1+1*3*triangle_num] != row[irow+1*3*triangle_num] ) { irow = irow + 1; continue; } side1 = row[irow-1+2*3*triangle_num]; tri1 = row[irow-1+3*3*triangle_num]; side2 = row[irow +2*3*triangle_num]; tri2 = row[irow +3*3*triangle_num]; triangle_neighbor[side1-1+(tri1-1)*3] = tri2; triangle_neighbor[side2-1+(tri2-1)*3] = tri1; irow = irow + 2; } delete [] row; return; } //****************************************************************************80 void triangulation_order6_neighbor_triangles ( int triangle_num, int triangle_node[], int triangle_neighbor[] ) //****************************************************************************80 // // Purpose: // // TRIANGULATION_ORDER6_NEIGHBOR_TRIANGLES determines triangle neighbors. // // Discussion: // // A triangulation of a set of nodes can be completely described by // the coordinates of the nodes, and the list of nodes that make up // each triangle. However, in some cases, it is necessary to know // triangle adjacency information, that is, which triangle, if any, // is adjacent to a given triangle on a particular side. // // This routine creates a data structure recording this information. // // The primary amount of work occurs in sorting a list of 3 * TRIANGLE_NUM // data items. // // Modified: // // 13 June 2005 // // Author: // // John Burkardt // // Parameters: // // Input, int TRIANGLE_NUM, the number of triangles. // // Input, int TRIANGLE_NODE[6*TRIANGLE_NUM], the nodes that make up each triangle. // // Output, int TRIANGLE_NEIGHBOR[3*TRIANGLE_NUM], the three triangles that are direct // neighbors of a given triangle. TRIANGLE_NEIGHBOR(1,I) is the index of the triangle // which touches side 1, defined by nodes 2 and 3, and so on. TRIANGLE_NEIGHBOR(1,I) // is negative if there is no neighbor on that side. In this case, that // side of the triangle lies on the boundary of the triangulation. // { int i; int irow; int j; int k; int *row; int side1; int side2; int tri; int tri1; int tri2; row = new int[3*triangle_num*4]; // // Step 1. // From the list of nodes for triangle T, of the form: (I,J,K) // construct the three neighbor relations: // // (I,J,1,T) or (J,I,1,T), // (J,K,2,T) or (K,J,2,T), // (K,I,3,T) or (I,K,3,T) // // where we choose (I,J,1,T) if I < J, or else (J,I,1,T) // for ( tri = 0; tri < triangle_num; tri++ ) { i = triangle_node[0+tri*6]; j = triangle_node[1+tri*6]; k = triangle_node[2+tri*6]; if ( i < j ) { row[3*tri+0*3*triangle_num] = i; row[3*tri+1*3*triangle_num] = j; row[3*tri+2*3*triangle_num] = 1; row[3*tri+3*3*triangle_num] = tri + 1; } else { row[3*tri+0*3*triangle_num] = j; row[3*tri+1*3*triangle_num] = i; row[3*tri+2*3*triangle_num] = 1; row[3*tri+3*3*triangle_num] = tri + 1; } if ( j < k ) { row[3*tri+1+0*3*triangle_num] = j; row[3*tri+1+1*3*triangle_num] = k; row[3*tri+1+2*3*triangle_num] = 2; row[3*tri+1+3*3*triangle_num] = tri + 1; } else { row[3*tri+1+0*3*triangle_num] = k; row[3*tri+1+1*3*triangle_num] = j; row[3*tri+1+2*3*triangle_num] = 2; row[3*tri+1+3*3*triangle_num] = tri + 1; } if ( k < i ) { row[3*tri+2+0*3*triangle_num] = k; row[3*tri+2+1*3*triangle_num] = i; row[3*tri+2+2*3*triangle_num] = 3; row[3*tri+2+3*3*triangle_num] = tri + 1; } else { row[3*tri+2+0*3*triangle_num] = i; row[3*tri+2+1*3*triangle_num] = k; row[3*tri+2+2*3*triangle_num] = 3; row[3*tri+2+3*3*triangle_num] = tri + 1; } } // // Step 2. Perform an ascending dictionary sort on the neighbor relations. // We only intend to sort on columns 1 and 2; the routine we call here // sorts on columns 1 through 4 but that won't hurt us. // // What we need is to find cases where two triangles share an edge. // Say they share an edge defined by the nodes I and J. Then there are // two rows of ROW that start out ( I, J, ?, ? ). By sorting ROW, // we make sure that these two rows occur consecutively. That will // make it easy to notice that the triangles are neighbors. // i4row_sort_a ( 3*triangle_num, 4, row ); // // Step 3. Neighboring triangles show up as consecutive rows with // identical first two entries. Whenever you spot this happening, // make the appropriate entries in TRIANGLE_NEIGHBOR. // for ( j = 0; j < triangle_num; j++ ) { for ( i = 0; i < 3; i++ ) { triangle_neighbor[i+j*3] = -1; } } irow = 1; for ( ; ; ) { if ( 3 * triangle_num <= irow ) { break; } if ( row[irow-1+0*3*triangle_num] != row[irow+0*3*triangle_num] || row[irow-1+1*3*triangle_num] != row[irow+1*3*triangle_num] ) { irow = irow + 1; continue; } side1 = row[irow-1+2*3*triangle_num]; tri1 = row[irow-1+3*3*triangle_num]; side2 = row[irow +2*3*triangle_num]; tri2 = row[irow +3*3*triangle_num]; triangle_neighbor[side1-1+(tri1-1)*3] = tri2; triangle_neighbor[side2-1+(tri2-1)*3] = tri1; irow = irow + 2; } delete [] row; return; } //****************************************************************************80 void triangulation_search ( int node_num, double node_xy[], int triangle_order, int triangle_num, int triangle_node[], int triangle_neighbor[], double p[2], int *triangle, int *edge ) //****************************************************************************80 // // Purpose: // // TRIANGULATION_SEARCH searches a triangulation for a point. // // Discussion: // // The algorithm "walks" from one triangle to its neighboring triangle, // and so on, until a triangle is found containing point P, or P is found // to be outside the convex hull. // // The algorithm computes the barycentric coordinates of the point with // respect to the current triangle. If all three quantities are positive, // the point is contained in the triangle. If the I-th coordinate is // negative, then (X,Y) lies on the far side of edge I, which is opposite // from vertex I. This gives a hint as to where to search next. // // For a Delaunay triangulation, the search is guaranteed to terminate. // For other triangulations, a cycle may occur. // // Note the surprising fact that, even for a Delaunay triangulation of // a set of nodes, the nearest point to (X,Y) need not be one of the // vertices of the triangle containing (X,Y). // // The code can be called for triangulations of any order, but only // the first three nodes in each triangle are considered. Thus, if // higher order triangles are used, and the extra nodes are intended // to give the triangle a polygonal shape, these will have no effect, // and the results obtained here might be misleading. // // Modified: // // 27 September 2006 // // Author: // // Barry Joe, // Department of Computing Science, // University of Alberta, // Edmonton, Alberta, Canada T6G 2H1 // // 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 NODE_NUM, the number of nodes. // // Input, double NODE_XY[2*NODE_NUM], the coordinates of the nodes. // // Input, int TRIANGLE_ORDER, the order of the triangles. // // Input, int TRIANGLE_NUM, the number of triangles in the triangulation. // // Input, int TRIANGLE_NODE[TRIANGLE_ORDER*TRIANGLE_NUM], // the nodes of each triangle. // // Input, int TRIANGLE_NEIGHBOR[3*TRIANGLE_NUM], the triangle neighbor list. // // Input, double P[2], the coordinates of a point. // // Output, int *TRIANGLE, the index of the triangle where the search ended. // If a cycle occurred, then TRIANGLE = -1. // // Output, int *EDGE, indicates the position of the point (X,Y) in // triangle TRIANGLE: // 0, the interior or boundary of the triangle; // -1, outside the convex hull of the triangulation, past edge 1; // -2, outside the convex hull of the triangulation, past edge 2; // -3, outside the convex hull of the triangulation, past edge 3. // { int a; double alpha; int b; double beta; int c; int count; double det; double dxp; double dxa; double dxb; double dyp; double dya; double dyb; double gamma; int seed; count = 0; edge = 0; seed = get_seed ( ); *triangle = i4_uniform ( 1, triangle_num, &seed ); for ( ; ; ) { count = count + 1; if ( triangle_num < count ) { cout << "\n"; cout << "TRIANGULATION_SEARCH - Fatal error!\n"; cout << " The algorithm seems to be cycling.\n"; cout << " Current triangle is " << *triangle << "\n"; *triangle = -1; *edge = -1; return; } // // Get the vertices of triangle TRIANGLE. // a = triangle_node[0+(*triangle-1)*triangle_order] - 1; b = triangle_node[1+(*triangle-1)*triangle_order] - 1; c = triangle_node[2+(*triangle-1)*triangle_order] - 1; // // Using vertex C as a base, compute the distances to vertices A and B, // and the point (X,Y). // dxa = node_xy[0+a*2] - node_xy[0+c*2]; dya = node_xy[1+a*2] - node_xy[1+c*2]; dxb = node_xy[0+b*2] - node_xy[0+c*2]; dyb = node_xy[1+b*2] - node_xy[1+c*2]; dxp = p[0] - node_xy[0+c*2]; dyp = p[1] - node_xy[1+c*2]; det = dxa * dyb - dya * dxb; // // Compute the barycentric coordinates of the point (X,Y) with respect // to this triangle. // alpha = ( dxp * dyb - dyp * dxb ) / det; beta = ( dxa * dyp - dya * dxp ) / det; gamma = 1.0 - alpha - beta; // // If the barycentric coordinates are all positive, then the point // is inside the triangle and we're done. // if ( 0.0 <= alpha && 0.0 <= beta && 0.0 <= gamma ) { break; } // // At least one barycentric coordinate is negative. // // If there is a negative barycentric coordinate for which there exists // an opposing triangle neighbor closer to the point, move to that triangle. // // (Two coordinates could be negative, in which case we could go for the // most negative one, or the most negative one normalized by the actual // distance it represents). // if ( alpha < 0.0 && 0 <= triangle_neighbor[1+(*triangle-1)*3] ) { *triangle = triangle_neighbor[1+(*triangle-1)*3]; continue; } else if ( beta < 0.0 && 0 <= triangle_neighbor[2+(*triangle-1)*3] ) { *triangle = triangle_neighbor[2+(*triangle-1)*3]; continue; } else if ( gamma < 0.0 && 0 <= triangle_neighbor[0+(*triangle-1)*3] ) { *triangle = triangle_neighbor[0+(*triangle-1)*3]; continue; } // // All negative barycentric coordinates correspond to vertices opposite // sides on the convex hull. // // Note the edge and exit. // if ( alpha < 0.0 ) { *edge = -2; break; } else if ( beta < 0.0 ) { *edge = -3; break; } else if ( gamma < 0.0 ) { *edge = -1; break; } else { cout << "\n"; cout << "TRIANGULATION_SEARCH - Fatal error!\n"; cout << " The algorithm seems to have reached a dead end\n"; cout << " after " << count << " steps.\n"; *triangle = -1; *edge = -1; return; } } return; }