program main !*****************************************************************************80 ! !! MAIN is the main program for TABLE_VORONOI. ! ! Discussion: ! ! TABLE_VORONOI uses GEOMPACK to get some Voronoi diagram information. ! ! Usage: ! ! table_voronoi file_name ! ! Licensing: ! ! This code is distributed under the GNU LGPL license. ! ! Modified: ! ! 21 October 2004 ! ! Author: ! ! John Burkardt ! ! Parameters: ! ! Command line argument FILE_NAME, is the name of a file in the ! TABLE format containing the coordinates of a point set to be analyzed. ! implicit none integer ( kind = 4 ) arg_num character ( len = 255 ) file_name integer ( kind = 4 ) iarg integer ( kind = 4 ) iargc call timestamp ( ) write ( *, '(a)' ) ' ' write ( *, '(a)' ) 'TABLE_VORONOI' write ( *, '(a)' ) ' FORTRAN90 version:' write ( *, '(a)' ) ' This program is given the coordinates of a set of' write ( *, '(a)' ) ' points in the plane, calls GEOMPACK to determine the' write ( *, '(a)' ) ' Delaunay triangulation of those points, and then' write ( *, '(a)' ) ' digests that data to produce information defining' write ( *, '(a)' ) ' the Voronoi diagram.' write ( *, '(a)' ) ' ' write ( *, '(a)' ) ' The input file contains the following data:' write ( *, '(a)' ) ' ' write ( *, '(a)' ) ' G_NUM: the number of generators;' write ( *, '(a)' ) ' G_XY: the (X,Y) coordinates of the generators.' write ( *, '(a)' ) ' ' write ( *, '(a)' ) ' The computed Voronoi information includes:' write ( *, '(a)' ) ' ' write ( *, '(a)' ) ' G_DEGREE: the degree of each Voronoi cell;' write ( *, '(a)' ) ' G_START: the index of the first Voronoi vertex;' write ( *, '(a)' ) ' G_FACE: the list of all Voronoi vertices;' write ( *, '(a)' ) ' ' write ( *, '(a)' ) ' V_NUM: the number of (finite) Voronoi vertices;' write ( *, '(a)' ) ' V_XY: the (X,Y) coordinates of the ' // & 'Voronoi vertices;' write ( *, '(a)' ) ' ' write ( *, '(a)' ) ' I_NUM: the number of Voronoi vertices at infinity;' write ( *, '(a)' ) ' I_XY: the directions associated with the Voronoi' write ( *, '(a)' ) ' vertices at infinity.' ! ! Get the number of command line arguments. ! arg_num = iargc ( ) ! ! If at least one command line argument, it's the input file name. ! if ( 1 <= arg_num ) then iarg = 1 call getarg ( iarg, file_name ) else write ( *, '(a)' ) ' ' write ( *, '(a)' ) ' Enter the generator file name.' read ( *, '(a)' ) file_name end if call handle_file ( file_name ) ! ! Terminate. ! write ( *, '(a)' ) ' ' write ( *, '(a)' ) 'TABLE_VORONOI' write ( *, '(a)' ) ' Normal end of execution.' write ( *, '(a)' ) ' ' call timestamp ( ) stop end subroutine ch_cap ( c ) !*****************************************************************************80 ! !! CH_CAP capitalizes a single character. ! ! Licensing: ! ! This code is distributed under the GNU LGPL license. ! ! Modified: ! ! 19 July 1998 ! ! Author: ! ! John Burkardt ! ! Parameters: ! ! Input/output, character C, the character to capitalize. ! implicit none character c integer ( kind = 4 ) itemp itemp = ichar ( c ) if ( 97 <= itemp .and. itemp <= 122 ) then c = char ( itemp - 32 ) end if return end function ch_eqi ( c1, c2 ) !*****************************************************************************80 ! !! CH_EQI is a case insensitive comparison of two characters for equality. ! ! Example: ! ! CH_EQI ( 'A', 'a' ) is .TRUE. ! ! Licensing: ! ! This code is distributed under the GNU LGPL license. ! ! Modified: ! ! 28 July 2000 ! ! Author: ! ! John Burkardt ! ! Parameters: ! ! Input, character C1, C2, the characters to compare. ! ! Output, logical CH_EQI, the result of the comparison. ! implicit none character c1 character c1_cap character c2 character c2_cap logical ch_eqi c1_cap = c1 c2_cap = c2 call ch_cap ( c1_cap ) call ch_cap ( c2_cap ) if ( c1_cap == c2_cap ) then ch_eqi = .true. else ch_eqi = .false. end if return end subroutine ch_to_digit ( c, digit ) !*****************************************************************************80 ! !! CH_TO_DIGIT returns the value of a base 10 digit. ! ! Example: ! ! C DIGIT ! --- ----- ! '0' 0 ! '1' 1 ! ... ... ! '9' 9 ! ' ' 0 ! 'X' -1 ! ! Licensing: ! ! This code is distributed under the GNU LGPL license. ! ! Modified: ! ! 04 August 1999 ! ! Author: ! ! John Burkardt ! ! Parameters: ! ! Input, character C, the decimal digit, '0' through '9' or blank ! are legal. ! ! Output, integer ( kind = 4 ) DIGIT, the corresponding value. ! If C was 'illegal', then DIGIT is -1. ! implicit none character c integer ( kind = 4 ) digit if ( lge ( c, '0' ) .and. lle ( c, '9' ) ) then digit = ichar ( c ) - 48 else if ( c == ' ' ) then digit = 0 else digit = -1 end if return end function diaedg ( x0, y0, x1, y1, x2, y2, x3, y3 ) !*****************************************************************************80 ! !! DIAEDG chooses a diagonal edge. ! ! Discussion: ! ! The routine determines whether 0--2 or 1--3 is the diagonal edge ! that should be chosen, based on the circumcircle criterion, where ! (X0,Y0), (X1,Y1), (X2,Y2), (X3,Y3) are the vertices of a simple ! quadrilateral in counterclockwise order. ! ! Licensing: ! ! This code is distributed under the GNU LGPL license. ! ! Modified: ! ! 19 February 2001 ! ! Author: ! ! Original FORTRAN77 version by Barry Joe. ! FORTRAN90 version by John Burkardt. ! ! Reference: ! ! Barry Joe, ! GEOMPACK - a software package for the generation of meshes ! using geometric algorithms, ! Advances in Engineering Software, ! Volume 13, pages 325-331, 1991. ! ! Parameters: ! ! Input, real ( kind = 8 ) X0, Y0, X1, Y1, X2, Y2, X3, Y3, the ! coordinates of the vertices of a quadrilateral, given in ! counter clockwise order. ! ! Output, integer ( kind = 4 ) DIAEDG, chooses a diagonal: ! +1, if diagonal edge 02 is chosen; ! -1, if diagonal edge 13 is chosen; ! 0, if the four vertices are cocircular. ! implicit none real ( kind = 8 ) ca real ( kind = 8 ) cb integer ( kind = 4 ) diaedg real ( kind = 8 ) dx10 real ( kind = 8 ) dx12 real ( kind = 8 ) dx30 real ( kind = 8 ) dx32 real ( kind = 8 ) dy10 real ( kind = 8 ) dy12 real ( kind = 8 ) dy30 real ( kind = 8 ) dy32 real ( kind = 8 ) s real ( kind = 8 ) tol real ( kind = 8 ) tola real ( kind = 8 ) tolb real ( kind = 8 ) x0 real ( kind = 8 ) x1 real ( kind = 8 ) x2 real ( kind = 8 ) x3 real ( kind = 8 ) y0 real ( kind = 8 ) y1 real ( kind = 8 ) y2 real ( kind = 8 ) y3 tol = 100.0D+00 * epsilon ( tol ) dx10 = x1 - x0 dy10 = y1 - y0 dx12 = x1 - x2 dy12 = y1 - y2 dx30 = x3 - x0 dy30 = y3 - y0 dx32 = x3 - x2 dy32 = y3 - y2 tola = tol * max ( abs ( dx10 ), abs ( dy10 ), abs ( dx30 ), abs ( dy30 ) ) tolb = tol * max ( abs ( dx12 ), abs ( dy12 ), abs ( dx32 ), abs ( dy32 ) ) ca = dx10 * dx30 + dy10 * dy30 cb = dx12 * dx32 + dy12 * dy32 if ( tola < ca .and. tolb < cb ) then diaedg = -1 else if ( ca < -tola .and. cb < -tolb ) then diaedg = 1 else tola = max ( tola, tolb ) s = ( dx10 * dy30 - dx30 * dy10 ) * cb + ( dx32 * dy12 - dx12 * dy32 ) * ca if ( tola < s ) then diaedg = -1 else if ( s < -tola ) then diaedg = 1 else diaedg = 0 end if end if return end subroutine dtris2 ( point_num, point_xy, tri_num, tri_vert, tri_nabe ) !*****************************************************************************80 ! !! DTRIS2 constructs a Delaunay triangulation of 2D vertices. ! ! Discussion: ! ! The routine constructs the Delaunay triangulation of a set of 2D vertices ! using an incremental approach and diagonal edge swaps. Vertices are ! first sorted in lexicographically increasing (X,Y) order, and ! then are inserted one at a time from outside the convex hull. ! ! Licensing: ! ! This code is distributed under the GNU LGPL license. ! ! Modified: ! ! 15 January 2004 ! ! Author: ! ! Original FORTRAN77 version by Barry Joe. ! FORTRAN90 version by John Burkardt. ! ! Reference: ! ! Barry Joe, ! GEOMPACK - a software package for the generation of meshes ! using geometric algorithms, ! Advances in Engineering Software, ! Volume 13, pages 325-331, 1991. ! ! Parameters: ! ! Input, integer ( kind = 4 ) POINT_NUM, the number of vertices. ! ! Input/output, real ( kind = 8 ) POINT_XY(2,POINT_NUM), the vertices. ! On output, the vertices have been sorted into dictionary order. ! ! Output, integer ( kind = 4 ) TRI_NUM, the number of triangles in the ! triangulation; TRI_NUM is equal to 2*POINT_NUM - NB - 2, where NB is the ! number of boundary vertices. ! ! Output, integer ( kind = 4 ) TRI_VERT(3,TRI_NUM), the nodes that make up ! each triangle. The elements are indices of POINT_XY. The vertices of the ! triangles are in counter clockwise order. ! ! Output, integer ( kind = 4 ) TRI_NABE(3,TRI_NUM), the triangle neighbor ! list. Positive elements are indices of TIL; negative elements are used ! for links of a counter clockwise linked list of boundary edges; ! LINK = -(3*I + J-1) where I, J = triangle, edge index; TRI_NABE(J,I) refers ! to the neighbor along edge from vertex J to J+1 (mod 3). ! implicit none integer ( kind = 4 ), parameter :: dim_num = 2 integer ( kind = 4 ) point_num real ( kind = 8 ) cmax integer ( kind = 4 ) e integer ( kind = 4 ) i integer ( kind = 4 ) ierr integer ( kind = 4 ) indx(point_num) integer ( kind = 4 ) j integer ( kind = 4 ) k integer ( kind = 4 ) l integer ( kind = 4 ) ledg integer ( kind = 4 ) lr integer ( kind = 4 ) lrline integer ( kind = 4 ) ltri integer ( kind = 4 ) m integer ( kind = 4 ) m1 integer ( kind = 4 ) m2 integer ( kind = 4 ) n real ( kind = 8 ) point_xy(dim_num,point_num) integer ( kind = 4 ) redg integer ( kind = 4 ) rtri integer ( kind = 4 ) stack(point_num) integer ( kind = 4 ) t real ( kind = 8 ) tol integer ( kind = 4 ) top integer ( kind = 4 ) tri_nabe(3,point_num*2) integer ( kind = 4 ) tri_num integer ( kind = 4 ) tri_vert(3,point_num*2) tol = 100.0D+00 * epsilon ( tol ) ierr = 0 ! ! Sort the vertices by increasing (x,y). ! call r82vec_sort_heap_index_a ( point_num, point_xy, indx ) call r82vec_permute ( point_num, point_xy, indx ) ! ! Make sure that the data points are "reasonably" distinct. ! m1 = 1 do i = 2, point_num m = m1 m1 = i k = 0 do j = 1, 2 cmax = max ( abs ( point_xy(j,m) ), abs ( point_xy(j,m1) ) ) if ( tol * ( cmax + 1.0D+00 ) & < abs ( point_xy(j,m) - point_xy(j,m1) ) ) then k = j exit end if end do if ( k == 0 ) then write ( *, '(a)' ) ' ' write ( *, '(a)' ) 'DTRIS2 - Fatal error!' write ( *, '(a,i8)' ) ' Fails for point number I = ', i write ( *, '(a,i8)' ) ' M = ', m write ( *, '(a,i8)' ) ' M1 = ', m1 write ( *, '(a,2g14.6)' ) ' X,Y(M) = ', point_xy(1,m), point_xy(2,m) write ( *, '(a,2g14.6)' ) ' X,Y(M1) = ', point_xy(1,m1), point_xy(2,m1) ierr = 224 return end if end do ! ! Starting from points M1 and M2, search for a third point M that ! makes a "healthy" triangle (M1,M2,M) ! m1 = 1 m2 = 2 j = 3 do if ( point_num < j ) then write ( *, '(a)' ) ' ' write ( *, '(a)' ) 'DTRIS2 - Fatal error!' ierr = 225 return end if m = j lr = lrline ( point_xy(1,m), point_xy(2,m), point_xy(1,m1), & point_xy(2,m1), point_xy(1,m2), point_xy(2,m2), 0.0D+00 ) if ( lr /= 0 ) then exit end if j = j + 1 end do ! ! Set up the triangle information for (M1,M2,M), and for any other ! triangles you created because points were collinear with M1, M2. ! tri_num = j - 2 if ( lr == -1 ) then tri_vert(1,1) = m1 tri_vert(2,1) = m2 tri_vert(3,1) = m tri_nabe(3,1) = -3 do i = 2, tri_num m1 = m2 m2 = i+1 tri_vert(1,i) = m1 tri_vert(2,i) = m2 tri_vert(3,i) = m tri_nabe(1,i-1) = -3 * i tri_nabe(2,i-1) = i tri_nabe(3,i) = i - 1 end do tri_nabe(1,tri_num) = -3 * tri_num - 1 tri_nabe(2,tri_num) = -5 ledg = 2 ltri = tri_num else tri_vert(1,1) = m2 tri_vert(2,1) = m1 tri_vert(3,1) = m tri_nabe(1,1) = -4 do i = 2, tri_num m1 = m2 m2 = i+1 tri_vert(1,i) = m2 tri_vert(2,i) = m1 tri_vert(3,i) = m tri_nabe(3,i-1) = i tri_nabe(1,i) = -3 * i - 3 tri_nabe(2,i) = i - 1 end do tri_nabe(3,tri_num) = -3 * tri_num tri_nabe(2,1) = -3 * tri_num - 2 ledg = 2 ltri = 1 end if ! ! Insert the vertices one at a time from outside the convex hull, ! determine visible boundary edges, and apply diagonal edge swaps until ! Delaunay triangulation of vertices (so far) is obtained. ! top = 0 do i = j+1, point_num m = i m1 = tri_vert(ledg,ltri) if ( ledg <= 2 ) then m2 = tri_vert(ledg+1,ltri) else m2 = tri_vert(1,ltri) end if lr = lrline ( point_xy(1,m), point_xy(2,m), point_xy(1,m1), & point_xy(2,m1), point_xy(1,m2), point_xy(2,m2), 0.0D+00 ) if ( 0 < lr ) then rtri = ltri redg = ledg ltri = 0 else l = -tri_nabe(ledg,ltri) rtri = l / 3 redg = mod(l,3) + 1 end if call vbedg ( point_xy(1,m), point_xy(2,m), point_num, point_xy, tri_num, & tri_vert, tri_nabe, ltri, ledg, rtri, redg ) n = tri_num + 1 l = -tri_nabe(ledg,ltri) do t = l / 3 e = mod ( l, 3 ) + 1 l = -tri_nabe(e,t) m2 = tri_vert(e,t) if ( e <= 2 ) then m1 = tri_vert(e+1,t) else m1 = tri_vert(1,t) end if tri_num = tri_num + 1 tri_nabe(e,t) = tri_num tri_vert(1,tri_num) = m1 tri_vert(2,tri_num) = m2 tri_vert(3,tri_num) = m tri_nabe(1,tri_num) = t tri_nabe(2,tri_num) = tri_num - 1 tri_nabe(3,tri_num) = tri_num + 1 top = top + 1 if ( point_num < top ) then ierr = 8 write ( *, '(a)' ) ' ' write ( *, '(a)' ) 'DTRIS2 - Fatal error!' write ( *, '(a)' ) ' Stack overflow.' return end if stack(top) = tri_num if ( t == rtri .and. e == redg ) then exit end if end do tri_nabe(ledg,ltri) = -3 * n - 1 tri_nabe(2,n) = -3 * tri_num - 2 tri_nabe(3,tri_num) = -l ltri = n ledg = 2 call swapec ( m, top, ltri, ledg, point_num, point_xy, tri_num, & tri_vert, tri_nabe, stack, ierr ) if ( ierr /= 0 ) then write ( *, '(a)' ) ' ' write ( *, '(a)' ) 'DTRIS2 - Fatal error!' write ( *, '(a)' ) ' Error return from SWAPEC.' return end if end do ! ! Now account for the sorting that we did. ! do i = 1, 3 do j = 1, tri_num tri_vert(i,j) = indx ( tri_vert(i,j) ) end do end do call perm_inv ( point_num, indx ) call r82vec_permute ( point_num, point_xy, indx ) return end subroutine file_column_count ( input_filename, column_num ) !*****************************************************************************80 ! !! FILE_COLUMN_COUNT counts the number of columns in the first line of a file. ! ! Discussion: ! ! The file is assumed to be a simple text file. ! ! Most lines of the file is presumed to consist of COLUMN_NUM words, ! separated by spaces. There may also be some blank lines, and some ! comment lines, ! which have a "#" in column 1. ! ! The routine tries to find the first non-comment non-blank line and ! counts the number of words in that line. ! ! If all lines are blanks or comments, it goes back and tries to analyze ! a comment line. ! ! Licensing: ! ! This code is distributed under the GNU LGPL license. ! ! Modified: ! ! 21 June 2001 ! ! Author: ! ! John Burkardt ! ! Parameters: ! ! Input, character ( len = * ) INPUT_FILENAME, the name of the file. ! ! Output, integer ( kind = 4 ) COLUMN_NUM, the number of columns in the file. ! implicit none integer ( kind = 4 ) column_num logical got_one character ( len = * ) input_filename integer ( kind = 4 ) input_unit integer ( kind = 4 ) ios character ( len = 255 ) line ! ! Open the file. ! call get_unit ( input_unit ) open ( unit = input_unit, file = input_filename, status = 'old', & form = 'formatted', access = 'sequential', iostat = ios ) if ( ios /= 0 ) then column_num = -1 write ( *, '(a)' ) ' ' write ( *, '(a)' ) 'FILE_COLUMN_COUNT - Fatal error!' write ( *, '(a)' ) ' Could not open the file:' write ( *, '(a)' ) ' ' // trim ( input_filename ) return end if ! ! Read one line, but skip blank lines and comment lines. ! got_one = .false. do read ( input_unit, '(a)', iostat = ios ) line if ( ios /= 0 ) then exit end if if ( len_trim ( line ) == 0 ) then cycle end if if ( line(1:1) == '#' ) then cycle end if got_one = .true. exit end do if ( .not. got_one ) then rewind ( input_unit ) do read ( input_unit, '(a)', iostat = ios ) line if ( ios /= 0 ) then exit end if if ( len_trim ( line ) == 0 ) then cycle end if got_one = .true. exit end do end if close ( unit = input_unit ) if ( .not. got_one ) then write ( *, '(a)' ) ' ' write ( *, '(a)' ) 'FILE_COLUMN_COUNT - Warning!' write ( *, '(a)' ) ' The file does not seem to contain any data.' column_num = -1 return end if call s_word_count ( line, column_num ) return end subroutine file_row_count ( input_filename, row_num ) !*****************************************************************************80 ! !! 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: ! ! 06 March 2003 ! ! Author: ! ! John Burkardt ! ! Parameters: ! ! Input, character ( len = * ) INPUT_FILENAME, the name of the input file. ! ! Output, integer ROW_NUM, the number of rows found. ! implicit none integer ( kind = 4 ) bad_num integer ( kind = 4 ) comment_num integer ( kind = 4 ) ierror character ( len = * ) input_filename integer ( kind = 4 ) input_unit integer ( kind = 4 ) ios character ( len = 255 ) line integer ( kind = 4 ) record_num integer ( kind = 4 ) row_num call get_unit ( input_unit ) open ( unit = input_unit, file = input_filename, status = 'old', & iostat = ios ) if ( ios /= 0 ) then row_num = -1; ierror = 1 write ( *, '(a)' ) ' ' write ( *, '(a)' ) 'FILE_ROW_COUNT - Fatal error!' write ( *, '(a)' ) ' Could not open the input file: ' // & trim ( input_filename ) stop end if comment_num = 0 row_num = 0 record_num = 0 bad_num = 0 do read ( input_unit, '(a)', iostat = ios ) line if ( ios /= 0 ) then ierror = record_num exit end if record_num = record_num + 1 if ( line(1:1) == '#' ) then comment_num = comment_num + 1 cycle end if if ( len_trim ( line ) == 0 ) then comment_num = comment_num + 1 cycle end if row_num = row_num + 1 end do close ( unit = input_unit ) return end subroutine get_unit ( iunit ) !*****************************************************************************80 ! !! GET_UNIT returns a free FORTRAN unit number. ! ! Discussion: ! ! A "free" FORTRAN unit number is an integer between 1 and 99 which ! is not currently associated with an I/O device. A free FORTRAN unit ! number is needed in order to open a file with the OPEN command. ! ! If IUNIT = 0, then no free FORTRAN unit could be found, although ! all 99 units were checked (except for units 5, 6 and 9, which ! are commonly reserved for console I/O). ! ! Otherwise, IUNIT is an integer between 1 and 99, representing a ! free FORTRAN unit. Note that GET_UNIT assumes that units 5 and 6 ! are special, and will never return those values. ! ! Licensing: ! ! This code is distributed under the GNU LGPL license. ! ! Modified: ! ! 18 September 2005 ! ! Author: ! ! John Burkardt ! ! Parameters: ! ! Output, integer ( kind = 4 ) IUNIT, the free unit number. ! implicit none integer ( kind = 4 ) i integer ( kind = 4 ) ios integer ( kind = 4 ) iunit logical lopen iunit = 0 do i = 1, 99 if ( i /= 5 .and. i /= 6 .and. i /= 9 ) then inquire ( unit = i, opened = lopen, iostat = ios ) if ( ios == 0 ) then if ( .not. lopen ) then iunit = i return end if end if end if end do return end subroutine handle_file ( file_name ) !*****************************************************************************80 ! !! HANDLE_FILE computes Voronoi information for a given set of data. ! ! Licensing: ! ! This code is distributed under the GNU LGPL license. ! ! Modified: ! ! 17 August 2003 ! ! Author: ! ! John Burkardt ! ! Parameters: ! ! Input, character FILE_NAME, the name of an XY file whose data ! is to be read and processed. ! ! Local Parameters: ! ! Local, integer ( kind = 4 ) G_DEGREE(G_NUM), the degree of each Voronoi ! cell. ! ! Local, integer ( kind = 4 ) G_FACE(6*G_NUM), the sequence of vertices to ! be used in a traversal of the boundary of the cell associated with each ! point. ! ! Local, integer ( kind = 4 ) G_NUM, the number of generators. ! ! Local, integer ( kind = 4 ) G_START(G_NUM), the index in G_FACE of the ! first vertex at which to begin a traversal of the boundary of the ! cell associated with each point. ! ! Local, real ( kind = 8 ) G_XY(2,G_NUM), the point coordinates. ! ! Local, integer ( kind = 4 ) I_NUM, the number of vertices at infinity of ! the Voronoi diagram. ! ! Local, real ( kind = 8 ) I_XY(2,I_NUM), the direction of the ! vertices at infinity. ! ! Local, integer ( kind = 4 ) V_NUM, the number of vertices of the Voronoi ! diagram. ! ! Local, real ( kind = 8 ) V_XY(2,V_NUM), the coordinates of the vertices ! of the Voronoi diagram. ! implicit none integer ( kind = 4 ), parameter :: dim_num = 2 character ( len = * ) :: file_name integer ( kind = 4 ) g_num integer ( kind = 4 ), allocatable, dimension ( : ) :: g_degree integer ( kind = 4 ), allocatable, dimension ( : ) :: g_face integer ( kind = 4 ), allocatable, dimension ( : ) :: g_start real ( kind = 8 ), allocatable, dimension ( :, : ) :: g_xy integer ( kind = 4 ) i integer ( kind = 4 ) i_num real ( kind = 8 ), allocatable, dimension ( :, : ) :: i_xy integer ( kind = 4 ) j integer ( kind = 4 ) k integer ( kind = 4 ) m integer ( kind = 4 ) v_num real ( kind = 8 ), allocatable, dimension ( :, : ) :: v_xy write ( *, '(a)' ) ' ' write ( *, '(a)' ) 'HANDLE_FILE' write ( *, '(a)' ) ' Read the TABLE file "' // trim ( file_name ) // '".' call r8mat_header_read ( file_name, m, g_num ) write ( *, '(a)' ) ' ' write ( *, '(a)' ) ' The header has been read.' write ( *, '(a)' ) ' ' write ( *, '(a,i8)' ) ' The spatial dimension of the data M = ', m write ( *, '(a,i8)' ) ' The number of generators, G_NUM = ', g_num if ( m < dim_num ) then write ( *, '(a)' ) ' ' write ( *, '(a)' ) 'HANDLE - Fatal error!' write ( *, '(a)' ) ' The input spatial dimension is less than 2.' stop else if ( dim_num < m ) then write ( *, '(a)' ) ' ' write ( *, '(a)' ) 'HANDLE - Warning:' write ( *, '(a)' ) ' The input spatial dimension is greater than 2.' write ( *, '(a)' ) ' Only the first two coordinates will be considered.' end if allocate ( i_xy(dim_num,g_num) ) allocate ( g_degree(g_num) ) allocate ( g_face(6*g_num) ) allocate ( g_start(g_num) ) allocate ( g_xy(dim_num,g_num) ) allocate ( v_xy(dim_num,2*g_num) ) call r8mat_data_read ( file_name, dim_num, g_num, g_xy ) write ( *, '(a)' ) ' ' write ( *, '(a)' ) ' The data has been read.' call r8mat_transpose_print ( dim_num, g_num, g_xy, ' The generators' ) call voronoi_data ( g_num, g_xy, g_degree, g_start, g_face, v_num, v_xy, & i_num, i_xy ) write ( *, '(a)' ) ' ' write ( *, '(a,i8)' ) ' G_START: The index of the first Voronoi vertex' write ( *, '(a,i8)' ) ' G_FACE: The Voronoi vertices' write ( *, '(a)' ) ' ' write ( *, '(a)' ) ' G G_START G_FACE' do j = 1, g_num k = g_start(j) write ( *, '(a)' ) ' ' write ( *, '(i8,5x,i8,2x,i8)' ) j, k, g_face(k) do i = 2, g_degree(j) k = k + 1 write ( *, '(8x,5x,8x,2x,i8)' ) g_face(k) end do end do write ( *, '(a)' ) ' ' write ( *, '(a,i8)' ) ' V_NUM: Number of Voronoi vertices = ', v_num call r8mat_transpose_print ( dim_num, v_num, v_xy, ' Voronoi vertices:' ) write ( *, '(a)' ) ' ' write ( *, '(a,i8)' ) & ' I_NUM: Number of Voronoi vertices at infinity = ', i_num call r8mat_transpose_print ( dim_num, i_num, i_xy, & ' Directions at infinity:' ) deallocate ( i_xy ) deallocate ( g_degree ) deallocate ( g_face ) deallocate ( g_start ) deallocate ( g_xy ) deallocate ( v_xy ) return end function i4_modp ( i, j ) !*****************************************************************************80 ! !! I4_MODP returns the nonnegative remainder of I4 division. ! ! Discussion: ! ! If ! NREM = I4_MODP ( I, J ) ! NMULT = ( I - NREM ) / J ! then ! I = J * NMULT + NREM ! where NREM is always nonnegative. ! ! The MOD function computes a result with the same sign as the ! quantity being divided. Thus, suppose you had an angle A, ! and you wanted to ensure that it was between 0 and 360. ! Then mod(A,360) would do, if A was positive, but if A ! was negative, your result would be between -360 and 0. ! ! On the other hand, I4_MODP(A,360) is between 0 and 360, always. ! ! Example: ! ! I J MOD I4_MODP Factorization ! ! 107 50 7 7 107 = 2 * 50 + 7 ! 107 -50 7 7 107 = -2 * -50 + 7 ! -107 50 -7 43 -107 = -3 * 50 + 43 ! -107 -50 -7 43 -107 = 3 * -50 + 43 ! ! Licensing: ! ! This code is distributed under the GNU LGPL license. ! ! Modified: ! ! 02 March 1999 ! ! Author: ! ! John Burkardt ! ! Parameters: ! ! Input, integer ( kind = 4 ) I, the number to be divided. ! ! Input, integer ( kind = 4 ) J, the number that divides I. ! ! Output, integer ( kind = 4 ) I4_MODP, the nonnegative remainder when I is ! divided by J. ! implicit none integer ( kind = 4 ) i integer ( kind = 4 ) i4_modp integer ( kind = 4 ) j if ( j == 0 ) then write ( *, '(a)' ) ' ' write ( *, '(a)' ) 'I4_MODP - Fatal error!' write ( *, '(a,i8)' ) ' I4_MODP ( I, J ) called with J = ', j stop end if i4_modp = mod ( i, j ) if ( i4_modp < 0 ) then i4_modp = i4_modp + abs ( j ) end if return end function i4_wrap ( ival, ilo, ihi ) !*****************************************************************************80 ! !! I4_WRAP forces an I4 to lie between given limits by wrapping. ! ! Example: ! ! ILO = 4, IHI = 8 ! ! I I4_WRAP ! ! -2 8 ! -1 4 ! 0 5 ! 1 6 ! 2 7 ! 3 8 ! 4 4 ! 5 5 ! 6 6 ! 7 7 ! 8 8 ! 9 4 ! 10 5 ! 11 6 ! 12 7 ! 13 8 ! 14 4 ! ! Licensing: ! ! This code is distributed under the GNU LGPL license. ! ! Modified: ! ! 15 July 2000 ! ! Author: ! ! John Burkardt ! ! Parameters: ! ! Input, integer ( kind = 4 ) IVAL, an integer value. ! ! Input, integer ( kind = 4 ) ILO, IHI, the desired bounds. ! ! Output, integer ( kind = 4 ) I4_WRAP, a "wrapped" version of IVAL. ! implicit none integer ( kind = 4 ) i4_modp integer ( kind = 4 ) i4_wrap integer ( kind = 4 ) ihi integer ( kind = 4 ) ilo integer ( kind = 4 ) ival integer ( kind = 4 ) wide wide = ihi + 1 - ilo if ( wide == 0 ) then i4_wrap = ilo else i4_wrap = ilo + i4_modp ( ival - ilo, wide ) end if return end subroutine i4mat_transpose_print ( m, n, a, title ) !*****************************************************************************80 ! !! I4MAT_TRANSPOSE_PRINT prints an I4MAT, transposed. ! ! Licensing: ! ! This code is distributed under the GNU LGPL license. ! ! Modified: ! ! 28 December 2004 ! ! Author: ! ! John Burkardt ! ! Parameters: ! ! Input, integer ( kind = 4 ) M, N, the number of rows and columns. ! ! Input, integer ( kind = 4 ) A(M,N), an M by N matrix to be printed. ! ! Input, character ( len = * ) TITLE, an optional title. ! implicit none integer ( kind = 4 ) m integer ( kind = 4 ) n integer ( kind = 4 ) a(m,n) character ( len = * ) title call i4mat_transpose_print_some ( m, n, a, 1, 1, m, n, title ) return end subroutine i4mat_transpose_print_some ( m, n, a, ilo, jlo, ihi, jhi, title ) !*****************************************************************************80 ! !! I4MAT_TRANSPOSE_PRINT_SOME prints some of the transpose of an I4MAT. ! ! Licensing: ! ! This code is distributed under the GNU LGPL license. ! ! Modified: ! ! 28 December 2004 ! ! Author: ! ! John Burkardt ! ! Parameters: ! ! Input, integer ( kind = 4 ) M, N, the number of rows and columns. ! ! Input, integer ( kind = 4 ) A(M,N), an M by N matrix to be printed. ! ! Input, integer ( kind = 4 ) ILO, JLO, the first row and column to print. ! ! Input, integer ( kind = 4 ) IHI, JHI, the last row and column to print. ! ! Input, character ( kind = 4 ) ( len = * ) TITLE, an optional title. ! implicit none integer ( kind = 4 ), parameter :: incx = 10 integer ( kind = 4 ) m integer ( kind = 4 ) n integer ( kind = 4 ) a(m,n) character ( len = 7 ) ctemp(incx) integer ( kind = 4 ) i integer ( kind = 4 ) i2 integer ( kind = 4 ) i2hi integer ( kind = 4 ) i2lo integer ( kind = 4 ) ihi integer ( kind = 4 ) ilo integer ( kind = 4 ) inc integer ( kind = 4 ) j integer ( kind = 4 ) j2hi integer ( kind = 4 ) j2lo integer ( kind = 4 ) jhi integer ( kind = 4 ) jlo character ( len = * ) title write ( *, '(a)' ) ' ' write ( *, '(a)' ) trim ( title ) do i2lo = max ( ilo, 1 ), min ( ihi, m ), incx i2hi = i2lo + incx - 1 i2hi = min ( i2hi, m ) i2hi = min ( i2hi, ihi ) inc = i2hi + 1 - i2lo write ( *, '(a)' ) ' ' do i = i2lo, i2hi i2 = i + 1 - i2lo write ( ctemp(i2), '(i7)') i end do write ( *, '('' Row '',10a7)' ) ctemp(1:inc) write ( *, '(a)' ) ' Col' write ( *, '(a)' ) ' ' j2lo = max ( jlo, 1 ) j2hi = min ( jhi, n ) do j = j2lo, j2hi do i2 = 1, inc i = i2lo - 1 + i2 write ( ctemp(i2), '(i7)' ) a(i,j) end do write ( *, '(i5,1x,10a7)' ) j, ( ctemp(i), i = 1, inc ) end do end do return end subroutine i4vec_indicator ( n, a ) !*****************************************************************************80 ! !! I4VEC_INDICATOR sets an I4VEC to the indicator vector. ! ! Licensing: ! ! This code is distributed under the GNU LGPL license. ! ! Modified: ! ! 09 November 2000 ! ! Author: ! ! John Burkardt ! ! Parameters: ! ! Input, integer ( kind = 4 ) N, the number of elements of A. ! ! Output, integer ( kind = 4 ) A(N), the array to be initialized. ! implicit none integer ( kind = 4 ) n integer ( kind = 4 ) a(n) integer ( kind = 4 ) i do i = 1, n a(i) = i end do return end subroutine i4vec_print ( n, a, title ) !*****************************************************************************80 ! !! I4VEC_PRINT prints an I4VEC. ! ! Licensing: ! ! This code is distributed under the GNU LGPL license. ! ! Modified: ! ! 28 November 2000 ! ! Author: ! ! John Burkardt ! ! Parameters: ! ! Input, integer ( kind = 4 ) N, the number of components of the vector. ! ! Input, integer ( kind = 4 ) A(N), the vector to be printed. ! ! Input, character ( len = * ) TITLE, a title to be printed first. ! TITLE may be blank. ! implicit none integer ( kind = 4 ) n integer ( kind = 4 ) a(n) integer ( kind = 4 ) big integer ( kind = 4 ) i character ( len = * ) title write ( *, '(a)' ) ' ' write ( *, '(a)' ) trim ( title ) big = maxval ( abs ( a(1:n) ) ) write ( *, '(a)' ) ' ' if ( big < 1000 ) then do i = 1, n write ( *, '(2x,i8,2x,i4)' ) i, a(i) end do else if ( big < 1000000 ) then do i = 1, n write ( *, '(2x,i8,2x,i7)' ) i, a(i) end do else do i = 1, n write ( *, '(2x,i8,2x,i12)' ) i, a(i) end do end if return end subroutine line_exp_normal_2d ( p1, p2, normal ) !*****************************************************************************80 ! !! LINE_EXP_NORMAL_2D computes the unit normal vector to a line in 2D. ! ! Discussion: ! ! The explicit form of a line in 2D is: ! ! ( P1, P2 ) = ( (X1,Y1), (X2,Y2) ). ! ! Licensing: ! ! This code is distributed under the GNU LGPL license. ! ! Modified: ! ! 01 January 2005 ! ! Author: ! ! John Burkardt ! ! Parameters: ! ! Input, real ( kind = 8 ) P1(2), P2(2), two points on the line. ! ! Output, real ( kind = 8 ) NORMAL(2), a unit normal vector to the line. ! implicit none integer ( kind = 4 ), parameter :: dim_num = 2 real ( kind = 8 ) norm real ( kind = 8 ) normal(dim_num) real ( kind = 8 ) p1(dim_num) real ( kind = 8 ) p2(dim_num) norm = sqrt ( ( p2(1) - p1(1) )**2 + ( p2(2) - p1(2) )**2 ) if ( norm == 0.0D+00 ) then normal(1:dim_num) = 0.0D+00 return end if normal(1) = ( p2(2) - p1(2) ) / norm normal(2) = - ( p2(1) - p1(1) ) / norm return end function lrline ( xu, yu, xv1, yv1, xv2, yv2, dv ) !*****************************************************************************80 ! !! LRLINE determines if a point is left of, right or, or on a directed line. ! ! Discussion: ! ! The directed line is parallel to, and at a signed distance DV from ! a directed base line from (XV1,YV1) to (XV2,YV2). ! ! Licensing: ! ! This code is distributed under the GNU LGPL license. ! ! Modified: ! ! 14 July 2001 ! ! Author: ! ! Original FORTRAN77 version by Barry Joe. ! FORTRAN90 version by John Burkardt. ! ! Reference: ! ! Barry Joe, ! GEOMPACK - a software package for the generation of meshes ! using geometric algorithms, ! Advances in Engineering Software, ! Volume 13, pages 325-331, 1991. ! ! Parameters: ! ! Input, real ( kind = 8 ) XU, YU, the coordinates of the point whose ! position relative to the directed line is to be determined. ! ! Input, real ( kind = 8 ) XV1, YV1, XV2, YV2, the coordinates of two points ! that determine the directed base line. ! ! Input, real ( kind = 8 ) DV, the signed distance of the directed line ! from the directed base line through the points (XV1,YV1) and (XV2,YV2). ! DV is positive for a line to the left of the base line. ! ! Output, integer ( kind = 4 ) LRLINE, the result: ! +1, the point is to the right of the directed line; ! 0, the point is on the directed line; ! -1, the point is to the left of the directed line. ! implicit none real ( kind = 8 ) dv real ( kind = 8 ) dx real ( kind = 8 ) dxu real ( kind = 8 ) dy real ( kind = 8 ) dyu integer ( kind = 4 ) lrline real ( kind = 8 ) t real ( kind = 8 ) tol real ( kind = 8 ) tolabs real ( kind = 8 ) xu real ( kind = 8 ) xv1 real ( kind = 8 ) xv2 real ( kind = 8 ) yu real ( kind = 8 ) yv1 real ( kind = 8 ) yv2 tol = 100.0D+00 * epsilon ( tol ) dx = xv2 - xv1 dy = yv2 - yv1 dxu = xu - xv1 dyu = yu - yv1 tolabs = tol * max ( abs ( dx ), abs ( dy ), abs ( dxu ), & abs ( dyu ), abs ( dv ) ) t = dy * dxu - dx * dyu + dv * sqrt ( dx * dx + dy * dy ) if ( tolabs < t ) then lrline = 1 else if ( -tolabs <= t ) then lrline = 0 else lrline = -1 end if return end subroutine perm_inv ( n, p ) !*****************************************************************************80 ! !! PERM_INV inverts a permutation "in place". ! ! Licensing: ! ! This code is distributed under the GNU LGPL license. ! ! Modified: ! ! 25 July 2000 ! ! Author: ! ! John Burkardt ! ! Parameters: ! ! Input, integer ( kind = 4 ) N, the number of objects being permuted. ! ! Input/output, integer ( kind = 4 ) P(N), the permutation, in standard ! index form. On output, P describes the inverse permutation ! implicit none integer ( kind = 4 ) n integer ( kind = 4 ) i integer ( kind = 4 ) i0 integer ( kind = 4 ) i1 integer ( kind = 4 ) i2 integer ( kind = 4 ) is integer ( kind = 4 ) p(n) if ( n <= 0 ) then write ( *, '(a)' ) ' ' write ( *, '(a)' ) 'PERM_INV - Fatal error!' write ( *, '(a,i8)' ) ' Input value of N = ', n stop end if is = 1 do i = 1, n i1 = p(i) do while ( i < i1 ) i2 = p(i1) p(i1) = -i2 i1 = i2 end do is = -sign ( 1, p(i) ) p(i) = sign ( p(i), is ) end do do i = 1, n i1 = -p(i) if ( 0 <= i1 ) then i0 = i do i2 = p(i1) p(i1) = i0 if ( i2 < 0 ) then exit end if i0 = i1 i1 = i2 end do end if end do return end subroutine r82vec_permute ( n, a, p ) !*****************************************************************************80 ! !! R82VEC_PERMUTE permutes an R82VEC in place. ! ! Discussion: ! ! This routine permutes an array of real "objects", but the same ! logic can be used to permute an array of objects of any arithmetic ! type, or an array of objects of any complexity. The only temporary ! storage required is enough to store a single object. The number ! of data movements made is N + the number of cycles of order 2 or more, ! which is never more than N + N/2. ! ! Example: ! ! Input: ! ! N = 5 ! P = ( 2, 4, 5, 1, 3 ) ! A = ( 1.0, 2.0, 3.0, 4.0, 5.0 ) ! (11.0, 22.0, 33.0, 44.0, 55.0 ) ! ! Output: ! ! A = ( 2.0, 4.0, 5.0, 1.0, 3.0 ) ! ( 22.0, 44.0, 55.0, 11.0, 33.0 ). ! ! Licensing: ! ! This code is distributed under the GNU LGPL license. ! ! Modified: ! ! 11 January 2004 ! ! Author: ! ! John Burkardt ! ! Parameters: ! ! Input, integer ( kind = 4 ) N, the number of objects. ! ! Input/output, real ( kind = 8 ) A(2,N), the array to be permuted. ! ! Input, integer ( kind = 4 ) P(N), the permutation. P(I) = J means ! that the I-th element of the output array should be the J-th ! element of the input array. P must be a legal permutation ! of the integers from 1 to N, otherwise the algorithm will ! fail catastrophically. ! implicit none integer ( kind = 4 ) n real ( kind = 8 ) a(2,n) real ( kind = 8 ) a_temp(2) integer ( kind = 4 ) iget integer ( kind = 4 ) iput integer ( kind = 4 ) istart integer ( kind = 4 ) p(n) ! ! Search for the next element of the permutation that has not been used. ! do istart = 1, n if ( p(istart) < 0 ) then cycle else if ( p(istart) == istart ) then p(istart) = - p(istart) cycle else a_temp(1:2) = a(1:2,istart) iget = istart ! ! Copy the new value into the vacated entry. ! do iput = iget iget = p(iget) p(iput) = - p(iput) if ( iget < 1 .or. n < iget ) then write ( *, '(a)' ) ' ' write ( *, '(a)' ) 'R82VEC_PERMUTE - Fatal error!' stop end if if ( iget == istart ) then a(1:2,iput) = a_temp(1:2) exit end if a(1:2,iput) = a(1:2,iget) end do end if end do ! ! Restore the signs of the entries. ! p(1:n) = -p(1:n) return end subroutine r82vec_sort_heap_index_a ( n, a, indx ) !*****************************************************************************80 ! !! R82VEC_SORT_HEAP_INDEX_A does an indexed heap ascending sort of an R82VEC. ! ! Discussion: ! ! The sorting is not actually carried out. Rather an index array is ! created which defines the sorting. This array may be used to sort ! or index the array, or to sort or index related arrays keyed on the ! original array. ! ! Once the index array is computed, the sorting can be carried out ! "implicitly: ! ! A(1:2,INDX(I)), I = 1 to N is sorted, ! ! or explicitly, by the call ! ! call R82VEC_PERMUTE ( N, A, INDX ) ! ! after which A(1:2,I), I = 1 to N is sorted. ! ! Licensing: ! ! This code is distributed under the GNU LGPL license. ! ! Modified: ! ! 11 January 2004 ! ! Author: ! ! John Burkardt ! ! Parameters: ! ! Input, integer ( kind = 4 ) N, the number of entries in the array. ! ! Input, real ( kind = 8 ) A(2,N), an array to be index-sorted. ! ! Output, integer ( kind = 4 ) INDX(N), the sort index. The ! I-th element of the sorted array is A(1:2,INDX(I)). ! implicit none integer ( kind = 4 ) n real ( kind = 8 ) a(2,n) real ( kind = 8 ) aval(2) integer ( kind = 4 ) i integer ( kind = 4 ) indx(n) integer ( kind = 4 ) indxt integer ( kind = 4 ) ir integer ( kind = 4 ) j integer ( kind = 4 ) l if ( n < 1 ) then return end if if ( n == 1 ) then indx(1) = 1 return end if call i4vec_indicator ( n, indx ) l = n / 2 + 1 ir = n do if ( 1 < l ) then l = l - 1 indxt = indx(l) aval(1:2) = a(1:2,indxt) else indxt = indx(ir) aval(1:2) = a(1:2,indxt) indx(ir) = indx(1) ir = ir - 1 if ( ir == 1 ) then indx(1) = indxt exit end if end if i = l j = l + l do while ( j <= ir ) if ( j < ir ) then if ( a(1,indx(j)) < a(1,indx(j+1)) .or. & ( a(1,indx(j)) == a(1,indx(j+1)) .and. & a(2,indx(j)) < a(2,indx(j+1)) ) ) then j = j + 1 end if end if if ( aval(1) < a(1,indx(j)) .or. & ( aval(1) == a(1,indx(j)) .and. & aval(2) < a(2,indx(j)) ) ) then indx(i) = indx(j) i = j j = j + j else j = ir + 1 end if end do indx(i) = indxt end do return end subroutine r8mat_data_read ( input_filename, m, n, table ) !*****************************************************************************80 ! !! R8MAT_DATA_READ reads data from an R8MAT file. ! ! Discussion: ! ! An R8MAT is an array of R8 values. ! ! Discussion: ! ! The file may contain more than N points, but this routine will ! return after reading N of them. ! ! Licensing: ! ! This code is distributed under the GNU LGPL license. ! ! Modified: ! ! 18 October 2008 ! ! Author: ! ! John Burkardt ! ! Parameters: ! ! Input, character ( len = * ) INPUT_FILENAME, the name of the input file. ! ! Input, integer ( kind = 4 ) M, the spatial dimension. ! ! Input, integer ( kind = 4 ) N, the number of points. ! ! Output, real ( kind = 8 ) TABLE(M,N), the table data. ! implicit none integer ( kind = 4 ) m integer ( kind = 4 ) n integer ( kind = 4 ) ierror character ( len = * ) input_filename integer ( kind = 4 ) input_status integer ( kind = 4 ) input_unit integer ( kind = 4 ) j character ( len = 255 ) line real ( kind = 8 ) table(m,n) real ( kind = 8 ) x(m) ierror = 0 call get_unit ( input_unit ) open ( unit = input_unit, file = input_filename, status = 'old', & iostat = input_status ) if ( input_status /= 0 ) then ierror = 1 write ( *, '(a)' ) ' ' write ( *, '(a)' ) 'R8MAT_DATA_READ - Fatal error!' write ( *, '(a,i8)' ) ' Could not open the input file "' // & trim ( input_filename ) // '" on unit ', input_unit stop end if j = 0 do while ( j < n ) read ( input_unit, '(a)', iostat = input_status ) line if ( input_status /= 0 ) then write ( *, '(a)' ) ' ' write ( *, '(a)' ) 'R8MAT_DATA_READ - Fatal error!' write ( *, '(a)' ) ' Error while reading lines of data.' write ( *, '(a,i8)' ) ' Number of values expected per line M = ', m write ( *, '(a,i8)' ) ' Number of data lines read, J = ', j write ( *, '(a,i8)' ) ' Number of data lines needed, N = ', n stop end if if ( line(1:1) == '#' .or. len_trim ( line ) == 0 ) then cycle end if call s_to_r8vec ( line, m, x, ierror ) if ( ierror /= 0 ) then cycle end if j = j + 1 table(1:m,j) = x(1:m) end do close ( unit = input_unit ) return end subroutine r8mat_header_read ( input_filename, m, n ) !*****************************************************************************80 ! !! R8MAT_HEADER_READ reads the header from an R8MAT file. ! ! Discussion: ! ! An R8MAT is an array of R8 values. ! ! Licensing: ! ! This code is distributed under the GNU LGPL license. ! ! Modified: ! ! 07 September 2004 ! ! Author: ! ! John Burkardt ! ! Parameters: ! ! Input, character ( len = * ) INPUT_FILENAME, the name of the input file. ! ! Output, integer ( kind = 4 ) M, spatial dimension. ! ! Output, integer ( kind = 4 ) N, the number of points. ! implicit none character ( len = * ) input_filename integer ( kind = 4 ) m integer ( kind = 4 ) n call file_column_count ( input_filename, m ) if ( m <= 0 ) then write ( *, '(a)' ) ' ' write ( *, '(a)' ) 'R8MAT_HEADER_READ - Fatal error!' write ( *, '(a)' ) ' There was some kind of I/O problem while trying' write ( *, '(a)' ) ' to count the number of data columns in' write ( *, '(a)' ) ' the file "' // trim ( input_filename ) // '".' stop end if call file_row_count ( input_filename, n ) if ( n <= 0 ) then write ( *, '(a)' ) ' ' write ( *, '(a)' ) 'R8MAT_HEADER_READ - Fatal error!' write ( *, '(a)' ) ' There was some kind of I/O problem while trying' write ( *, '(a)' ) ' to count the number of data rows in' write ( *, '(a)' ) ' the file "' // trim ( input_filename ) // '".' stop end if return end subroutine r8mat_transpose_print ( m, n, a, title ) !*****************************************************************************80 ! !! R8MAT_TRANSPOSE_PRINT prints an R8MAT, transposed. ! ! Licensing: ! ! This code is distributed under the GNU LGPL license. ! ! Modified: ! ! 14 June 2004 ! ! Author: ! ! John Burkardt ! ! Parameters: ! ! Input, integer ( kind = 4 ) M, N, the number of rows and columns. ! ! Input, real ( kind = 8 ) A(M,N), an M by N matrix to be printed. ! ! Input, character ( len = * ) TITLE, an optional title. ! implicit none integer ( kind = 4 ) m integer ( kind = 4 ) n real ( kind = 8 ) a(m,n) character ( len = * ) title call r8mat_transpose_print_some ( m, n, a, 1, 1, m, n, title ) return end subroutine r8mat_transpose_print_some ( m, n, a, ilo, jlo, ihi, jhi, title ) !*****************************************************************************80 ! !! R8MAT_TRANSPOSE_PRINT_SOME prints some of an R8MAT, transposed. ! ! Licensing: ! ! This code is distributed under the GNU LGPL license. ! ! Modified: ! ! 14 June 2004 ! ! Author: ! ! John Burkardt ! ! Parameters: ! ! Input, integer ( kind = 4 ) M, N, the number of rows and columns. ! ! Input, real ( kind = 8 ) A(M,N), an M by N matrix to be printed. ! ! Input, integer ( kind = 4 ) ILO, JLO, the first row and column to print. ! ! Input, integer ( kind = 4 ) IHI, JHI, the last row and column to print. ! ! Input, character ( len = * ) TITLE, an optional title. ! implicit none integer ( kind = 4 ), parameter :: incx = 5 integer ( kind = 4 ) m integer ( kind = 4 ) n real ( kind = 8 ) a(m,n) character ( len = 14 ) ctemp(incx) integer ( kind = 4 ) i integer ( kind = 4 ) i2 integer ( kind = 4 ) i2hi integer ( kind = 4 ) i2lo integer ( kind = 4 ) ihi integer ( kind = 4 ) ilo integer ( kind = 4 ) inc integer ( kind = 4 ) j integer ( kind = 4 ) j2hi integer ( kind = 4 ) j2lo integer ( kind = 4 ) jhi integer ( kind = 4 ) jlo character ( len = * ) title if ( 0 < len_trim ( title ) ) then write ( *, '(a)' ) ' ' write ( *, '(a)' ) trim ( title ) end if do i2lo = max ( ilo, 1 ), min ( ihi, m ), incx i2hi = i2lo + incx - 1 i2hi = min ( i2hi, m ) i2hi = min ( i2hi, ihi ) inc = i2hi + 1 - i2lo write ( *, '(a)' ) ' ' do i = i2lo, i2hi i2 = i + 1 - i2lo write ( ctemp(i2), '(i7,7x)') i end do write ( *, '('' Row '',5a14)' ) ctemp(1:inc) write ( *, '(a)' ) ' Col' write ( *, '(a)' ) ' ' j2lo = max ( jlo, 1 ) j2hi = min ( jhi, n ) do j = j2lo, j2hi do i2 = 1, inc i = i2lo - 1 + i2 write ( ctemp(i2), '(g14.6)' ) a(i,j) end do write ( *, '(i5,1x,5a14)' ) j, ( ctemp(i), i = 1, inc ) end do end do return end subroutine s_to_r8 ( s, r, ierror, length ) !*****************************************************************************80 ! !! 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 ! '-1D2ABCD' -100.0 ! '-1X2ABCD' -1.0 ! ' 2D-1' 0.2 ! '23.45' 23.45 ! '-4.2D+2' -420.0 ! '17d2' 1700.0 ! '-14d-2' -0.14 ! 'd2' 100.0 ! '-12.73d-9.23' -12.73 * 10.0^(-9.23) ! ! Licensing: ! ! This code is distributed under the GNU LGPL license. ! ! Modified: ! ! 12 February 2001 ! ! Author: ! ! John Burkardt ! ! Parameters: ! ! Input, character ( len = * ) 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, real ( kind = 8 ) R, the real value that was read from the string. ! ! Output, integer ( kind = 4 ) IERROR, error flag. ! ! 0, no errors occurred. ! ! 1, 2, 6 or 7, the input number was garbled. The ! value of IERROR is the last type of input successfully ! read. For instance, 1 means initial blanks, 2 means ! a plus or minus sign, and so on. ! ! Output, integer LENGTH, the number of characters read from ! the string to form the number, including any terminating ! characters such as a trailing comma or blanks. ! implicit none character c logical ch_eqi integer ( kind = 4 ) ierror integer ( kind = 4 ) ihave integer ( kind = 4 ) isgn integer ( kind = 4 ) iterm integer ( kind = 4 ) jbot integer ( kind = 4 ) jsgn integer ( kind = 4 ) jtop integer ( kind = 4 ) length integer ( kind = 4 ) nchar integer ( kind = 4 ) ndig real ( kind = 8 ) r real ( kind = 8 ) rbot real ( kind = 8 ) rexp real ( kind = 8 ) rtop character ( len = * ) s character, parameter :: TAB = char ( 9 ) nchar = len_trim ( s ) ierror = 0 r = 0.0D+00 length = - 1 isgn = 1 rtop = 0.0D+00 rbot = 1.0D+00 jsgn = 1 jtop = 0 jbot = 1 ihave = 1 iterm = 0 do length = length + 1 c = s(length+1:length+1) ! ! Blank or TAB character. ! if ( c == ' ' .or. c == TAB ) then if ( ihave == 2 ) then else if ( ihave == 6 .or. ihave == 7 ) then iterm = 1 else if ( 1 < ihave ) then ihave = 11 end if ! ! Comma. ! else if ( c == ',' .or. c == ';' ) then if ( ihave /= 1 ) then iterm = 1 ihave = 12 length = length + 1 end if ! ! Minus sign. ! else if ( c == '-' ) then if ( ihave == 1 ) then ihave = 2 isgn = - 1 else if ( ihave == 6 ) then ihave = 7 jsgn = - 1 else iterm = 1 end if ! ! Plus sign. ! else if ( c == '+' ) then if ( ihave == 1 ) then ihave = 2 else if ( ihave == 6 ) then ihave = 7 else iterm = 1 end if ! ! Decimal point. ! else if ( c == '.' ) then if ( ihave < 4 ) then ihave = 4 else if ( 6 <= ihave .and. ihave <= 8 ) then ihave = 9 else iterm = 1 end if ! ! Exponent marker. ! else if ( ch_eqi ( c, 'E' ) .or. ch_eqi ( c, 'D' ) ) then if ( ihave < 6 ) then ihave = 6 else iterm = 1 end if ! ! Digit. ! else if ( ihave < 11 .and. lge ( c, '0' ) .and. lle ( c, '9' ) ) then if ( ihave <= 2 ) then ihave = 3 else if ( ihave == 4 ) then ihave = 5 else if ( ihave == 6 .or. ihave == 7 ) then ihave = 8 else if ( ihave == 9 ) then ihave = 10 end if call ch_to_digit ( c, ndig ) if ( ihave == 3 ) then rtop = 10.0D+00 * rtop + real ( ndig, kind = 8 ) else if ( ihave == 5 ) then rtop = 10.0D+00 * rtop + real ( ndig, kind = 8 ) rbot = 10.0D+00 * rbot else if ( ihave == 8 ) then jtop = 10 * jtop + ndig else if ( ihave == 10 ) then jtop = 10 * jtop + ndig jbot = 10 * jbot end if ! ! Anything else is regarded as a terminator. ! else iterm = 1 end if ! ! If we haven't seen a terminator, and we haven't examined the ! entire string, go get the next character. ! if ( iterm == 1 .or. nchar <= length+1 ) then exit end if end do ! ! If we haven't seen a terminator, and we have examined the ! entire string, then we're done, and LENGTH is equal to NCHAR. ! if ( iterm /= 1 .and. length+1 == nchar ) then length = nchar end if ! ! 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 .or. ihave == 2 .or. ihave == 6 .or. ihave == 7 ) then ierror = ihave return end if ! ! Number seems OK. Form it. ! if ( jtop == 0 ) then rexp = 1.0D+00 else if ( jbot == 1 ) then rexp = 10.0D+00**( jsgn * jtop ) else rexp = jsgn * jtop rexp = rexp / jbot rexp = 10.0D+00**rexp end if end if r = isgn * rexp * rtop / rbot return end subroutine s_to_r8vec ( s, n, rvec, ierror ) !*****************************************************************************80 ! !! S_TO_R8VEC reads an R8VEC from a string. ! ! Licensing: ! ! This code is distributed under the GNU LGPL license. ! ! Modified: ! ! 19 February 2001 ! ! Author: ! ! John Burkardt ! ! Parameters: ! ! Input, character ( len = * ) S, the string to be read. ! ! Input, integer ( kind = 4 ) N, the number of values expected. ! ! Output, real ( kind = 8 ) RVEC(N), the values read from the string. ! ! Output, integer ( kind = 4 ) IERROR, error flag. ! 0, no errors occurred. ! -K, could not read data for entries -K through N. ! implicit none integer ( kind = 4 ) n integer ( kind = 4 ) i integer ( kind = 4 ) ierror integer ( kind = 4 ) ilo integer ( kind = 4 ) length real ( kind = 8 ) rvec(n) character ( len = * ) s i = 0 ierror = 0 ilo = 1 do while ( i < n ) i = i + 1 call s_to_r8 ( s(ilo:), rvec(i), ierror, length ) if ( ierror /= 0 ) then ierror = -i exit end if ilo = ilo + length end do return end subroutine s_word_count ( s, nword ) !*****************************************************************************80 ! !! S_WORD_COUNT counts the number of "words" in a string. ! ! Licensing: ! ! This code is distributed under the GNU LGPL license. ! ! Modified: ! ! 14 April 1999 ! ! Author: ! ! John Burkardt ! ! Parameters: ! ! Input, character ( len = * ) S, the string to be examined. ! ! Output, integer ( kind = 4 ) NWORD, the number of "words" in the string. ! Words are presumed to be separated by one or more blanks. ! implicit none logical blank integer ( kind = 4 ) i integer ( kind = 4 ) lens integer ( kind = 4 ) nword character ( len = * ) s nword = 0 lens = len ( s ) if ( lens <= 0 ) then return end if blank = .true. do i = 1, lens if ( s(i:i) == ' ' ) then blank = .true. else if ( blank ) then nword = nword + 1 blank = .false. end if end do return end subroutine swapec ( i, top, btri, bedg, point_num, point_xy, tri_num, & tri_vert, tri_nabe, stack, ierr ) !*****************************************************************************80 ! !! SWAPEC swaps diagonal edges until all triangles are Delaunay. ! ! Discussion: ! ! The routine swaps diagonal edges in a 2D triangulation, based on ! the empty circumcircle criterion, until all triangles are Delaunay, ! given that I is the index of the new vertex added to the triangulation. ! ! Licensing: ! ! This code is distributed under the GNU LGPL license. ! ! Modified: ! ! 14 July 2001 ! ! Author: ! ! Original FORTRAN77 version by Barry Joe. ! FORTRAN90 version by John Burkardt. ! ! Reference: ! ! Barry Joe, ! GEOMPACK - a software package for the generation of meshes ! using geometric algorithms, ! Advances in Engineering Software, ! Volume 13, pages 325-331, 1991. ! ! Parameters: ! ! Input, integer ( kind = 4 ) I, the index of the new vertex. ! ! Input/output, integer ( kind = 4 ) TOP, the index of the top of the stack. ! On output, TOP is zero. ! ! Input/output, integer ( kind = 4 ) BTRI, BEDG; on input, if positive, are ! the triangle and edge indices of a boundary edge whose updated indices ! must be recorded. On output, these may be updated because of swaps. ! ! Input, integer ( kind = 4 ) POINT_NUM, the number of points. ! ! Input, real ( kind = 8 ) POINT_XY(2,POINT_NUM), the coordinates ! of the points. ! ! Input, integer ( kind = 4 ) TRI_NUM, the number of triangles. ! ! Input/output, integer ( kind = 4 ) TRI_VERT(3,TRI_NUM), the triangle ! incidence list. May be updated on output because of swaps. ! ! Input/output, integer ( kind = 4 ) TRI_NABE(3,TRI_NUM), the triangle ! neighbor list; negative values are used for links of the counter-clockwise ! linked list of boundary edges; May be updated on output because of swaps. ! ! LINK = -(3*I + J-1) where I, J = triangle, edge index. ! ! Workspace, integer ( kind = 4 ) STACK(MAXST); on input, entries 1 through ! TOP contain the indices of initial triangles (involving vertex I) ! put in stack; the edges opposite I should be in interior; entries ! TOP+1 through MAXST are used as a stack. ! ! Output, integer IERR is set to 8 for abnormal return. ! implicit none integer ( kind = 4 ) point_num integer ( kind = 4 ) tri_num integer ( kind = 4 ) a integer ( kind = 4 ) b integer ( kind = 4 ) bedg integer ( kind = 4 ) btri integer ( kind = 4 ) c integer ( kind = 4 ) diaedg integer ( kind = 4 ) e integer ( kind = 4 ) ee integer ( kind = 4 ) em1 integer ( kind = 4 ) ep1 integer ( kind = 4 ) f integer ( kind = 4 ) fm1 integer ( kind = 4 ) fp1 integer ( kind = 4 ) i integer ( kind = 4 ) ierr integer ( kind = 4 ) i4_wrap integer ( kind = 4 ) l real ( kind = 8 ) point_xy(2,point_num) integer ( kind = 4 ) r integer ( kind = 4 ) s integer ( kind = 4 ) stack(point_num) integer ( kind = 4 ) swap integer ( kind = 4 ) t integer ( kind = 4 ) top integer ( kind = 4 ) tri_nabe(3,tri_num) integer ( kind = 4 ) tri_vert(3,tri_num) integer ( kind = 4 ) tt integer ( kind = 4 ) u real ( kind = 8 ) x real ( kind = 8 ) y ! ! Determine whether triangles in stack are Delaunay, and swap ! diagonal edge of convex quadrilateral if not. ! x = point_xy(1,i) y = point_xy(2,i) do if ( top <= 0 ) then exit end if t = stack(top) top = top - 1 if ( tri_vert(1,t) == i ) then e = 2 b = tri_vert(3,t) else if ( tri_vert(2,t) == i ) then e = 3 b = tri_vert(1,t) else e = 1 b = tri_vert(2,t) end if a = tri_vert(e,t) u = tri_nabe(e,t) if ( tri_nabe(1,u) == t ) then f = 1 c = tri_vert(3,u) else if ( tri_nabe(2,u) == t ) then f = 2 c = tri_vert(1,u) else f = 3 c = tri_vert(2,u) end if swap = diaedg ( x, y, point_xy(1,a), point_xy(2,a), point_xy(1,c), & point_xy(2,c), point_xy(1,b), point_xy(2,b) ) if ( swap == 1 ) then em1 = i4_wrap ( e - 1, 1, 3 ) ep1 = i4_wrap ( e + 1, 1, 3 ) fm1 = i4_wrap ( f - 1, 1, 3 ) fp1 = i4_wrap ( f + 1, 1, 3 ) tri_vert(ep1,t) = c tri_vert(fp1,u) = i r = tri_nabe(ep1,t) s = tri_nabe(fp1,u) tri_nabe(ep1,t) = u tri_nabe(fp1,u) = t tri_nabe(e,t) = s tri_nabe(f,u) = r if ( 0 < tri_nabe(fm1,u) ) then top = top + 1 stack(top) = u end if if ( 0 < s ) then if ( tri_nabe(1,s) == u ) then tri_nabe(1,s) = t else if ( tri_nabe(2,s) == u ) then tri_nabe(2,s) = t else tri_nabe(3,s) = t end if top = top + 1 if ( point_num < top ) then ierr = 8 return end if stack(top) = t else if ( u == btri .and. fp1 == bedg ) then btri = t bedg = e end if l = - ( 3 * t + e - 1 ) tt = t ee = em1 do while ( 0 < tri_nabe(ee,tt) ) tt = tri_nabe(ee,tt) if ( tri_vert(1,tt) == a ) then ee = 3 else if ( tri_vert(2,tt) == a ) then ee = 1 else ee = 2 end if end do tri_nabe(ee,tt) = l end if if ( 0 < r ) then if ( tri_nabe(1,r) == t ) then tri_nabe(1,r) = u else if ( tri_nabe(2,r) == t ) then tri_nabe(2,r) = u else tri_nabe(3,r) = u end if else if ( t == btri .and. ep1 == bedg ) then btri = u bedg = f end if l = - ( 3 * u + f - 1 ) tt = u ee = fm1 do while ( 0 < tri_nabe(ee,tt) ) tt = tri_nabe(ee,tt) if ( tri_vert(1,tt) == b ) then ee = 3 else if ( tri_vert(2,tt) == b ) then ee = 1 else ee = 2 end if end do tri_nabe(ee,tt) = l end if end if end do return end subroutine timestamp ( ) !*****************************************************************************80 ! !! TIMESTAMP prints the current YMDHMS date as a time stamp. ! ! Example: ! ! 31 May 2001 9:45:54.872 AM ! ! Licensing: ! ! This code is distributed under the GNU LGPL license. ! ! Modified: ! ! 18 May 2013 ! ! Author: ! ! John Burkardt ! ! Parameters: ! ! None ! implicit none character ( len = 8 ) ampm integer ( kind = 4 ) d integer ( kind = 4 ) h integer ( kind = 4 ) m integer ( kind = 4 ) mm character ( len = 9 ), parameter, dimension(12) :: month = (/ & 'January ', 'February ', 'March ', 'April ', & 'May ', 'June ', 'July ', 'August ', & 'September', 'October ', 'November ', 'December ' /) integer ( kind = 4 ) n integer ( kind = 4 ) s integer ( kind = 4 ) values(8) integer ( kind = 4 ) y call date_and_time ( values = values ) y = values(1) m = values(2) d = values(3) h = values(5) n = values(6) s = values(7) mm = values(8) if ( h < 12 ) then ampm = 'AM' else if ( h == 12 ) then if ( n == 0 .and. s == 0 ) then ampm = 'Noon' else ampm = 'PM' end if else h = h - 12 if ( h < 12 ) then ampm = 'PM' else if ( h == 12 ) then if ( n == 0 .and. s == 0 ) then ampm = 'Midnight' else ampm = 'AM' end if end if end if write ( *, '(i2,1x,a,1x,i4,2x,i2,a1,i2.2,a1,i2.2,a1,i3.3,1x,a)' ) & d, trim ( month(m) ), y, h, ':', n, ':', s, '.', mm, trim ( ampm ) return end subroutine tri_augment ( v_num, nodtri, v_inf ) !*****************************************************************************80 ! !! TRI_AUGMENT augments the triangle data using vertices at infinity. ! ! Discussion: ! ! The algorithm simply looks at the list of triangle edges stored ! in NODTRI, and determines which edges, of the form (P1,P2), do ! not have a matching (P2,P1) occurrence. These correspond to ! boundary edges of the convex hull. To simplify our computations, ! we adjust the NODTRI array to accommodate an extra triangle with ! one vertex at infinity for each such unmatched edge. ! ! The algorithm used here is ruinously inefficient for large V_NUM. ! Assuming that this data structure modification is the way to go, ! the routine should be rewritten to determine the boundary edges ! more efficiently. ! ! The fictitious vertices at infinity show up in the augmenting ! rows of the NODTRI array with negative indices. ! ! Licensing: ! ! This code is distributed under the GNU LGPL license. ! ! Modified: ! ! 25 August 2003 ! ! Author: ! ! John Burkardt ! ! Parameters: ! ! Input, integer ( kind = 4 ) V_NUM, the number of Voronoi vertices. ! ! Input/output, integer ( kind = 4 ) NODTRI(3,*), the list of nodes that ! comprise each Delaunay triangle. On input, there are V_NUM ! sets of this data. On output, for every pair of nodes (P1,P2) ! for which the pair (P2,P1) does not occur, an augmenting triangle ! has been created with exactly this edge (plus a vertex at infinity). ! On output, there are V_NUM + V_INF sets of data. ! ! Output, integer ( kind = 4 ) V_INF, the number of augmenting triangles and ! vertices at infinity that were created. ! implicit none logical found integer ( kind = 4 ) i integer ( kind = 4 ) i4_wrap integer ( kind = 4 ) i2 integer ( kind = 4 ) ip1 integer ( kind = 4 ) nodtri(3,*) integer ( kind = 4 ) s integer ( kind = 4 ) s2 integer ( kind = 4 ) t integer ( kind = 4 ) t2 integer ( kind = 4 ) v integer ( kind = 4 ) v_inf integer ( kind = 4 ) v_num integer ( kind = 4 ) v2 v_inf = 0 do v = 1, v_num do i = 1, 3 s = nodtri(i,v) ip1 = i4_wrap ( i + 1, 1, 3 ) t = nodtri(ip1,v) found = .false. do v2 = 1, v_num do i2 = 1, 3 s2 = nodtri(i2,v2) ip1 = i4_wrap ( i2 + 1, 1, 3 ) t2 = nodtri(ip1,v2) if ( s == t2 .and. t == s2 ) then found = .true. exit end if end do if ( found ) then exit end if end do if ( .not. found ) then v_inf = v_inf + 1 nodtri(1:3,v_num+v_inf) = (/ -v_inf, t, s /) end if end do end do write ( *, '(a)' ) ' ' write ( *, '(a)' ) 'TRI_AUGMENT:' write ( *, '(a,i8)' ) ' Number of boundary triangles = ', v_inf return end subroutine triangle_area_2d ( t, area ) !*****************************************************************************80 ! !! TRIANGLE_AREA_2D computes the area of a triangle in 2D. ! ! Discussion: ! ! If the triangle's vertices are given in counter clockwise order, ! the area will be positive. If the triangle's vertices are given ! in clockwise order, the area will be negative! ! ! An earlier version of this routine always returned the absolute ! value of the computed area. I am convinced now that that is ! a less useful result! For instance, by returning the signed ! area of a triangle, it is possible to easily compute the area ! of a nonconvex polygon as the sum of the (possibly negative) ! areas of triangles formed by node 1 and successive pairs of vertices. ! ! Licensing: ! ! This code is distributed under the GNU LGPL license. ! ! Modified: ! ! 17 October 2005 ! ! Author: ! ! John Burkardt ! ! Parameters: ! ! Input, real ( kind = 8 ) T(2,3), the triangle vertices. ! ! Output, real ( kind = 8 ) AREA, the area of the triangle. ! implicit none integer ( kind = 4 ), parameter :: dim_num = 2 real ( kind = 8 ) area real ( kind = 8 ) t(dim_num,3) area = 0.5D+00 * ( & t(1,1) * ( t(2,2) - t(2,3) ) & + t(1,2) * ( t(2,3) - t(2,1) ) & + t(1,3) * ( t(2,1) - t(2,2) ) ) return end subroutine triangle_circumcenter_2d ( t, center ) !*****************************************************************************80 ! !! TRIANGLE_CIRCUMCENTER_2D computes the circumcenter of a triangle in 2D. ! ! Discussion: ! ! The circumcenter of a triangle is the center of the circumcircle, the ! circle that passes through the three vertices of the triangle. ! ! The circumcircle contains the triangle, but it is not necessarily the ! smallest triangle to do so. ! ! If all angles of the triangle are no greater than 90 degrees, then ! the center of the circumscribed circle will lie inside the triangle. ! Otherwise, the center will lie outside the triangle. ! ! The circumcenter is the intersection of the perpendicular bisectors ! of the sides of the triangle. ! ! In geometry, the circumcenter of a triangle is often symbolized by "O". ! ! Licensing: ! ! This code is distributed under the GNU LGPL license. ! ! Modified: ! ! 18 February 2005 ! ! Author: ! ! John Burkardt ! ! Parameters: ! ! Input, real ( kind = 8 ) T(2,3), the triangle vertices. ! ! Output, real ( kind = 8 ) CENTER(2), the circumcenter of the triangle. ! implicit none integer ( kind = 4 ), parameter :: dim_num = 2 real ( kind = 8 ) asq real ( kind = 8 ) bot real ( kind = 8 ) center(dim_num) real ( kind = 8 ) csq real ( kind = 8 ) t(dim_num,3) real ( kind = 8 ) top(dim_num) asq = ( t(1,2) - t(1,1) )**2 + ( t(2,2) - t(2,1) )**2 csq = ( t(1,3) - t(1,1) )**2 + ( t(2,3) - t(2,1) )**2 top(1) = ( t(2,2) - t(2,1) ) * csq - ( t(2,3) - t(2,1) ) * asq top(2) = - ( t(1,2) - t(1,1) ) * csq + ( t(1,3) - t(1,1) ) * asq bot = ( t(2,2) - t(2,1) ) * ( t(1,3) - t(1,1) ) & - ( t(2,3) - t(2,1) ) * ( t(1,2) - t(1,1) ) center(1:2) = t(1:2,1) + 0.5D+00 * top(1:2) / bot return end subroutine vbedg ( x, y, point_num, point_xy, tri_num, tri_vert, tri_nabe, & ltri, ledg, rtri, redg ) !*****************************************************************************80 ! !! VBEDG determines which boundary edges are visible to a point. ! ! Discussion: ! ! The point (X,Y) is assumed to be outside the convex hull of the ! region covered by the 2D triangulation. ! ! Licensing: ! ! This code is distributed under the GNU LGPL license. ! ! Modified: ! ! 25 August 2001 ! ! Author: ! ! Original FORTRAN77 version by Barry Joe. ! FORTRAN90 version by John Burkardt. ! ! Reference: ! ! Barry Joe, ! GEOMPACK - a software package for the generation of meshes ! using geometric algorithms, ! Advances in Engineering Software, ! Volume 13, pages 325-331, 1991. ! ! Parameters: ! ! Input, real ( kind = 8 ) X, Y, a point outside the convex hull ! of the current triangulation. ! ! Input, integer ( kind = 4 ) POINT_NUM, the number of points. ! ! Input, real ( kind = 8 ) POINT_XY(2,POINT_NUM), the vertices. ! ! Input, integer ( kind = 4 ) TRI_NUM, the number of triangles. ! ! Input, integer ( kind = 4 ) TRI_VERT(3,TRI_NUM), the triangle ! incidence list. ! ! Input, integer ( kind = 4 ) TRI_NABE(3,TRI_NUM), the triangle neighbor ! list; negative values are used for links of a counter clockwise linked list ! of boundary edges; ! LINK = -(3*I + J-1) where I, J = triangle, edge index. ! ! Input/output, integer ( kind = 4 ) LTRI, LEDG. If LTRI /= 0 then these ! values are assumed to be already computed and are not changed, else they ! are updated. On output, LTRI is the index of boundary triangle to the left ! of the leftmost boundary triangle visible from (X,Y), and LEDG is the ! boundary edge of triangle LTRI to the left of the leftmost boundary edge ! visible from (X,Y). 1 <= LEDG <= 3. ! ! Input/output, integer ( kind = 4 ) RTRI. On input, the index of the ! boundary triangle to begin the search at. On output, the index of the ! rightmost boundary triangle visible from (X,Y). ! ! Input/output, integer ( kind = 4 ) REDG, the edge of triangle RTRI that is ! visible from (X,Y). 1 <= REDG <= 3. ! implicit none integer ( kind = 4 ) point_num integer ( kind = 4 ) tri_num integer ( kind = 4 ) a integer ( kind = 4 ) b integer ( kind = 4 ) e integer ( kind = 4 ) i4_wrap integer ( kind = 4 ) l logical ldone integer ( kind = 4 ) ledg integer ( kind = 4 ) lr integer ( kind = 4 ) lrline integer ( kind = 4 ) ltri real ( kind = 8 ) point_xy(2,point_num) integer ( kind = 4 ) redg integer ( kind = 4 ) rtri integer ( kind = 4 ) t integer ( kind = 4 ) tri_nabe(3,tri_num) integer ( kind = 4 ) tri_vert(3,tri_num) real ( kind = 8 ) x real ( kind = 8 ) y ! ! Find the rightmost visible boundary edge using links, then possibly ! leftmost visible boundary edge using triangle neighbor information. ! if ( ltri == 0 ) then ldone = .false. ltri = rtri ledg = redg else ldone = .true. end if do l = -tri_nabe(redg,rtri) t = l / 3 e = mod ( l, 3 ) + 1 a = tri_vert(e,t) if ( e <= 2 ) then b = tri_vert(e+1,t) else b = tri_vert(1,t) end if lr = lrline ( x, y, point_xy(1,a), point_xy(2,a), point_xy(1,b), & point_xy(2,b), 0.0D+00 ) if ( lr <= 0 ) then exit end if rtri = t redg = e end do if ( ldone ) then return end if t = ltri e = ledg do b = tri_vert(e,t) e = i4_wrap ( e-1, 1, 3 ) do while ( 0 < tri_nabe(e,t) ) t = tri_nabe(e,t) if ( tri_vert(1,t) == b ) then e = 3 else if ( tri_vert(2,t) == b ) then e = 1 else e = 2 end if end do a = tri_vert(e,t) lr = lrline ( x, y, point_xy(1,a), point_xy(2,a), point_xy(1,b), & point_xy(2,b), 0.0D+00 ) if ( lr <= 0 ) then exit end if end do ltri = t ledg = e return end subroutine voronoi_data ( g_num, g_xy, g_degree, g_start, g_face, v_num, & v_xy, i_num, i_xy ) !*****************************************************************************80 ! !! VORONOI_DATA returns data defining the Voronoi diagram. ! ! Discussion: ! ! The routine first determines the Delaunay triangulation. ! ! The Voronoi diagram is then determined from this information. ! ! In particular, the circumcenter of each Delaunay triangle ! is a vertex of a Voronoi polygon. ! ! Licensing: ! ! This code is distributed under the GNU LGPL license. ! ! Modified: ! ! 08 February 2005 ! ! Author: ! ! John Burkardt ! ! Parameters: ! ! Input, integer ( kind = 4 ) G_NUM, the number of generators. ! ! Input, real ( kind = 8 ) G_XY(2,G_NUM), the point coordinates. ! ! Output, integer ( kind = 4 ) G_DEGREE(G_NUM), the degree of each ! Voronoi cell. ! ! Output, integer ( kind = 4 ) G_START(G_NUM), the index in G_FACE of the ! first vertex at which to begin a traversal of the boundary of the ! cell associated with each point. ! ! Output, integer ( kind = 4 ) G_FACE(6*G_NUM), the sequence of vertices to ! be used in a traversal of the boundary of the cell associated with each ! point. ! ! Output, integer ( kind = 4 ) V_NUM, the number of vertices of the Voronoi ! diagram. ! ! Output, real ( kind = 8 ) V_XY(2,V_NUM), the coordinates of the vertices ! of the Voronoi diagram. ! ! Output, integer ( kind = 4 ) I_NUM, the number of vertices at infinity ! of the Voronoi diagram. ! ! Output, real ( kind = 8 ) I_XY(2,I_NUM), the direction of the ! vertices at infinity. ! implicit none integer ( kind = 4 ) g_num integer ( kind = 4 ), parameter :: dim_num = 2 real ( kind = 8 ) area integer ( kind = 4 ) count logical, parameter :: debug = .true. integer ( kind = 4 ) g integer ( kind = 4 ) g_degree(g_num) integer ( kind = 4 ) g_face(6*g_num) integer ( kind = 4 ) g_next integer ( kind = 4 ) g_start(g_num) real ( kind = 8 ) g_xy(dim_num,g_num) integer ( kind = 4 ) i integer ( kind = 4 ) i_num integer ( kind = 4 ) i4_wrap real ( kind = 8 ) i_xy(dim_num,g_num) integer ( kind = 4 ) i1 integer ( kind = 4 ) i2 integer ( kind = 4 ) i3 integer ( kind = 4 ) ix1 integer ( kind = 4 ) ix2 integer ( kind = 4 ) j integer ( kind = 4 ) jp1 integer ( kind = 4 ) k integer ( kind = 4 ) nodtri(3,2*g_num) integer ( kind = 4 ) s integer ( kind = 4 ) sp1 integer ( kind = 4 ) s_next real ( kind = 8 ) t(dim_num,3) integer ( kind = 4 ) tnbr(3,2*g_num) integer ( kind = 4 ) v integer ( kind = 4 ) v_inf integer ( kind = 4 ) v_next integer ( kind = 4 ) v_num integer ( kind = 4 ) v_old integer ( kind = 4 ) v_save real ( kind = 8 ) v_xy(dim_num,2*g_num) real ( kind = 8 ) x1 real ( kind = 8 ) x2 real ( kind = 8 ) y1 real ( kind = 8 ) y2 ! ! Compute the Delaunay triangulation. ! call dtris2 ( g_num, g_xy, v_num, nodtri, tnbr ) ! ! Compute and print the areas of the finite triangles. ! write ( *, '(a)' ) ' ' write ( *, '(a)' ) ' Triangle Area' write ( *, '(a)' ) ' ' do v = 1, v_num i1 = nodtri(1,v) i2 = nodtri(2,v) i3 = nodtri(3,v) t(1:dim_num,1:3) = reshape ( (/ & g_xy(1:2,i1), g_xy(1:2,i2), g_xy(1:2,i3) /), (/ dim_num, 3 /) ) call triangle_area_2d ( t, area ) write ( *, '(2x,i8,2x,g14.6)' ) v, area end do ! ! Extend the NODTRI data structure, adding fictitious vertices at infinity, ! so that the Delaunay triangulation can be regarded as covering the ! entire plane. ! call tri_augment ( v_num, nodtri, v_inf ) if ( debug ) then write ( *, '(a)' ) ' ' write ( *, '(a)' ) ' The generators that form each Delaunay triangle:' write ( *, '(a)' ) ' (Negative values are fictitious nodes at infinity.)' write ( *, '(a)' ) ' ' call i4mat_transpose_print ( 3, v_num+v_inf, nodtri, ' Triangle nodes:' ) end if ! ! Negative entries in TNBR indicate a semi-infinite Voronoi side. ! However, DTRIS2 uses a peculiar numbering. Renumber them. ! i_num = 0 do v = 1, v_num do i = 1, 3 if ( tnbr(i,v) < 0 ) then i_num = i_num + 1 tnbr(i,v) = -i_num end if end do end do if ( debug ) then write ( *, '(a)' ) ' ' write ( *, '(a)' ) ' Neighboring triangles of each Delaunay triangle:' write ( *, '(a)' ) ' Negative values indicate no finite neighbor.' write ( *, '(a)' ) ' ' call i4mat_transpose_print ( 3, v_num, tnbr, ' Neighbor triangles:' ) end if ! ! Determine the degree of each cell. ! g_degree(1:g_num) = 0 do j = 1, v_num + v_inf do i = 1, 3 k = nodtri(i,j) if ( 0 < k ) then g_degree(k) = g_degree(k) + 1 end if end do end do call i4vec_print ( g_num, g_degree, ' Voronoi cell degrees' ) ! ! Each (finite) Delaunay triangle contains a vertex of the Voronoi polygon, ! at the triangle's circumcenter. ! do v = 1, v_num i1 = nodtri(1,v) i2 = nodtri(2,v) i3 = nodtri(3,v) t(1:dim_num,1:3) = reshape ( (/ & g_xy(1:2,i1), g_xy(1:2,i2), g_xy(1:2,i3) /), (/ dim_num, 3 /) ) call triangle_circumcenter_2d ( t, v_xy(1:2,v) ) end do call r8mat_transpose_print ( dim_num, v_num, v_xy, ' The Voronoi vertices:' ) ! ! For each generator G: ! Determine if its region is infinite. ! Find a Delaunay triangle containing G. ! Seek another triangle containing the next node in that triangle. ! count = 0 g_start(1:g_num) = 0 do g = 1, g_num v_next = 0 do v = 1, v_num + v_inf do s = 1, 3 if ( nodtri(s,v) == g ) then v_next = v s_next = s exit end if end do if ( v_next /= 0 ) then exit end if end do v_save = v_next do s_next = i4_wrap ( s_next + 1, 1, 3 ) g_next = nodtri(s_next,v_next) if ( g_next == g ) then s_next = i4_wrap ( s_next + 1, 1, 3 ) g_next = nodtri(s_next,v_next) end if v_old = v_next v_next = 0 do v = 1, v_num + v_inf if ( v == v_old ) then cycle end if do s = 1, 3 if ( nodtri(s,v) == g ) then sp1 = i4_wrap ( s + 1, 1, 3 ) if ( nodtri(sp1,v) == g_next ) then v_next = v s_next = sp1 exit end if sp1 = i4_wrap ( s + 2, 1, 3 ) if ( nodtri(sp1,v) == g_next ) then v_next = v s_next = sp1 exit end if end if end do if ( v_next /= 0 ) then exit end if end do if ( v_next == v_save ) then exit end if if ( v_next == 0 ) then v_next = v_old exit end if end do ! ! Now, starting in the current triangle, V_NEXT, cycle again, ! and copy the list of nodes into the array. ! v_save = v_next count = count + 1 g_start(g) = count g_face(count) = v_next do s_next = i4_wrap ( s_next + 1, 1, 3 ) g_next = nodtri(s_next,v_next) if ( g_next == g ) then s_next = i4_wrap ( s_next + 1, 1, 3 ) g_next = nodtri(s_next,v_next) end if v_old = v_next v_next = 0 do v = 1, v_num + v_inf if ( v == v_old ) then cycle end if do s = 1, 3 if ( nodtri(s,v) == g ) then sp1 = i4_wrap ( s + 1, 1, 3 ) if ( nodtri(sp1,v) == g_next ) then v_next = v s_next = sp1 exit end if sp1 = i4_wrap ( s + 2, 1, 3 ) if ( nodtri(sp1,v) == g_next ) then v_next = v s_next = sp1 exit end if end if end do if ( v_next /= 0 ) then exit end if end do if ( v_next == v_save ) then exit end if if ( v_next == 0 ) then exit end if count = count + 1 g_face(count) = v_next end do end do ! ! Mark all the vertices at infinity with a negative sign, ! so that the data in G_FACE is easier to interpret. ! do i = 1, count if ( v_num < g_face(i) ) then g_face(i) = -g_face(i) end if end do ! ! For each (finite) Delaunay triangle, I ! For each side J, ! do i = 1, v_num do j = 1, 3 k = tnbr(j,i) ! If there is no neighboring triangle on that side, ! extend a line from the circumcenter of I in the direction of the ! outward normal to that side. This is an infinite edge of ! an infinite Voronoi polygon. ! if ( k < 0 ) then ix1 = nodtri(j,i) x1 = g_xy(1,ix1) y1 = g_xy(2,ix1) jp1 = i4_wrap ( j+1, 1, 3 ) ix2 = nodtri(jp1,i) x2 = g_xy(1,ix2) y2 = g_xy(2,ix2) ! ! Compute the direction I_XY(1:2,-K). ! call line_exp_normal_2d ( g_xy(1:2,ix1), g_xy(1:2,ix2), i_xy(1:2,-k) ) end if end do end do return end