# include # include # include # include # include # include # include using namespace std; int main ( int argc, char *argv[] ); char ch_cap ( char ch ); bool ch_eqi ( char ch1, char ch2 ); int ch_to_digit ( char ch ); void comp_next ( int n, int k, int a[], bool *more, int *h, int *t ); int file_column_count ( string input_filename ); int file_row_count ( string input_filename ); double *monomial_value ( int dim_num, int point_num, double x[], int expon[] ); double *r8mat_data_read ( string input_filename, int m, int n ); double r8mat_det_4d ( double a[4*4] ); void r8mat_header_read ( string input_filename, int *m, int *n ); double r8vec_dot ( int n, double a1[], double a2[] ); int s_len_trim ( string s ); int s_to_i4 ( string s, int *last, bool *error ); 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 ); double tet01_monomial_integral ( int dim_num, int expon[] ); double tet01_monomial_quadrature ( int dim_num, int expon[], int point_num, double x[], double weight[] ); void tetrahedron_order4_physical_to_reference ( double tetra[], int n, double phy[], double ref[] ); double tetrahedron_volume ( double tetra[3*4] ); void timestamp ( ); //****************************************************************************80 int main ( int argc, char *argv[] ) //****************************************************************************80 // // Purpose: // // MAIN is the main program for TETRAHEDRON_EXACTNESS. // // Discussion: // // This program investigates the polynomial exactness of a quadrature // rule for the tetrahedron. // // Licensing: // // This code is distributed under the GNU LGPL license. // // Modified: // // 07 October 2010 // // Author: // // John Burkardt // { int degree; int degree_max; int dim; int dim_num; int dim_num2; int dim_num3; bool error; int *expon; int h; int last; bool more; int point; int point_num; int point_num2; int point_num3; double quad_error; string quad_filename; string quad_r_filename; string quad_w_filename; string quad_x_filename; double *r; int t; double volume; double *w; double *x; double *x_ref; timestamp ( ); cout << "\n"; cout << "TETRAHEDRON_EXACTNESS\n"; cout << " C++ version\n"; cout << "\n"; cout << " Investigate the polynomial exactness of a quadrature\n"; cout << " rule for the tetrahedron by integrating all monomials\n"; cout << " of a given degree.\n"; cout << "\n"; cout << " The rule will be adjusted to the unit tetrahedron.\n"; // // Get the quadrature file root name: // if ( 1 < argc ) { quad_filename = argv[1]; } else { cout << "\n"; cout << "TETRAHEDRON_EXACTNESS:\n"; cout << " Enter the \"root\" name of the quadrature files.\n"; cin >> quad_filename; } // // Create the names of: // the quadrature X file; // the quadrature W file; // the quadrature R file; // quad_r_filename = quad_filename + "_r.txt"; quad_w_filename = quad_filename + "_w.txt"; quad_x_filename = quad_filename + "_x.txt"; // // The second command line argument is the maximum degree. // if ( 2 < argc ) { degree_max = s_to_i4 ( argv[2], &last, &error ); } else { cout << "\n"; cout << "TETRAHEDRON_EXACTNESS:\n"; cout << " Please enter the maximum total degree to check.\n"; cin >> degree_max; } // // Summarize the input. // cout << "\n"; cout << "TETRAHEDRON_EXACTNESS: User input:\n"; cout << " Quadrature rule X file = \"" << quad_x_filename << "\".\n"; cout << " Quadrature rule W file = \"" << quad_w_filename << "\".\n"; cout << " Quadrature rule R file = \"" << quad_r_filename << "\".\n"; cout << " Maximum total degree to check = " << degree_max << "\n"; // // Read the X file. // r8mat_header_read ( quad_x_filename, &dim_num, &point_num ); cout << "\n"; cout << " Spatial dimension = " << dim_num << "\n"; cout << " Number of points = " << point_num << "\n"; if ( dim_num != 3 ) { cerr << "\n"; cerr << "TETRAHEDRON_EXACTNESS - Fatal error!\n"; cerr << " The quadrature abscissas must be 3 dimensional.\n"; exit ( 1 ); } x = r8mat_data_read ( quad_x_filename, dim_num, point_num ); // // Read the W file. // r8mat_header_read ( quad_w_filename, &dim_num2, &point_num2 ); if ( dim_num2 != 1 ) { cerr << "\n"; cerr << "TETRAHEDRON_EXACTNESS - Fatal error!\n"; cerr << " The quadrature weight file should have exactly\n"; cerr << " one value on each line.\n"; exit ( 1 ); } if ( point_num2 != point_num ) { cerr << "\n"; cerr << "TETRAHEDRON_EXACTNESS - Fatal error!\n"; cerr << " The quadrature weight file should have exactly\n"; cerr << " the same number of lines as the abscissa file.\n"; exit ( 1 ); } w = r8mat_data_read ( quad_w_filename, 1, point_num ); // // Read the R file. // r8mat_header_read ( quad_r_filename, &dim_num3, &point_num3 ); if ( dim_num3 != dim_num ) { cerr << "\n"; cerr << "TETRAHEDRON_EXACTNESS - Fatal error!\n"; cerr << " The quadrature region file should have the same\n"; cerr << " number of values on each line as the abscissa file\n"; cerr << " does.\n"; exit ( 1 ); } if ( point_num3 != 4 ) { cerr << "\n"; cerr << "TETRAHEDRON_EXACTNESS - Fatal error!\n"; cerr << " The quadrature region file should have 4 lines.\n"; exit ( 1 ); } r = r8mat_data_read ( quad_r_filename, dim_num, 4 ); // // Rescale the weights. // volume = tetrahedron_volume ( r ); for ( point = 0; point < point_num; point++ ) { w[point] = ( 1.0 / 6.0 ) * w[point] / volume; } // // Translate the abscissas. // x_ref = new double[dim_num*point_num]; tetrahedron_order4_physical_to_reference ( r, point_num, x, x_ref ); // // Explore the monomials. // expon = new int[dim_num]; cout << "\n"; cout << " Error Degree Exponents\n"; cout << "\n"; for ( degree = 0; degree <= degree_max; degree++ ) { more = false; h = 0; t = 0; for ( ; ; ) { comp_next ( degree, dim_num, expon, &more, &h, &t ); quad_error = tet01_monomial_quadrature ( dim_num, expon, point_num, x_ref, w ); cout << " " << setw(12) << quad_error << " " << setw(2) << degree << " "; for ( dim = 0; dim < dim_num; dim++ ) { cout << setw(3) << expon[dim]; } cout << "\n"; if ( !more ) { break; } } cout << "\n"; } // // Free memory. // delete [] expon; delete [] r; delete [] w; delete [] x; delete [] x_ref; // // Terminate. // cout << "\n"; cout << "TETRAHEDRON_EXACTNESS:\n"; cout << " Normal end of execution.\n"; cout << "\n"; timestamp ( ); return 0; } //****************************************************************************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. // // Licensing: // // This code is distributed under the GNU LGPL license. // // 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. // // Licensing: // // This code is distributed under the GNU LGPL license. // // 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 // // Licensing: // // This code is distributed under the GNU LGPL license. // // 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 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 void comp_next ( int n, int k, int a[], bool *more, int *h, int *t ) //****************************************************************************80 // // Purpose: // // COMP_NEXT computes the compositions of the integer N into K parts. // // Discussion: // // A composition of the integer N into K parts is an ordered sequence // of K nonnegative integers which sum to N. The compositions (1,2,1) // and (1,1,2) are considered to be distinct. // // The routine computes one composition on each call until there are no more. // For instance, one composition of 6 into 3 parts is // 3+2+1, another would be 6+0+0. // // On the first call to this routine, set MORE = FALSE. The routine // will compute the first element in the sequence of compositions, and // return it, as well as setting MORE = TRUE. If more compositions // are desired, call again, and again. Each time, the routine will // return with a new composition. // // However, when the LAST composition in the sequence is computed // and returned, the routine will reset MORE to FALSE, signaling that // the end of the sequence has been reached. // // This routine originally used a SAVE statement to maintain the // variables H and T. I have decided that it is safer // to pass these variables as arguments, even though the user should // never alter them. This allows this routine to safely shuffle // between several ongoing calculations. // // // There are 28 compositions of 6 into three parts. This routine will // produce those compositions in the following order: // // I A // - --------- // 1 6 0 0 // 2 5 1 0 // 3 4 2 0 // 4 3 3 0 // 5 2 4 0 // 6 1 5 0 // 7 0 6 0 // 8 5 0 1 // 9 4 1 1 // 10 3 2 1 // 11 2 3 1 // 12 1 4 1 // 13 0 5 1 // 14 4 0 2 // 15 3 1 2 // 16 2 2 2 // 17 1 3 2 // 18 0 4 2 // 19 3 0 3 // 20 2 1 3 // 21 1 2 3 // 22 0 3 3 // 23 2 0 4 // 24 1 1 4 // 25 0 2 4 // 26 1 0 5 // 27 0 1 5 // 28 0 0 6 // // Licensing: // // This code is distributed under the GNU LGPL license. // // Modified: // // 02 July 2008 // // Author: // // Original FORTRAN77 version by Albert Nijenhuis, Herbert Wilf. // C++ version by John Burkardt. // // Reference: // // Albert Nijenhuis, Herbert Wilf, // Combinatorial Algorithms for Computers and Calculators, // Second Edition, // Academic Press, 1978, // ISBN: 0-12-519260-6, // LC: QA164.N54. // // Parameters: // // Input, int N, the integer whose compositions are desired. // // Input, int K, the number of parts in the composition. // // Input/output, int A[K], the parts of the composition. // // Input/output, bool *MORE. // Set MORE = FALSE on first call. It will be reset to TRUE on return // with a new composition. Each new call returns another composition until // MORE is set to FALSE when the last composition has been computed // and returned. // // Input/output, int *H, *T, two internal parameters needed for the // computation. The user should allocate space for these in the calling // program, include them in the calling sequence, but never alter them! // { int i; if ( !( *more ) ) { *t = n; *h = 0; a[0] = n; for ( i = 1; i < k; i++ ) { a[i] = 0; } } else { if ( 1 < *t ) { *h = 0; } *h = *h + 1; *t = a[*h-1]; a[*h-1] = 0; a[0] = *t - 1; a[*h] = a[*h] + 1; } *more = ( a[k-1] != n ); return; } //****************************************************************************80 int file_column_count ( string filename ) //****************************************************************************80 // // Purpose: // // FILE_COLUMN_COUNT counts the columns in the first line of a file. // // Discussion: // // The file is assumed to be a simple text file. // // Most lines of the file are 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: // // 05 July 2009 // // Author: // // John Burkardt // // Parameters: // // Input, string 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; string text; // // Open the file. // input.open ( 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 << " \"" << filename << "\"\n"; exit ( 1 ); } // // Read one line, but skip blank lines and comment lines. // got_one = false; for ( ; ; ) { getline ( input, text ); if ( input.eof ( ) ) { break; } if ( s_len_trim ( text ) <= 0 ) { continue; } if ( text[0] == '#' ) { continue; } got_one = true; break; } if ( !got_one ) { input.close ( ); input.open ( filename.c_str ( ) ); for ( ; ; ) { input >> text; if ( input.eof ( ) ) { break; } if ( s_len_trim ( text ) == 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 ( text ); 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 bad_num; int comment_num; ifstream input; int i; string line; int record_num; int row_num; row_num = 0; comment_num = 0; record_num = 0; bad_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"; exit ( 1 ); } for ( ; ; ) { getline ( input, 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 double *monomial_value ( int dim_num, int point_num, double x[], int expon[] ) //****************************************************************************80 // // Purpose: // // MONOMIAL_VALUE evaluates a monomial. // // Discussion: // // This routine evaluates a monomial of the form // // product ( 1 <= dim <= dim_num ) x(dim)^expon(dim) // // where the exponents are nonnegative integers. Note that // if the combination 0^0 is encountered, it should be treated // as 1. // // Licensing: // // This code is distributed under the GNU LGPL license. // // Modified: // // 05 May 2007 // // Author: // // John Burkardt // // Parameters: // // Input, int DIM_NUM, the spatial dimension. // // Input, int POINT_NUM, the number of points at which the // monomial is to be evaluated. // // Input, double X[DIM_NUM*POINT_NUM], the point coordinates. // // Input, int EXPON[DIM_NUM], the exponents. // // Output, double MONOMIAL_VALUE[POINT_NUM], the value of the monomial. // { int dim; int point; double *value; value = new double[point_num]; for ( point = 0; point < point_num; point++ ) { value[point] = 1.0; } for ( dim = 0; dim < dim_num; dim++ ) { if ( 0 != expon[dim] ) { for ( point = 0; point < point_num; point++ ) { value[point] = value[point] * pow ( x[dim+point*dim_num], expon[dim] ); } } } return value; } //****************************************************************************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"; exit ( 1 ); } 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 double r8mat_det_4d ( double a[4*4] ) //****************************************************************************80 // // Purpose: // // R8MAT_DET_4D computes the determinant of a 4 by 4 R8MAT. // // Discussion: // // The two dimensional array is stored as a one dimensional vector, // by COLUMNS. // // Licensing: // // This code is distributed under the GNU LGPL license. // // Modified: // // 10 September 2003 // // Author: // // John Burkardt // // Parameters: // // Input, double A[4*4], the matrix whose determinant is desired. // // Output, double R8MAT_DET_4D, the determinant of the matrix. // { double det; det = a[0+0*4] * ( a[1+1*4] * ( a[2+2*4] * a[3+3*4] - a[2+3*4] * a[3+2*4] ) - a[1+2*4] * ( a[2+1*4] * a[3+3*4] - a[2+3*4] * a[3+1*4] ) + a[1+3*4] * ( a[2+1*4] * a[3+2*4] - a[2+2*4] * a[3+1*4] ) ) - a[0+1*4] * ( a[1+0*4] * ( a[2+2*4] * a[3+3*4] - a[2+3*4] * a[3+2*4] ) - a[1+2*4] * ( a[2+0*4] * a[3+3*4] - a[2+3*4] * a[3+0*4] ) + a[1+3*4] * ( a[2+0*4] * a[3+2*4] - a[2+2*4] * a[3+0*4] ) ) + a[0+2*4] * ( a[1+0*4] * ( a[2+1*4] * a[3+3*4] - a[2+3*4] * a[3+1*4] ) - a[1+1*4] * ( a[2+0*4] * a[3+3*4] - a[2+3*4] * a[3+0*4] ) + a[1+3*4] * ( a[2+0*4] * a[3+1*4] - a[2+1*4] * a[3+0*4] ) ) - a[0+3*4] * ( a[1+0*4] * ( a[2+1*4] * a[3+2*4] - a[2+2*4] * a[3+1*4] ) - a[1+1*4] * ( a[2+0*4] * a[3+2*4] - a[2+2*4] * a[3+0*4] ) + a[1+2*4] * ( a[2+0*4] * a[3+1*4] - a[2+1*4] * a[3+0*4] ) ); return det; } //****************************************************************************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"; exit ( 1 ); } *n = file_row_count ( input_filename ); if ( *n <= 0 ) { cerr << "\n"; cerr << "R8MAT_HEADER_READ - Fatal error!\n"; cerr << " FILE_ROW_COUNT failed.\n"; exit ( 1 ); } return; } //****************************************************************************80 double r8vec_dot ( int n, double a1[], double a2[] ) //****************************************************************************80 // // Purpose: // // R8VEC_DOT computes the dot product of a pair of R8VEC's. // // Licensing: // // This code is distributed under the GNU LGPL license. // // Modified: // // 03 July 2005 // // Author: // // John Burkardt // // Parameters: // // Input, int N, the number of entries in the vectors. // // Input, double A1[N], A2[N], the two vectors to be considered. // // Output, double R8VEC_DOT, the dot product of the vectors. // { int i; double value; value = 0.0; for ( i = 0; i < n; i++ ) { value = value + a1[i] * a2[i]; } return value; } //****************************************************************************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 int s_to_i4 ( string s, int *last, bool *error ) //****************************************************************************80 // // Purpose: // // S_TO_I4 reads an I4 from a string. // // Licensing: // // This code is distributed under the GNU LGPL license. // // Modified: // // 05 July 2009 // // Author: // // John Burkardt // // Parameters: // // Input, string 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; for ( ; ; ) { 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 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 double tet01_monomial_integral ( int dim_num, int expon[] ) //****************************************************************************80 // // Purpose: // // TET01_MONOMIAL_INTEGRAL integrates a monomial over the unit tetrahedron. // // Discussion: // // This routine evaluates a monomial of the form // // product ( 1 <= dim <= dim_num ) x(dim)^expon(dim) // // where the exponents are nonnegative integers. Note that // if the combination 0^0 is encountered, it should be treated // as 1. // // Integral ( over unit tetrahedron ) x^l y^m z^n dx dy // = l! * m! * n! / ( l + m + n + 3 )! // // Licensing: // // This code is distributed under the GNU LGPL license. // // Modified: // // 04 July 2007 // // Author: // // John Burkardt // // Parameters: // // Input, int DIM_NUM, the spatial dimension. // // Input, int EXPON[DIM_NUM], the exponents. // // Output, double TET01_MONOMIAL_INTEGRAL, the value of the integral of the // monomial. // { int i; int k; double value; // // The first computation ends with VALUE = 1.0; // value = 1.0; k = 0; for ( i = 1; i <= expon[0]; i++ ) { k = k + 1; // value = value * ( double ) ( i ) / ( double ) ( k ); } for ( i = 1; i <= expon[1]; i++ ) { k = k + 1; value = value * ( double ) ( i ) / ( double ) ( k ); } for ( i = 1; i <= expon[2]; i++ ) { k = k + 1; value = value * ( double ) ( i ) / ( double ) ( k ); } k = k + 1; value = value / ( double ) ( k ); k = k + 1; value = value / ( double ) ( k ); k = k + 1; value = value / ( double ) ( k ); return value; } //****************************************************************************80 double tet01_monomial_quadrature ( int dim_num, int expon[], int point_num, double x[], double weight[] ) //****************************************************************************80 // // Purpose: // // TET01_MONOMIAL_QUADRATURE applies quadrature to a monomial in a tetrahedron. // // Licensing: // // This code is distributed under the GNU LGPL license. // // Modified: // // 05 July 2007 // // Author: // // John Burkardt // // Parameters: // // Input, int DIM_NUM, the spatial dimension. // // Input, int EXPON[DIM_NUM], the exponents. // // Input, int POINT_NUM, the number of points in the rule. // // Input, double X[DIM_NUM*POINT_NUM], the quadrature points. // // Input, double WEIGHT[POINT_NUM], the quadrature weights. // // Output, double MONOMIAL_QUADRATURE, the quadrature error. // { double exact; double quad; double quad_error; double scale; double *value; double volume; // // Get the exact value of the integral of the unscaled monomial. // scale = tet01_monomial_integral ( dim_num, expon ); // // Evaluate the monomial at the quadrature points. // value = monomial_value ( dim_num, point_num, x, expon ); // // Compute the weighted sum and divide by the exact value. // volume = 1.0 / 6.0; quad = volume * r8vec_dot ( point_num, weight, value ) / scale; // // Error: // exact = 1.0; quad_error = fabs ( quad - exact ); delete [] value; return quad_error; } //****************************************************************************80 void tetrahedron_order4_physical_to_reference ( double tetra[], int n, double phy[], double ref[] ) //****************************************************************************80 // // Purpose: // // TETRAHEDRON_ORDER4_PHYSICAL_TO_REFERENCE maps physical points to reference points. // // Discussion: // // Given the vertices of an order 4 physical tetrahedron and a point // (X,Y,Z) in the physical tetrahedron, the routine computes the value // of the corresponding image point (R,S,T) in reference space. // // This routine may be appropriate for an order 10 tetrahedron, // if the mapping between reference and physical space is linear. // This implies, in particular, that the edges of the image tetrahedron // are straight, the faces are flat, and the "midside" nodes in the // physical tetrahedron are halfway along the sides of the physical // tetrahedron. // // Licensing: // // This code is distributed under the GNU LGPL license. // // Modified: // // 06 December 2006 // // Author: // // John Burkardt // // Parameters: // // Input, double TETRA[3*4], the coordinates of the vertices. // The vertices are assumed to be the images of // (0,0,0), (1,0,0), (0,1,0) and (0,0,1) respectively. // // Input, int N, the number of points to transform. // // Input, double PHY[3*N], the coordinates of physical points // to be transformed. // // Output, double REF[3*N], the coordinates of the corresponding // points in the reference space. // { double a[3*3]; double det; int i; int j; // // Set up the matrix. // for ( i = 0; i < 3; i++ ) { a[i+0*3] = tetra[i+1*3] - tetra[i+0*3]; a[i+1*3] = tetra[i+2*3] - tetra[i+0*3]; a[i+2*3] = tetra[i+3*3] - tetra[i+0*3]; } // // Compute the determinant. // det = a[0+0*3] * ( a[1+1*3] * a[2+2*3] - a[1+2*3] * a[2+1*3] ) + a[0+1*3] * ( a[1+2*3] * a[2+0*3] - a[1+0*3] * a[2+2*3] ) + a[0+2*3] * ( a[1+0*3] * a[2+1*3] - a[1+1*3] * a[2+0*3] ); // // If the determinant is zero, bail out. // if ( det == 0.0 ) { for ( j = 0; j < n; j++ ) { for ( i = 0; i < 3; i++ ) { ref[i+j*3] = 0.0; } } return; } // // Compute the solution. // for ( j = 0; j < n; j++ ) { ref[0+j*3] = ( ( a[1+1*3] * a[2+2*3] - a[1+2*3] * a[2+1*3] ) * ( phy[0+j*3] - tetra[0+0*3] ) - ( a[0+1*3] * a[2+2*3] - a[0+2*3] * a[2+1*3] ) * ( phy[1+j*3] - tetra[1+0*3] ) + ( a[0+1*3] * a[1+2*3] - a[0+2*3] * a[1+1*3] ) * ( phy[2+j*3] - tetra[2+0*3] ) ) / det; ref[1+j*3] = ( - ( a[1+0*3] * a[2+2*3] - a[1+2*3] * a[2+0*3] ) * ( phy[0+j*3] - tetra[0+0*3] ) + ( a[0+0*3] * a[2+2*3] - a[0+2*3] * a[2+0*3] ) * ( phy[1+j*3] - tetra[1+0*3] ) - ( a[0+0*3] * a[1+2*3] - a[0+2*3] * a[1+0*3] ) * ( phy[2+j*3] - tetra[2+0*3] ) ) / det; ref[2+j*3] = ( ( a[1+0*3] * a[2+1*3] - a[1+1*3] * a[2+0*3] ) * ( phy[0+j*3] - tetra[0+0*3] ) - ( a[0+0*3] * a[2+1*3] - a[0+1*3] * a[2+0*3] ) * ( phy[1+j*3] - tetra[1+0*3] ) + ( a[0+0*3] * a[1+1*3] - a[0+1*3] * a[1+0*3] ) * ( phy[2+j*3] - tetra[2+0*3] ) ) / det; } return; } //****************************************************************************80 double tetrahedron_volume ( double tetra[3*4] ) //****************************************************************************80 // // Purpose: // // TETRAHEDRON_VOLUME computes the volume of a tetrahedron. // // Licensing: // // This code is distributed under the GNU LGPL license. // // Modified: // // 06 August 2005 // // Author: // // John Burkardt // // Parameters: // // Input, double TETRA[3*4], the vertices of the tetrahedron. // // Output, double TETRAHEDRON_VOLUME, the volume of the tetrahedron. // { double a[4*4]; int i; int j; double volume; for ( i = 0; i < 3; i++ ) { for ( j = 0; j < 4; j++ ) { a[i+j*4] = tetra[i+j*3]; } } i = 3; for ( j = 0; j < 4; j++ ) { a[i+j*4] = 1.0; } volume = fabs ( r8mat_det_4d ( a ) ) / 6.0; return volume; } //****************************************************************************80 void timestamp ( ) //****************************************************************************80 // // Purpose: // // TIMESTAMP prints the current YMDHMS date as a time stamp. // // Example: // // 31 May 2001 09:45:54 AM // // Licensing: // // This code is distributed under the GNU LGPL license. // // Modified: // // 08 July 2009 // // Author: // // John Burkardt // // Parameters: // // None // { # define TIME_SIZE 40 static char time_buffer[TIME_SIZE]; const struct std::tm *tm_ptr; size_t len; std::time_t now; now = std::time ( NULL ); tm_ptr = std::localtime ( &now ); len = std::strftime ( time_buffer, TIME_SIZE, "%d %B %Y %I:%M:%S %p", tm_ptr ); std::cout << time_buffer << "\n"; return; # undef TIME_SIZE }