program main !*****************************************************************************80 ! !! QUAD_MESH_RCM applies the RCM reordering to a mesh of quadrilaterals. ! ! Discussion: ! ! The user supplies a node file and a element file, containing ! the coordinates of the nodes, and the indices of the nodes that ! make up each element. ! ! The program reads the data, computes the adjacency information, ! carries out the RCM algorithm to get the permutation, applies ! the permutation to the nodes and elements, and writes out ! new node and element files that correspond to the RCM permutation. ! ! Note that node data is normally two dimensional, that is, ! each node has an X and Y coordinate. In some applications, it ! may be desirable to specify more information. This program ! will accept node data that includes DIM_NUM entries on each line, ! as long as DIM_NUM is the same for each entry. ! ! Usage: ! ! quad_mesh_rcm prefix ! ! where 'prefix' is the common filename prefix: ! ! * prefix_nodes.txt contains the node coordinates, ! * prefix_elements.txt contains the element definitions. ! * prefix_rcm_nodes.txt will contain the RCM node coordinates, ! * prefix_rcm_elements.txt will contain the RCM element definitions. ! ! Licensing: ! ! This code is distributed under the MIT license. ! ! Modified: ! ! 30 September 2009 ! ! Author: ! ! John Burkardt ! implicit none integer, parameter :: rk = kind ( 1.0D+00 ) integer, allocatable, dimension ( : ) :: adj integer adj_bandwidth integer adj_num integer adj_perm_bandwidth integer, allocatable, dimension ( : ) :: adj_row integer arg_num integer bandwidth integer dim_num character ( len = 255 ) :: element_filename = ' ' integer element_num integer element_order integer, allocatable, dimension ( :, : ) :: element_neighbor integer, allocatable, dimension ( :, : ) :: element_node character ( len = 255 ) :: element_rcm_filename = ' ' integer i integer iarg integer iargc integer j integer node character ( len = 255 ) :: node_filename = ' ' integer node_num character ( len = 255 ) :: node_rcm_filename = ' ' real ( kind = rk ), allocatable, dimension ( :, : ) :: node_xy integer, allocatable, dimension ( : ) :: perm integer, allocatable, dimension ( : ) :: perm_inv character ( len = 255 ) prefix call timestamp ( ) write ( *, '(a)' ) ' ' write ( *, '(a)' ) 'QUAD_MESH_RCM:' write ( *, '(a)' ) ' FORTRAN90 version' write ( *, '(a)' ) ' ' write ( *, '(a)' ) ' Read a node dataset of NODE_NUM points in 2 dimensions.' write ( *, '(a)' ) ' Read an associated quad mesh dataset of ELEMENT_NUM ' write ( *, '(a)' ) ' elements using ELEMENT_ORDER nodes.' write ( *, '(a)' ) ' ' write ( *, '(a)' ) ' Apply the RCM reordering (Reverse Cuthill-McKee).' write ( *, '(a)' ) ' ' write ( *, '(a)' ) ' Reorder the data and write it out to files.' ! ! Get the number of command line arguments. ! arg_num = iargc ( ) ! ! The command line argument is the common filename prefix. ! if ( 1 <= arg_num ) then iarg = 1 call getarg ( iarg, prefix ) else write ( *, '(a)' ) ' ' write ( *, '(a)' ) 'QUAD_MESH_RCM:' write ( *, '(a)' ) ' Please enter the filename prefix.' read ( *, '(a)' ) prefix end if ! ! Create the filenames. ! node_filename = trim ( prefix ) // '_nodes.txt' element_filename = trim ( prefix ) // '_elements.txt' node_rcm_filename = trim ( prefix ) // '_rcm_nodes.txt' element_rcm_filename = trim ( prefix ) // '_rcm_elements.txt' ! ! Read the node data. ! call r8mat_header_read ( node_filename, dim_num, node_num ) write ( *, '(a)' ) ' ' write ( *, '(a)' ) ' Read the header of "' & // trim ( node_filename ) //'".' write ( *, '(a)' ) ' ' write ( *, '(a,i8)' ) ' Spatial dimension DIM_NUM = ', dim_num write ( *, '(a,i8)' ) ' Number of nodes NODE_NUM = ', node_num allocate ( node_xy(1:dim_num,1:node_num) ) call r8mat_data_read ( node_filename, dim_num, node_num, node_xy ) write ( *, '(a)' ) ' ' write ( *, '(a)' ) ' Read the data in "' & // trim ( node_filename ) //'".' call r8mat_transpose_print_some ( dim_num, node_num, node_xy, 1, 1, & dim_num, 5, ' Coordinates of first 5 nodes:' ) ! ! Read the mesh data. ! call i4mat_header_read ( element_filename, element_order, & element_num ) write ( *, '(a)' ) ' ' write ( *, '(a)' ) ' Read the header of "' & // trim ( element_filename ) //'".' write ( *, '(a)' ) ' ' write ( *, '(a,i8)' ) ' Element order ELEMENT_ORDER = ', element_order write ( *, '(a,i8)' ) ' Number of elements ELEMENT_NUM = ', element_num if ( element_order /= 4 ) then write ( *, '(a)' ) ' ' write ( *, '(a)' ) 'QUAD_MESH_RCM - Fatal error!' write ( *, '(a)' ) ' This program can only handle quadrilaterals' write ( *, '(a)' ) ' of order 4.' write ( *, '(a)' ) ' ' write ( *, '(a)' ) ' The input mesh seems to have' write ( *, '(a,i8)' ) ' order = ', element_order stop end if allocate ( element_node(1:element_order,1:element_num) ) call i4mat_data_read ( element_filename, element_order, & element_num, element_node ) write ( *, '(a)' ) ' ' write ( *, '(a)' ) ' Read the data in "' & // trim ( element_filename ) //'".' call i4mat_transpose_print_some ( element_order, element_num, & element_node, 1, 1, element_order, 5, ' First 5 elements:' ) ! ! Create the element neighbor array. ! allocate ( element_neighbor(1:4,1:element_num) ) call neighbor_elements_q4_mesh ( element_num, element_node, element_neighbor ) ! ! Count the number of adjacencies, and set up the ADJ_ROW ! adjacency pointer array. ! allocate ( adj_row(1:node_num+1) ) call adj_size_q4_mesh ( node_num, element_num, element_node, & element_neighbor, adj_num, adj_row ) ! ! Set up the ADJ adjacency array. ! allocate ( adj(1:adj_num) ) call adj_set_q4_mesh ( node_num, element_num, element_node, & element_neighbor, adj_num, adj_row, adj ) ! ! Determine the bandwidth. ! bandwidth = adj_bandwidth ( node_num, adj_num, adj_row, adj ) write ( *, '(a)' ) ' ' write ( *, '(a,i8)' ) ' ADJ bandwidth = ', bandwidth ! ! Compute the RCM permutation. ! allocate ( perm(1:node_num) ) call genrcm ( node_num, adj_num, adj_row, adj, perm ) allocate ( perm_inv(1:node_num) ) call perm_inverse3 ( node_num, perm, perm_inv ) bandwidth = adj_perm_bandwidth ( node_num, adj_num, adj_row, adj, & perm, perm_inv ) write ( *, '(a)' ) ' ' write ( *, '(a,i8)' ) ' Permuted ADJ bandwidth = ', bandwidth ! ! Permute the nodes according to the permutation vector. ! call r8col_permute ( dim_num, node_num, perm, node_xy ) ! ! Permute the node indices in the element array. ! do j = 1, element_num do i = 1, element_order node = element_node(i,j) element_node(i,j) = perm_inv ( node ) end do end do ! ! Write the nodes. ! call r8mat_write ( node_rcm_filename, dim_num, node_num, node_xy ) write ( *, '(a)' ) ' ' write ( *, '(a)' ) ' Created the node file "' & // trim ( node_rcm_filename ) //'".' call i4mat_write ( element_rcm_filename, element_order, & element_num, element_node ) write ( *, '(a)' ) ' ' write ( *, '(a)' ) ' Created the triangulation file "' & // trim ( element_rcm_filename ) //'".' ! ! Free memory. ! deallocate ( adj ) deallocate ( adj_row ) deallocate ( element_neighbor ) deallocate ( element_node ) deallocate ( node_xy ) deallocate ( perm ) deallocate ( perm_inv ) ! ! Terminate. ! write ( *, '(a)' ) ' ' write ( *, '(a)' ) 'QUAD_MESH_RCM' write ( *, '(a)' ) ' Normal end of execution.' write ( *, '(a)' ) ' ' call timestamp ( ) stop end function adj_bandwidth ( node_num, adj_num, adj_row, adj ) !*****************************************************************************80 ! !! ADJ_BANDWIDTH computes the bandwidth of an adjacency matrix. ! ! Licensing: ! ! This code is distributed under the MIT license. ! ! Modified: ! ! 11 March 2005 ! ! Author: ! ! John Burkardt ! ! Author: ! ! Original FORTRAN77 version by Alan George, Joseph Liu. ! FORTRAN90 version by John Burkardt. ! ! Reference: ! ! Alan George, Joseph Liu, ! Computer Solution of Large Sparse Positive Definite Systems, ! Prentice Hall, 1981. ! ! Parameters: ! ! Input, integer NODE_NUM, the number of nodes. ! ! Input, integer ADJ_NUM, the number of adjacency entries. ! ! Input, integer ADJ_ROW(NODE_NUM+1). Information about ! row I is stored in entries ADJ_ROW(I) through ADJ_ROW(I+1)-1 of ADJ. ! ! Input, integer ADJ(ADJ_NUM), the adjacency structure. ! For each row, it contains the column indices of the nonzero entries. ! ! Output, integer ADJ_BANDWIDTH, the bandwidth of the adjacency ! matrix. ! implicit none integer adj_num integer node_num integer adj(adj_num) integer adj_bandwidth integer adj_row(node_num+1) integer band_hi integer band_lo integer col integer i integer j band_lo = 0 band_hi = 0 do i = 1, node_num do j = adj_row(i), adj_row(i+1) - 1 col = adj(j) band_lo = max ( band_lo, i - col ) band_hi = max ( band_hi, col - i ) end do end do adj_bandwidth = band_lo + 1 + band_hi return end function adj_perm_bandwidth ( node_num, adj_num, adj_row, adj, perm, perm_inv ) !*****************************************************************************80 ! !! ADJ_PERM_BANDWIDTH computes the bandwidth of a permuted adjacency matrix. ! ! Discussion: ! ! The matrix is defined by the adjacency information and a permutation. ! ! The routine also computes the bandwidth and the size of the envelope. ! ! Licensing: ! ! This code is distributed under the MIT license. ! ! Modified: ! ! 11 March 2005 ! ! Author: ! ! John Burkardt ! ! Reference: ! ! Alan George, Joseph Liu, ! Computer Solution of Large Sparse Positive Definite Systems, ! Prentice Hall, 1981. ! ! Parameters: ! ! Input, integer NODE_NUM, the number of nodes. ! ! Input, integer ADJ_NUM, the number of adjacency entries. ! ! Input, integer ADJ_ROW(NODE_NUM+1). Information about ! row I is stored in entries ADJ_ROW(I) through ADJ_ROW(I+1)-1 of ADJ. ! ! Input, integer ADJ(ADJ_NUM), the adjacency structure. ! For each row, it contains the column indices of the nonzero entries. ! ! Input, integer PERM(NODE_NUM), PERM_INV(NODE_NUM), the ! permutation and inverse permutation. ! ! Output, integer ADJ_PERM_BANDWIDTH, the bandwidth of the ! permuted adjacency matrix. ! implicit none integer adj_num integer node_num integer adj(adj_num) integer adj_perm_bandwidth integer adj_row(node_num+1) integer band_hi integer band_lo integer col integer i integer j integer perm(node_num) integer perm_inv(node_num) band_lo = 0 band_hi = 0 do i = 1, node_num do j = adj_row(perm(i)), adj_row(perm(i)+1) - 1 col = perm_inv(adj(j)) band_lo = max ( band_lo, i - col ) band_hi = max ( band_hi, col - i ) end do end do adj_perm_bandwidth = band_lo + 1 + band_hi return end subroutine adj_set_q4_mesh ( node_num, element_num, element_node, & element_neighbor, adj_num, adj_row, adj ) !*****************************************************************************80 ! !! ADJ_SET_Q4_MESH sets adjacencies in a Q4 mesh. ! ! Discussion: ! ! This routine is called to set the adjacency values, after the ! appropriate amount of memory has been set aside for storage. ! ! The mesh is assumed to involve 4-node quadrilaterals. ! ! Two nodes are "adjacent" if they are both nodes in some element. ! Also, a node is considered to be adjacent to itself. ! ! This routine can be used to create the compressed column storage ! for a linear element finite element discretization of ! Poisson's equation in two dimensions. ! ! Diagram: ! ! side 3 ! 4-------3 ! s | | s ! i | | i ! d | | d ! e | | e ! | | ! 4 | | 2 ! | | ! 1-------2 ! ! side 1 ! ! The local node numbering ! ! ! 20-21-22-23-24 ! | | | | | ! | | | | | ! 15-16-17-18-19 ! | | | | | ! | | | | | ! 10-11-12-13-14 ! | | | | | ! | | | | | ! 5--6--7--8--9 ! | | | | | ! | | | | | ! 0--1--2--3--4 ! ! A sample grid. ! ! Licensing: ! ! This code is distributed under the MIT license. ! ! Modified: ! ! 29 September 2009 ! ! Author: ! ! John Burkardt ! ! Parameters ! ! Input, integer NODE_NUM, the number of nodes. ! ! Input, integer ELEMENT_NUM, the number of elements. ! ! Input, integer ELEMENT_NODE(4,ELEMENT_NUM), lists the nodes ! that make up each element in counterclockwise order. ! ! Input, integer ELEMENT_NEIGHBOR(4,ELEMENT_NUM), for each ! side of an element, lists the neighboring element, or -1 if there is ! no neighbor. ! ! Input, integer ADJ_NUM, the number of adjacencies. ! ! Input, integer ADJ_ROW(NODE_NUM+1). Information about ! column J is stored in entries ADJ_ROW(J) through ADJ_ROW(J+1)-1 of ADJ. ! ! Output, integer ADJ(ADJ_NUM), the adjacency information. ! implicit none integer adj_num integer node_num integer element_num integer, parameter :: element_order = 4 integer adj(adj_num) integer adj_row(node_num+1) integer adj_copy(node_num) integer k1 integer k2 integer n1 integer n2 integer n3 integer n4 integer node integer number integer element integer element2 integer element_neighbor(4,element_num) integer element_node(element_order,element_num) adj(1:adj_num) = -1 adj_copy(1:node_num) = adj_row(1:node_num) ! ! Set every node to be adjacent to itself. ! do node = 1, node_num adj(adj_copy(node)) = node adj_copy(node) = adj_copy(node) + 1 end do ! ! Examine each element. ! do element = 1, element_num n1 = element_node(1,element) n2 = element_node(2,element) n3 = element_node(3,element) n4 = element_node(4,element) ! ! Add edges (1,3) and (2,4). There is no need to check for redundancy, ! since this is the only case when these nodes can share an element. ! adj(adj_copy(n1)) = n3 adj_copy(n1) = adj_copy(n1) + 1 adj(adj_copy(n3)) = n1 adj_copy(n3) = adj_copy(n3) + 1 adj(adj_copy(n2)) = n4 adj_copy(n2) = adj_copy(n2) + 1 adj(adj_copy(n4)) = n2 adj_copy(n4) = adj_copy(n4) + 1 ! ! Add edge (1,2) if this is the first occurrence, ! that is, if the edge (1,2) is on a boundary (ELEMENT2 <= 0) ! or if this element is the first of the pair in which the edge ! occurs (ELEMENT < ELEMENT2). ! element2 = element_neighbor(1,element) if ( element2 < 0 .or. element < element2 ) then adj(adj_copy(n1)) = n2 adj_copy(n1) = adj_copy(n1) + 1 adj(adj_copy(n2)) = n1 adj_copy(n2) = adj_copy(n2) + 1 end if ! ! Add edge (2,3). ! element2 = element_neighbor(2,element) if ( element2 < 0 .or. element < element2 ) then adj(adj_copy(n2)) = n3 adj_copy(n2) = adj_copy(n2) + 1 adj(adj_copy(n3)) = n2 adj_copy(n3) = adj_copy(n3) + 1 end if ! ! Add edge (3,4). ! element2 = element_neighbor(3,element) if ( element2 < 0 .or. element < element2 ) then adj(adj_copy(n4)) = n3 adj_copy(n4) = adj_copy(n4) + 1 adj(adj_copy(n3)) = n4 adj_copy(n3) = adj_copy(n3) + 1 end if ! ! Add edge (4,1). ! element2 = element_neighbor(4,element) if ( element2 < 0 .or. element < element2 ) then adj(adj_copy(n1)) = n4 adj_copy(n1) = adj_copy(n1) + 1 adj(adj_copy(n4)) = n1 adj_copy(n4) = adj_copy(n4) + 1 end if end do ! ! Ascending sort the entries for each node. ! do node = 1, node_num k1 = adj_row(node) k2 = adj_row(node+1)-1 number = k2 + 1 - k1 call i4vec_sort_heap_a ( number, adj(k1:k2) ) end do return end subroutine adj_size_q4_mesh ( node_num, element_num, element_node, & element_neighbor, adj_num, adj_row ) !*****************************************************************************80 ! !! ADJ_SIZE_Q4_MESH counts adjacencies in a Q4 mesh. ! ! Discussion: ! ! This routine is called to count the adjacencies, so that the ! appropriate amount of memory can be set aside for storage when ! the adjacency structure is created. ! ! The mesh is assumed to involve 4-node quadrilaterals. ! ! Two nodes are "adjacent" if they are both nodes in some quadrilateral. ! Also, a node is considered to be adjacent to itself. ! ! Diagram: ! ! side 3 ! 4-------3 ! s | | s ! i | | i ! d | | d ! e | | e ! | | ! 4 | | 2 ! | | ! 1-------2 ! ! side 1 ! ! The local node numbering ! ! ! 20-21-22-23-24 ! | | | | | ! | | | | | ! 15-16-17-18-19 ! | | | | | ! | | | | | ! 10-11-12-13-14 ! | | | | | ! | | | | | ! 5--6--7--8--9 ! | | | | | ! | | | | | ! 0--1--2--3--4 ! ! A sample grid. ! ! Licensing: ! ! This code is distributed under the MIT license. ! ! Modified: ! ! 29 September 2009 ! ! Author: ! ! John Burkardt ! ! Parameters ! ! Input, integer NODE_NUM, the number of nodes. ! ! Input, integer ELEMENT_NUM, the number of elements. ! ! Input, integer ELEMENT_NODE(4,ELEMENT_NUM), lists the ! nodes that make up each element, in counterclockwise order. ! ! Input, integer ELEMENT_NEIGHBOR(4,ELEMENT_NUM), for each ! side of a element, lists the neighboring elment, or -1 if there is ! no neighbor. ! ! Output, integer ADJ_NUM, the number of adjacencies. ! ! Output, integer ADJ_ROW(NODE_NUM+1), Information about ! column J is stored in entries ADJ_ROW(J) through ADJ_ROW(J+1)-1 of ADJ. ! implicit none integer element_num integer, parameter :: element_order = 4 integer node_num integer adj_row(node_num+1) integer adj_num integer element integer element_neighbor(4,element_num) integer element_node(element_order,element_num) integer element2 integer n1 integer n2 integer n3 integer n4 integer node adj_num = 0 ! ! Set every node to be adjacent to itself. ! adj_row(1:node_num) = 1 ! ! Examine each element. ! do element = 1, element_num n1 = element_node(1,element) n2 = element_node(2,element) n3 = element_node(3,element) n4 = element_node(4,element) ! ! Add edge (1,3). ! adj_row(n1) = adj_row(n1) + 1 adj_row(n3) = adj_row(n3) + 1 ! ! Add edge (2,4). ! adj_row(n2) = adj_row(n2) + 1 adj_row(n4) = adj_row(n4) + 1 ! ! Add edge (1,2) if this is the first occurrence, ! that is, if the edge (1,2) is on a boundary (ELEMENT2 <= 0) ! or if this element is the first of the pair in which the edge ! occurs (ELEMENT < ELEMENT2). ! element2 = element_neighbor(1,element) if ( element2 < 0 .or. element < element2 ) then adj_row(n1) = adj_row(n1) + 1 adj_row(n2) = adj_row(n2) + 1 end if ! ! Add edge (2,3). ! element2 = element_neighbor(2,element) if ( element2 < 0 .or. element < element2 ) then adj_row(n2) = adj_row(n2) + 1 adj_row(n3) = adj_row(n3) + 1 end if ! ! Add edge (3,4). ! element2 = element_neighbor(3,element) if ( element2 < 0 .or. element < element2 ) then adj_row(n3) = adj_row(n3) + 1 adj_row(n4) = adj_row(n4) + 1 end if ! ! Add edge (4,1). ! element2 = element_neighbor(4,element) if ( element2 < 0 .or. element < element2 ) then adj_row(n4) = adj_row(n4) + 1 adj_row(n1) = adj_row(n1) + 1 end if end do ! ! We used ADJ_ROW to count the number of entries in each column. ! Convert it to pointers into the ADJ array. ! do node = node_num + 1, 2, -1 adj_row(node) = adj_row(node-1) end do adj_row(1) = 1 do node = 2, node_num + 1 adj_row(node) = adj_row(node) + adj_row(node-1) end do ! ! Finally, record the total number of adjacencies. ! adj_num = adj_row(node_num+1) - 1 return end subroutine bandwidth ( element_order, element_num, element_node, ml, mu, m ) !*****************************************************************************80 ! !! BANDWIDTH determines the bandwidth associated with a finite element mesh. ! ! Discussion: ! ! The quantity computed here is the "geometric" bandwidth determined ! by the finite element mesh alone. ! ! If a single finite element variable is associated with each node ! of the mesh, and if the nodes and variables are numbered in the ! same way, then the geometric bandwidth is the same as the bandwidth ! of a typical finite element matrix. ! ! The bandwidth M is defined in terms of the lower and upper bandwidths: ! ! M = ML + 1 + MU ! ! where ! ! ML = maximum distance from any diagonal entry to a nonzero ! entry in the same row, but earlier column, ! ! MU = maximum distance from any diagonal entry to a nonzero ! entry in the same row, but later column. ! ! Because the finite element node adjacency relationship is symmetric, ! we are guaranteed that ML = MU. ! ! Licensing: ! ! This code is distributed under the MIT license. ! ! Modified: ! ! 02 September 2006 ! ! Author: ! ! John Burkardt ! ! Parameters: ! ! Input, integer ELEMENT_ORDER, the order of the elements. ! ! Input, integer ELEMENT_NUM, the number of elements. ! ! Input, integer ELEMENT_NODE(ELEMENT_ORDER,ELEMENT_NUM); ! ELEMENT_NODE(I,J) is the global index of local node I in element J. ! ! Output, integer ML, MU, the lower and upper bandwidths ! of the matrix. ! ! Output, integer M, the bandwidth of the matrix. ! implicit none integer element_num integer element_order integer element integer element_node(element_order,element_num) integer global_i integer global_j integer local_i integer local_j integer m integer ml integer mu ml = 0 mu = 0 do element = 1, element_num do local_i = 1, element_order global_i = element_node(local_i,element) do local_j = 1, element_order global_j = element_node(local_j,element) mu = max ( mu, global_j - global_i ) ml = max ( ml, global_i - global_j ) end do end do end do m = ml + 1 + mu return end subroutine ch_cap ( c ) !*****************************************************************************80 ! !! CH_CAP capitalizes a single character. ! ! Licensing: ! ! This code is distributed under the MIT license. ! ! Modified: ! ! 19 July 1998 ! ! Author: ! ! John Burkardt ! ! Parameters: ! ! Input/output, character C, the character to capitalize. ! implicit none character c integer 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 MIT 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 logical ch_eqi character c1 character c1_cap character c2 character c2_cap 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 MIT license. ! ! Modified: ! ! 04 August 1999 ! ! Author: ! ! John Burkardt ! ! Parameters: ! ! Input, character C, the decimal digit, '0' through '9' or blank ! are legal. ! ! Output, integer DIGIT, the corresponding value. ! If C was 'illegal', then DIGIT is -1. ! implicit none character c integer 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 subroutine degree ( root, adj_num, adj_row, adj, mask, deg, iccsze, ls, & node_num ) !*****************************************************************************80 ! !! DEGREE computes the degrees of the nodes in the connected component. ! ! Discussion: ! ! The connected component is specified by MASK and ROOT. ! Nodes for which MASK is zero are ignored. ! ! Licensing: ! ! This code is distributed under the MIT license. ! ! Modified: ! ! 05 January 2003 ! ! Author: ! ! Original FORTRAN77 version by Alan George, Joseph Liu. ! FORTRAN90 version by John Burkardt ! ! Reference: ! ! Alan George, Joseph Liu, ! Computer Solution of Large Sparse Positive Definite Systems, ! Prentice Hall, 1981. ! ! Parameters: ! ! Input, integer ROOT, the node that defines the connected ! component. ! ! Input, integer ADJ_NUM, the number of adjacency entries. ! ! Input, integer ADJ_ROW(NODE_NUM+1). Information about ! row I is stored in entries ADJ_ROW(I) through ADJ_ROW(I+1)-1 of ADJ. ! ! Input, integer ADJ(ADJ_NUM), the adjacency structure. ! For each row, it contains the column indices of the nonzero entries. ! ! Input, integer MASK(NODE_NUM), is nonzero for those nodes ! which are to be considered. ! ! Output, integer DEG(NODE_NUM), contains, for each node in ! the connected component, its degree. ! ! Output, integer ICCSIZE, the number of nodes in the ! connected component. ! ! Output, integer LS(NODE_NUM), stores in entries 1 through ! ICCSIZE the nodes in the connected component, starting with ROOT, and ! proceeding by levels. ! ! Input, integer NODE_NUM, the number of nodes. ! implicit none integer adj_num integer node_num integer adj(adj_num) integer adj_row(node_num+1) integer deg(node_num) integer i integer iccsze integer ideg integer j integer jstop integer jstrt integer lbegin integer ls(node_num) integer lvlend integer lvsize integer mask(node_num) integer nbr integer node integer root ! ! The sign of ADJ_ROW(I) is used to indicate if node I has been considered. ! ls(1) = root adj_row(root) = -adj_row(root) lvlend = 0 iccsze = 1 ! ! LBEGIN is the pointer to the beginning of the current level, and ! LVLEND points to the end of this level. ! do lbegin = lvlend + 1 lvlend = iccsze ! ! Find the degrees of nodes in the current level, ! and at the same time, generate the next level. ! do i = lbegin, lvlend node = ls(i) jstrt = -adj_row(node) jstop = abs ( adj_row(node+1) ) - 1 ideg = 0 do j = jstrt, jstop nbr = adj(j) if ( mask(nbr) /= 0 ) then ideg = ideg + 1 if ( 0 <= adj_row(nbr) ) then adj_row(nbr) = -adj_row(nbr) iccsze = iccsze + 1 ls(iccsze) = nbr end if end if end do deg(node) = ideg end do ! ! Compute the current level width. ! lvsize = iccsze - lvlend ! ! If the current level width is nonzero, generate another level. ! if ( lvsize == 0 ) then exit end if end do ! ! Reset ADJ_ROW to its correct sign and return. ! do i = 1, iccsze node = ls(i) adj_row(node) = -adj_row(node) end do return end subroutine file_column_count ( input_file_name, 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 MIT license. ! ! Modified: ! ! 21 June 2001 ! ! Author: ! ! John Burkardt ! ! Parameters: ! ! Input, character ( len = * ) INPUT_FILE_NAME, the name of the file. ! ! Output, integer COLUMN_NUM, the number of columns in the file. ! implicit none integer column_num logical got_one character ( len = * ) input_file_name integer input_status integer input_unit character ( len = 255 ) line ! ! Open the file. ! call get_unit ( input_unit ) open ( unit = input_unit, file = input_file_name, status = 'old', & form = 'formatted', access = 'sequential', iostat = input_status ) if ( input_status /= 0 ) then column_num = -1 write ( *, '(a)' ) ' ' write ( *, '(a)' ) 'FILE_COLUMN_COUNT - Fatal error!' write ( *, '(a,i8)' ) ' Could not open the input file "' & // trim ( input_file_name ) // '" on unit ', input_unit return end if ! ! Read one line, but skip blank lines and comment lines. ! got_one = .false. do read ( input_unit, '(a)', iostat = input_status ) line if ( input_status /= 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 = input_status ) line if ( input_status /= 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_file_name, 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 MIT license. ! ! Modified: ! ! 06 March 2003 ! ! Author: ! ! John Burkardt ! ! Parameters: ! ! Input, character ( len = * ) INPUT_FILE_NAME, the name of the input file. ! ! Output, integer ROW_NUM, the number of rows found. ! implicit none integer bad_num integer comment_num integer ierror character ( len = * ) input_file_name integer input_status integer input_unit character ( len = 255 ) line integer record_num integer row_num call get_unit ( input_unit ) open ( unit = input_unit, file = input_file_name, status = 'old', & iostat = input_status ) if ( input_status /= 0 ) then row_num = -1; ierror = 1 write ( *, '(a)' ) ' ' write ( *, '(a)' ) 'FILE_ROW_COUNT - Fatal error!' write ( *, '(a,i8)' ) ' Could not open the input file "' // & trim ( input_file_name ) // '" on unit ', input_unit stop end if comment_num = 0 row_num = 0 record_num = 0 bad_num = 0 do read ( input_unit, '(a)', iostat = input_status ) line if ( input_status /= 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 genrcm ( node_num, adj_num, adj_row, adj, perm ) !*****************************************************************************80 ! !! GENRCM finds the reverse Cuthill-Mckee ordering for a general graph. ! ! Discussion: ! ! For each connected component in the graph, the routine obtains ! an ordering by calling RCM. ! ! Licensing: ! ! This code is distributed under the MIT license. ! ! Modified: ! ! 04 January 2003 ! ! Author: ! ! Original FORTRAN77 version by Alan George, Joseph Liu. ! FORTRAN90 version by John Burkardt ! ! Reference: ! ! Alan George, Joseph Liu, ! Computer Solution of Large Sparse Positive Definite Systems, ! Prentice Hall, 1981. ! ! Parameters: ! ! Input, integer NODE_NUM, the number of nodes. ! ! Input, integer ADJ_NUM, the number of adjacency entries. ! ! Input, integer ADJ_ROW(NODE_NUM+1). Information about ! row I is stored in entries ADJ_ROW(I) through ADJ_ROW(I+1)-1 of ADJ. ! ! Input, integer ADJ(ADJ_NUM), the adjacency structure. ! For each row, it contains the column indices of the nonzero entries. ! ! Output, integer PERM(NODE_NUM), the RCM ordering. ! ! Local Parameters: ! ! Local, integer LEVEL_ROW(NODE_NUM+1), the index vector for a level ! structure. The level structure is stored in the currently unused ! spaces in the permutation vector PERM. ! ! Local, integer MASK(NODE_NUM), marks variables that have been numbered. ! implicit none integer adj_num integer node_num integer adj(adj_num) integer adj_row(node_num+1) integer i integer iccsze integer mask(node_num) integer level_num integer level_row(node_num+1) integer num integer perm(node_num) integer root mask(1:node_num) = 1 num = 1 do i = 1, node_num ! ! For each masked connected component... ! if ( mask(i) /= 0 ) then root = i ! ! Find a pseudo-peripheral node ROOT. The level structure found by ! ROOT_FIND is stored starting at PERM(NUM). ! call root_find ( root, adj_num, adj_row, adj, mask, level_num, & level_row, perm(num), node_num ) write ( *, * ) ' ROOT = ',root ! ! RCM orders the component using ROOT as the starting node. ! call rcm ( root, adj_num, adj_row, adj, mask, perm(num), iccsze, & node_num ) num = num + iccsze ! ! We can stop once every node is in one of the connected components. ! if ( node_num < num ) then return end if end if end do 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 MIT license. ! ! Modified: ! ! 18 September 2005 ! ! Author: ! ! John Burkardt ! ! Parameters: ! ! Output, integer IUNIT, the free unit number. ! implicit none integer i integer ios integer 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 i4col_compare ( m, n, a, i, j, isgn ) !*****************************************************************************80 ! !! I4COL_COMPARE compares columns I and J of an I4COL. ! ! Example: ! ! Input: ! ! M = 3, N = 4, I = 2, J = 4 ! ! A = ( ! 1 2 3 4 ! 5 6 7 8 ! 9 10 11 12 ) ! ! Output: ! ! ISGN = -1 ! ! Licensing: ! ! This code is distributed under the MIT license. ! ! Modified: ! ! 30 June 2000 ! ! Author: ! ! John Burkardt ! ! Parameters: ! ! Input, integer M, N, the number of rows and columns. ! ! Input, integer A(M,N), an array of N columns of vectors ! of length M. ! ! Input, integer I, J, the columns to be compared. ! I and J must be between 1 and N. ! ! Output, integer ISGN, the results of the comparison: ! -1, column I < column J, ! 0, column I = column J, ! +1, column J < column I. ! implicit none integer m integer n integer a(m,n) integer i integer isgn integer j integer k ! ! Check. ! if ( i < 1 .or. n < i ) then write ( *, '(a)' ) ' ' write ( *, '(a)' ) 'I4COL_COMPARE - Fatal error!' write ( *, '(a)' ) ' Column index I is out of bounds.' stop end if if ( j < 1 .or. n < j ) then write ( *, '(a)' ) ' ' write ( *, '(a)' ) 'I4COL_COMPARE - Fatal error!' write ( *, '(a)' ) ' Column index J is out of bounds.' stop end if isgn = 0 if ( i == j ) then return end if k = 1 do while ( k <= m ) if ( a(k,i) < a(k,j) ) then isgn = -1 return else if ( a(k,j) < a(k,i) ) then isgn = +1 return end if k = k + 1 end do return end subroutine i4col_sort_a ( m, n, a ) !*****************************************************************************80 ! !! I4COL_SORT_A ascending sorts an I4COL. ! ! Discussion: ! ! In lexicographic order, the statement "X < Y", applied to two real ! vectors X and Y of length M, means that there is some index I, with ! 1 <= I <= M, with the property that ! ! X(J) = Y(J) for J < I, ! and ! X(I) < Y(I). ! ! In other words, the first time they differ, X is smaller. ! ! Licensing: ! ! This code is distributed under the MIT license. ! ! Modified: ! ! 25 September 2001 ! ! Author: ! ! John Burkardt ! ! Parameters: ! ! Input, integer M, the number of rows of A, and the length of ! a vector of data. ! ! Input, integer N, the number of columns of A. ! ! Input/output, integer A(M,N). ! On input, the array of N columns of M-vectors. ! On output, the columns of A have been sorted in ascending ! lexicographic order. ! implicit none integer m integer n integer a(m,n) integer i integer indx integer isgn integer j if ( m <= 0 ) then return end if if ( n <= 1 ) then return end if ! ! Initialize. ! i = 0 indx = 0 isgn = 0 j = 0 ! ! Call the external heap sorter. ! do call sort_heap_external ( n, indx, i, j, isgn ) ! ! Interchange the I and J objects. ! if ( 0 < indx ) then call i4col_swap ( m, n, a, i, j ) ! ! Compare the I and J objects. ! else if ( indx < 0 ) then call i4col_compare ( m, n, a, i, j, isgn ) else if ( indx == 0 ) then exit end if end do return end subroutine i4col_swap ( m, n, a, i, j ) !*****************************************************************************80 ! !! I4COL_SWAP swaps columns I and J of an I4COL. ! ! Example: ! ! Input: ! ! M = 3, N = 4, I = 2, J = 4 ! ! A = ( ! 1 2 3 4 ! 5 6 7 8 ! 9 10 11 12 ) ! ! Output: ! ! A = ( ! 1 4 3 2 ! 5 8 7 6 ! 9 12 11 10 ) ! ! Licensing: ! ! This code is distributed under the MIT license. ! ! Modified: ! ! 04 April 2001 ! ! Author: ! ! John Burkardt ! ! Parameters: ! ! Input, integer M, N, the number of rows and columns ! in the array. ! ! Input/output, integer A(M,N), an array of N columns ! of length M. ! ! Input, integer I, J, the columns to be swapped. ! implicit none integer m integer n integer a(m,n) integer col(m) integer i integer j if ( i < 1 .or. n < i .or. j < 1 .or. n < j ) then write ( *, '(a)' ) ' ' write ( *, '(a)' ) 'I4COL_SWAP - Fatal error!' write ( *, '(a)' ) ' I or J is out of bounds.' write ( *, '(a,i8)' ) ' I = ', i write ( *, '(a,i8)' ) ' J = ', j write ( *, '(a,i8)' ) ' N = ', n stop end if if ( i == j ) then return end if col(1:m) = a(1:m,i) a(1:m,i) = a(1:m,j) a(1:m,j) = col(1:m) return end subroutine i4mat_data_read ( input_filename, m, n, table ) !*****************************************************************************80 ! !! I4MAT_DATA_READ reads data from an I4MAT file. ! ! Discussion: ! ! The file may contain more than N points, but this routine ! will return after reading N points. ! ! Licensing: ! ! This code is distributed under the MIT license. ! ! Modified: ! ! 27 January 2005 ! ! Author: ! ! John Burkardt ! ! Parameters: ! ! Input, character ( len = * ) INPUT_FILENAME, the name of the input file. ! ! Input, integer M, the spatial dimension. ! ! Input, integer N, the number of points. ! ! Output, integer TABLE(M,N), the table data. ! implicit none integer m integer n integer ierror character ( len = * ) input_filename integer input_status integer input_unit integer j character ( len = 255 ) line integer table(m,n) integer 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)' ) 'I4MAT_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 ierror = 2 write ( *, '(a)' ) ' ' write ( *, '(a)' ) 'I4MAT_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_i4vec ( 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 i4mat_header_read ( input_filename, m, n ) !*****************************************************************************80 ! !! I4MAT_HEADER_READ reads the header from an I4MAT. ! ! Licensing: ! ! This code is distributed under the MIT license. ! ! Modified: ! ! 04 June 2004 ! ! Author: ! ! John Burkardt ! ! Parameters: ! ! Input, character ( len = * ) INPUT_FILENAME, the name of the input file. ! ! Output, integer M, spatial dimension. ! ! Output, integer N, the number of points. ! implicit none character ( len = * ) input_filename integer m integer n call file_column_count ( input_filename, m ) if ( m <= 0 ) then write ( *, '(a)' ) ' ' write ( *, '(a)' ) 'I4MAT_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)' ) 'I4MAT_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 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. ! ! Discussion: ! ! An I4MAT is a rectangular array of I4 values. ! ! Licensing: ! ! This code is distributed under the MIT license. ! ! Modified: ! ! 09 February 2005 ! ! Author: ! ! John Burkardt ! ! Parameters: ! ! Input, integer M, N, the number of rows and columns. ! ! Input, integer A(M,N), an M by N matrix to be printed. ! ! Input, integer ILO, JLO, the first row and column to print. ! ! Input, integer IHI, JHI, the last row and column to print. ! ! Input, character ( len = * ) TITLE, a title. ! implicit none integer, parameter :: incx = 10 integer m integer n integer a(m,n) character ( len = 8 ) ctemp(incx) integer i integer i2 integer i2hi integer i2lo integer ihi integer ilo integer inc integer j integer j2hi integer j2lo integer jhi integer 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), '(i8)' ) i end do write ( *, '('' Row '',10a8)' ) 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), '(i8)' ) a(i,j) end do write ( *, '(i5,1x,10a8)' ) j, ( ctemp(i), i = 1, inc ) end do end do write ( *, '(a)' ) ' ' return end subroutine i4mat_write ( output_filename, m, n, table ) !*****************************************************************************80 ! !! I4MAT_WRITE writes an I4MAT file. ! ! Licensing: ! ! This code is distributed under the MIT license. ! ! Modified: ! ! 31 August 2009 ! ! Author: ! ! John Burkardt ! ! Parameters: ! ! Input, character ( len = * ) OUTPUT_FILENAME, the output file name. ! ! Input, integer M, the spatial dimension. ! ! Input, integer N, the number of points. ! ! Input, integer TABLE(M,N), the table data. ! implicit none integer m integer n integer j character ( len = * ) output_filename integer output_status integer output_unit character ( len = 30 ) string integer table(m,n) ! ! Open the file. ! call get_unit ( output_unit ) open ( unit = output_unit, file = output_filename, & status = 'replace', iostat = output_status ) if ( output_status /= 0 ) then write ( *, '(a)' ) ' ' write ( *, '(a)' ) 'I4MAT_WRITE - Fatal error!' write ( *, '(a,i8)' ) ' Could not open the output file "' // & trim ( output_filename ) // '" on unit ', output_unit output_unit = -1 stop end if ! ! Create a format string. ! if ( 0 < m .and. 0 < n ) then write ( string, '(a1,i8,a4)' ) '(', m, 'i10)' ! ! Write the data. ! do j = 1, n write ( output_unit, string ) table(1:m,j) end do end if ! ! Close the file. ! close ( unit = output_unit ) return end subroutine i4vec_heap_d ( n, a ) !*****************************************************************************80 ! !! I4VEC_HEAP_D reorders an I4VEC into an descending heap. ! ! Discussion: ! ! An I4VEC is a vector of integer values. ! ! A descending heap is an array A with the property that, for every index J, ! A(J) >= A(2*J) and A(J) >= A(2*J+1), (as long as the indices ! 2*J and 2*J+1 are legal). ! ! A(1) ! / \ ! A(2) A(3) ! / \ / \ ! A(4) A(5) A(6) A(7) ! / \ / \ ! A(8) A(9) A(10) A(11) ! ! Licensing: ! ! This code is distributed under the MIT license. ! ! Modified: ! ! 15 April 1999 ! ! Author: ! ! John Burkardt ! ! Reference: ! ! Albert Nijenhuis, Herbert Wilf, ! Combinatorial Algorithms, ! Academic Press, 1978, second edition, ! ISBN 0-12-519260-6. ! ! Parameters: ! ! Input, integer N, the size of the input array. ! ! Input/output, integer A(N). ! On input, an unsorted array. ! On output, the array has been reordered into a heap. ! implicit none integer n integer a(n) integer i integer ifree integer key integer m ! ! Only nodes N/2 down to 1 can be "parent" nodes. ! do i = n/2, 1, -1 ! ! Copy the value out of the parent node. ! Position IFREE is now "open". ! key = a(i) ifree = i do ! ! Positions 2*IFREE and 2*IFREE + 1 are the descendants of position ! IFREE. (One or both may not exist because they exceed N.) ! m = 2 * ifree ! ! Does the first position exist? ! if ( n < m ) then exit end if ! ! Does the second position exist? ! if ( m + 1 <= n ) then ! ! If both positions exist, take the larger of the two values, ! and update M if necessary. ! if ( a(m) < a(m+1) ) then m = m + 1 end if end if ! ! If the large descendant is larger than KEY, move it up, ! and update IFREE, the location of the free position, and ! consider the descendants of THIS position. ! if ( a(m) <= key ) then exit end if a(ifree) = a(m) ifree = m end do ! ! Once there is no more shifting to do, KEY moves into the free spot IFREE. ! a(ifree) = key end do return end subroutine i4vec_print ( n, a, title ) !*****************************************************************************80 ! !! I4VEC_PRINT prints an I4VEC. ! ! Discussion: ! ! An I4VEC is a vector of I4's. ! ! Licensing: ! ! This code is distributed under the MIT license. ! ! Modified: ! ! 28 November 2000 ! ! Author: ! ! John Burkardt ! ! Parameters: ! ! Input, integer N, the number of components of the vector. ! ! Input, integer A(N), the vector to be printed. ! ! Input, character ( len = * ) TITLE, a title. ! implicit none integer n integer a(n) integer i character ( len = * ) title if ( 0 < len_trim ( title ) ) then write ( *, '(a)' ) ' ' write ( *, '(a)' ) trim ( title ) end if write ( *, '(a)' ) ' ' do i = 1, n write ( *, '(2x,i8,2x,i12)' ) i, a(i) end do return end subroutine i4vec_reverse ( n, a ) !*****************************************************************************80 ! !! I4VEC_REVERSE reverses the elements of an I4VEC. ! ! Example: ! ! Input: ! ! N = 5, ! A = ( 11, 12, 13, 14, 15 ). ! ! Output: ! ! A = ( 15, 14, 13, 12, 11 ). ! ! Licensing: ! ! This code is distributed under the MIT license. ! ! Modified: ! ! 30 September 2009 ! ! Author: ! ! John Burkardt ! ! Parameters: ! ! Input, integer N, the number of entries in the array. ! ! Input/output, integer A(N), the array to be reversed. ! implicit none integer n integer a(n) a(1:n) = a(n:1:-1) return end subroutine i4vec_sort_heap_a ( n, a ) !*****************************************************************************80 ! !! I4VEC_SORT_HEAP_A ascending sorts an I4VEC using heap sort. ! ! Discussion: ! ! An I4VEC is a vector of integer values. ! ! Licensing: ! ! This code is distributed under the MIT license. ! ! Modified: ! ! 30 September 2009 ! ! Author: ! ! John Burkardt ! ! Reference: ! ! Albert Nijenhuis, Herbert Wilf, ! Combinatorial Algorithms, ! Academic Press, 1978, second edition, ! ISBN 0-12-519260-6. ! ! Parameters: ! ! Input, integer N, the number of entries in the array. ! ! Input/output, integer A(N). ! On input, the array to be sorted; ! On output, the array has been sorted. ! implicit none integer n integer a(n) integer n1 integer t if ( n <= 1 ) then return end if ! ! 1: Put A into descending heap form. ! call i4vec_heap_d ( n, a ) ! ! 2: Sort A. ! ! The largest object in the heap is in A(1). ! Move it to position A(N). ! t = a(1) a(1) = a(n) a(n) = t ! ! Consider the diminished heap of size N1. ! do n1 = n - 1, 2, -1 ! ! Restore the heap structure of A(1) through A(N1). ! call i4vec_heap_d ( n1, a ) ! ! Take the largest object from A(1) and move it to A(N1). ! t = a(1) a(1) = a(n1) a(n1) = t end do return end subroutine level_set ( root, adj_num, adj_row, adj, mask, level_num, & level_row, level, node_num ) !*****************************************************************************80 ! !! LEVEL_SET generates the connected level structure rooted at a given node. ! ! Discussion: ! ! Only nodes for which MASK is nonzero will be considered. ! ! The root node chosen by the user is assigned level 1, and masked. ! All (unmasked) nodes reachable from a node in level 1 are ! assigned level 2 and masked. The process continues until there ! are no unmasked nodes adjacent to any node in the current level. ! The number of levels may vary between 2 and NODE_NUM. ! ! Licensing: ! ! This code is distributed under the MIT license. ! ! Modified: ! ! 28 October 2003 ! ! Author: ! ! Original FORTRAN77 version by Alan George, Joseph Liu. ! FORTRAN90 version by John Burkardt ! ! Reference: ! ! Alan George, Joseph Liu, ! Computer Solution of Large Sparse Positive Definite Systems, ! Prentice Hall, 1981. ! ! Parameters: ! ! Input, integer ROOT, the node at which the level structure ! is to be rooted. ! ! Input, integer ADJ_NUM, the number of adjacency entries. ! ! Input, integer ADJ_ROW(NODE_NUM+1). Information about ! row I is stored in entries ADJ_ROW(I) through ADJ_ROW(I+1)-1 of ADJ. ! ! Input, integer ADJ(ADJ_NUM), the adjacency structure. ! For each row, it contains the column indices of the nonzero entries. ! ! Input/output, integer MASK(NODE_NUM). On input, only nodes ! with nonzero MASK are to be processed. On output, those nodes which were ! included in the level set have MASK set to 1. ! ! Output, integer LEVEL_NUM, the number of levels in the level ! structure. ROOT is in level 1. The neighbors of ROOT ! are in level 2, and so on. ! ! Output, integer LEVEL_ROW(NODE_NUM+1), LEVEL(NODE_NUM), ! the rooted level structure. ! ! Input, integer NODE_NUM, the number of nodes. ! implicit none integer adj_num integer node_num integer adj(adj_num) integer adj_row(node_num+1) integer i integer iccsze integer j integer jstop integer jstrt integer lbegin integer level_num integer level_row(node_num+1) integer level(node_num) integer lvlend integer lvsize integer mask(node_num) integer nbr integer node integer root mask(root) = 0 level(1) = root level_num = 0 lvlend = 0 iccsze = 1 ! ! LBEGIN is the pointer to the beginning of the current level, and ! LVLEND points to the end of this level. ! do lbegin = lvlend + 1 lvlend = iccsze level_num = level_num + 1 level_row(level_num) = lbegin ! ! Generate the next level by finding all the masked neighbors of nodes ! in the current level. ! do i = lbegin, lvlend node = level(i) jstrt = adj_row(node) jstop = adj_row(node+1) - 1 do j = jstrt, jstop nbr = adj(j) if ( mask(nbr) /= 0 ) then iccsze = iccsze + 1 level(iccsze) = nbr mask(nbr) = 0 end if end do end do ! ! Compute the current level width (the number of nodes encountered.) ! If it is positive, generate the next level. ! lvsize = iccsze - lvlend if ( lvsize <= 0 ) then exit end if end do level_row(level_num+1) = lvlend + 1 ! ! Reset MASK to 1 for the nodes in the level structure. ! mask(level(1:iccsze)) = 1 return end subroutine neighbor_elements_q4_mesh ( element_num, element_node, & element_neighbor ) !*****************************************************************************80 ! !! NEIGHBOR_ELEMENTS_Q4_MESH determines element neighbors in a Q4 mesh. ! ! Discussion: ! ! A quadrilateral mesh of a set of nodes can be completely described by ! the coordinates of the nodes, and the list of nodes that make up ! each element. However, in some cases, it is necessary to know ! element adjacency information, that is, which element, if any, ! is adjacent to a given element on a particular side. ! ! This routine creates a data structure recording this information. ! ! The primary amount of work occurs in sorting a list of 4 * ELEMENT_NUM ! data items. ! ! Note that COL is a work array allocated dynamically inside this ! routine. It is possible, for very large values of ELEMENT_NUM, ! that the necessary amount of memory will not be accessible, and the ! routine will fail. This is a limitation of the implementation of ! dynamic arrays in FORTRAN90. One way to get around this would be ! to require the user to declare ROW in the calling routine ! as an allocatable array, get the necessary memory explicitly with ! an ALLOCATE statement, and then pass ROW into this routine. ! ! Of course, the point of dynamic arrays was to make it easy to ! hide these sorts of temporary work arrays from the poor user! ! ! This routine was revised to store the edge data in a column ! array rather than a row array. ! ! Licensing: ! ! This code is distributed under the MIT license. ! ! Modified: ! ! 14 February 2006 ! ! Author: ! ! John Burkardt ! ! Parameters: ! ! Input, integer ELEMENT_NUM, the number of elements. ! ! Input, integer ELEMENT_NODE(4,ELEMENT_NUM), the nodes ! that make up each element. ! ! Output, integer ELEMENT_NEIGHBOR(4,ELEMENT_NUM), lists the ! neighboring element on each side of a given element, or -1 if there is ! no neighbor. ! implicit none integer element_num integer, parameter :: element_order = 4 integer, allocatable :: col(:,:) integer element_neighbor(4,element_num) integer element_node(element_order,element_num) integer i integer icol integer j integer k integer l integer q integer q1 integer q2 integer side1 integer side2 allocate ( col (4,4*element_num) ) ! ! Step 1. ! From the list of nodes for element Q, of the form: (I,J,K,L) ! construct the four neighbor relations: ! ! (I,J,1,Q) or (J,I,1,Q), ! (J,K,2,Q) or (K,J,2,Q), ! (K,L,3,Q) or (L,K,3,Q) ! (K,I,4,Q) or (I,K,4,Q) ! ! where we choose (I,J,1,Q) if I < J, or else (J,I,1,Q) ! do q = 1, element_num i = element_node(1,q) j = element_node(2,q) k = element_node(3,q) l = element_node(4,q) if ( i < j ) then col(1:4,4*(q-1)+1) = (/ i, j, 1, q /) else col(1:4,4*(q-1)+1) = (/ j, i, 1, q /) end if if ( j < k ) then col(1:4,4*(q-1)+2) = (/ j, k, 2, q /) else col(1:4,4*(q-1)+2) = (/ k, j, 2, q /) end if if ( k < l ) then col(1:4,4*(q-1)+3) = (/ k, l, 3, q /) else col(1:4,4*(q-1)+3) = (/ l, k, 3, q /) end if if ( l < i ) then col(1:4,4*(q-1)+4) = (/ l, i, 4, q /) else col(1:4,4*(q-1)+4) = (/ i, l, 4, q /) end if end do ! ! Step 2. Perform an ascending dictionary sort on the neighbor relations. ! We only intend to sort on rows 1 and 2; the routine we call here ! sorts on rows 1 through 4 but that won't hurt us. ! ! What we need is to find cases where two elements share an edge. ! Say they share an edge defined by the nodes I and J. Then there are ! two columns of COL that start out ( I, J, ?, ? ). By sorting COL, ! we make sure that these two columns occur consecutively. That will ! make it easy to notice that the elements are neighbors. ! call i4col_sort_a ( 4, 4*element_num, col ) ! ! Step 3. Neighboring elements show up as consecutive columns with ! identical first two entries. Whenever you spot this happening, ! make the appropriate entries in ELEMENT_NEIGHBOR. ! element_neighbor(1:4,1:element_num) = -1 icol = 1 do if ( 4 * element_num <= icol ) then exit end if if ( col(1,icol) /= col(1,icol+1) .or. col(2,icol) /= col(2,icol+1) ) then icol = icol + 1 cycle end if side1 = col(3,icol) q1 = col(4,icol) side2 = col(3,icol+1) q2 = col(4,icol+1) element_neighbor(side1,q1) = q2 element_neighbor(side2,q2) = q1 icol = icol + 2 end do deallocate ( col ) return end subroutine perm_check ( n, p, base, ierror ) !*****************************************************************************80 ! !! PERM_CHECK checks that a vector represents a permutation. ! ! Discussion: ! ! The routine verifies that each of the integers from BASE to ! to BASE+N-1 occurs among the N entries of the permutation. ! ! Set the input quantity BASE to 0, if P is a 0-based permutation, ! or to 1 if P is a 1-based permutation. ! ! Licensing: ! ! This code is distributed under the MIT license. ! ! Modified: ! ! 31 October 2008 ! ! Author: ! ! John Burkardt ! ! Parameters: ! ! Input, integer N, the number of entries. ! ! Input, integer P(N), the array to check. ! ! Input, integer BASE, the index base. ! ! Output, integer IERROR, error flag. ! 0, the array represents a permutation. ! nonzero, the array does not represent a permutation. The smallest ! missing value is equal to IERROR. ! implicit none integer n integer base integer find integer ierror integer p(n) integer seek ierror = 0 do seek = base, base + n - 1 ierror = 1 do find = 1, n if ( p(find) == seek ) then ierror = 0 exit end if end do if ( ierror /= 0 ) then write ( *, '(a)' ) ' ' write ( *, '(a)' ) 'PERM_CHECK - Fatal error!' write ( *, '(a)' ) ' The input array does not represent' write ( *, '(a)' ) ' a proper permutation.' stop end if end do return end subroutine perm_inverse3 ( n, perm, perm_inv ) !*****************************************************************************80 ! !! PERM_INVERSE3 produces the inverse of a given permutation. ! ! Licensing: ! ! This code is distributed under the MIT license. ! ! Modified: ! ! 28 October 2003 ! ! Author: ! ! John Burkardt ! ! Parameters: ! ! Input, integer N, the number of items permuted. ! ! Input, integer PERM(N), a permutation. ! ! Output, integer PERM_INV(N), the inverse permutation. ! implicit none integer n integer i integer perm(n) integer perm_inv(n) do i = 1, n perm_inv(perm(i)) = i end do return end subroutine r8col_permute ( m, n, p, a ) !*****************************************************************************80 ! !! R8COL_PERMUTE permutes an R8COL in place. ! ! Discussion: ! ! An R8COL is an M by N array of R8's, regarded as an array of N columns, ! each of length M. ! ! 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: ! ! M = 2 ! 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 MIT license. ! ! Modified: ! ! 09 December 2006 ! ! Author: ! ! John Burkardt ! ! Parameters: ! ! Input, integer M, the dimension of objects. ! ! Input, integer N, the number of objects. ! ! Input, integer P(N), the permutation. P(I) = J means ! that the I-th element of the output array should be the J-th ! element of the input array. ! ! Input/output, real ( kind = rk ) A(M,N), the array to be permuted. ! implicit none integer, parameter :: rk = kind ( 1.0D+00 ) integer m integer n real ( kind = rk ) a(m,n) real ( kind = rk ) a_temp(m) integer, parameter :: base = 1 integer ierror integer iget integer iput integer istart integer p(n) call perm_check ( n, p, base, ierror ) if ( ierror /= 0 ) then write ( *, '(a)' ) ' ' write ( *, '(a)' ) 'R8COL_PERMUTE - Fatal error!' write ( *, '(a)' ) ' PERM_CHECK rejects this permutation.' stop end if ! ! 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:m) = a(1:m,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)' ) 'R8COL_PERMUTE - Fatal error!' write ( *, '(a)' ) ' A permutation index is out of range.' write ( *, '(a,i8,a,i8)' ) ' P(', iput, ') = ', iget stop end if if ( iget == istart ) then a(1:m,iput) = a_temp(1:m) exit end if a(1:m,iput) = a(1:m,iget) end do end if end do ! ! Restore the signs of the entries. ! p(1:n) = - p(1:n) return end subroutine r8mat_data_read ( input_filename, m, n, table ) !*****************************************************************************80 ! !! R8MAT_DATA_READ reads data from an R8MAT file. ! ! 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 MIT license. ! ! Modified: ! ! 18 October 2008 ! ! Author: ! ! John Burkardt ! ! Parameters: ! ! Input, character ( len = * ) INPUT_FILENAME, the name of the input file. ! ! Input, integer M, the spatial dimension. ! ! Input, integer N, the number of points. ! ! Output, real ( kind = rk ) TABLE(M,N), the table data. ! implicit none integer, parameter :: rk = kind ( 1.0D+00 ) integer m integer n integer ierror character ( len = * ) input_filename integer input_status integer input_unit integer j character ( len = 255 ) line real ( kind = rk ) table(m,n) real ( kind = rk ) 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. ! ! Licensing: ! ! This code is distributed under the MIT license. ! ! Modified: ! ! 07 September 2004 ! ! Author: ! ! John Burkardt ! ! Parameters: ! ! Input, character ( len = * ) INPUT_FILENAME, the name of the input file. ! ! Output, integer M, spatial dimension. ! ! Output, integer N, the number of points. ! implicit none character ( len = * ) input_filename integer m integer 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_some ( m, n, a, ilo, jlo, ihi, jhi, title ) !*****************************************************************************80 ! !! R8MAT_TRANSPOSE_PRINT_SOME prints some of an R8MAT, transposed. ! ! Discussion: ! ! An R8MAT is an array of R8 values. ! ! Licensing: ! ! This code is distributed under the MIT license. ! ! Modified: ! ! 14 June 2004 ! ! Author: ! ! John Burkardt ! ! Parameters: ! ! Input, integer M, N, the number of rows and columns. ! ! Input, real ( kind = rk ) A(M,N), an M by N matrix to be printed. ! ! Input, integer ILO, JLO, the first row and column to print. ! ! Input, integer IHI, JHI, the last row and column to print. ! ! Input, character ( len = * ) TITLE, a title. ! implicit none integer, parameter :: rk = kind ( 1.0D+00 ) integer, parameter :: incx = 5 integer m integer n real ( kind = rk ) a(m,n) character ( len = 14 ) ctemp(incx) integer i integer i2 integer i2hi integer i2lo integer ihi integer ilo integer inc integer j integer j2hi integer j2lo integer jhi integer 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), '(i8,6x)' ) 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 write ( *, '(a)' ) ' ' return end subroutine r8mat_write ( output_filename, m, n, table ) !*****************************************************************************80 ! !! R8MAT_WRITE writes an R8MAT file. ! ! Licensing: ! ! This code is distributed under the MIT license. ! ! Modified: ! ! 31 May 2009 ! ! Author: ! ! John Burkardt ! ! Parameters: ! ! Input, character ( len = * ) OUTPUT_FILENAME, the output file name. ! ! Input, integer M, the spatial dimension. ! ! Input, integer N, the number of points. ! ! Input, real ( kind = rk ) TABLE(M,N), the table data. ! implicit none integer, parameter :: rk = kind ( 1.0D+00 ) integer m integer n integer j character ( len = * ) output_filename integer output_status integer output_unit character ( len = 30 ) string real ( kind = rk ) table(m,n) ! ! Open the file. ! call get_unit ( output_unit ) open ( unit = output_unit, file = output_filename, & status = 'replace', iostat = output_status ) if ( output_status /= 0 ) then write ( *, '(a)' ) ' ' write ( *, '(a)' ) 'R8MAT_WRITE - Fatal error!' write ( *, '(a,i8)' ) ' Could not open the output file "' // & trim ( output_filename ) // '" on unit ', output_unit output_unit = -1 stop end if ! ! Create a format string. ! ! For greater precision in the output file, try: ! ! '(', m, 'g', 24, '.', 16, ')' ! if ( 0 < m .and. 0 < n ) then write ( string, '(a1,i8,a1,i8,a1,i8,a1)' ) '(', m, 'g', 14, '.', 6, ')' ! ! Write the data. ! do j = 1, n write ( output_unit, string ) table(1:m,j) end do end if ! ! Close the file. ! close ( unit = output_unit ) return end subroutine rcm ( root, adj_num, adj_row, adj, mask, perm, iccsze, node_num ) !*****************************************************************************80 ! !! RCM renumbers a connected component by the reverse Cuthill McKee algorithm. ! ! Discussion: ! ! The connected component is specified by a node ROOT and a mask. ! The numbering starts at the root node. ! ! An outline of the algorithm is as follows: ! ! X(1) = ROOT. ! ! for ( I = 1 to N-1) ! Find all unlabeled neighbors of X(I), ! assign them the next available labels, in order of increasing degree. ! ! When done, reverse the ordering. ! ! Licensing: ! ! This code is distributed under the MIT license. ! ! Modified: ! ! 02 January 2007 ! ! Author: ! ! Original FORTRAN77 version by Alan George, Joseph Liu. ! FORTRAN90 version by John Burkardt ! ! Reference: ! ! Alan George, Joseph Liu, ! Computer Solution of Large Sparse Positive Definite Systems, ! Prentice Hall, 1981. ! ! Parameters: ! ! Input, integer ROOT, the node that defines the connected ! component. It is used as the starting point for the RCM ordering. ! ! Input, integer ADJ_NUM, the number of adjacency entries. ! ! Input, integer ADJ_ROW(NODE_NUM+1). Information about ! row I is stored in entries ADJ_ROW(I) through ADJ_ROW(I+1)-1 of ADJ. ! ! Input, integer ADJ(ADJ_NUM), the adjacency structure. ! For each row, it contains the column indices of the nonzero entries. ! ! Input/output, integer MASK(NODE_NUM), a mask for the nodes. ! Only those nodes with nonzero input mask values are considered by the ! routine. The nodes numbered by RCM will have their mask values ! set to zero. ! ! Output, integer PERM(NODE_NUM), the RCM ordering. ! ! Output, integer ICCSZE, the size of the connected component ! that has been numbered. ! ! Input, integer NODE_NUM, the number of nodes. ! ! Local Parameters: ! ! Workspace, integer DEG(NODE_NUM), a temporary vector used to hold ! the degree of the nodes in the section graph specified by mask and root. ! implicit none integer adj_num integer node_num integer adj(adj_num) integer adj_row(node_num+1) integer deg(node_num) integer fnbr integer i integer iccsze integer j integer jstop integer jstrt integer k integer l integer lbegin integer lnbr integer lperm integer lvlend integer mask(node_num) integer nbr integer node integer perm(node_num) integer root ! ! Find the degrees of the nodes in the component specified by MASK and ROOT. ! call degree ( root, adj_num, adj_row, adj, mask, deg, iccsze, perm, node_num ) mask(root) = 0 if ( iccsze <= 1 ) then return end if lvlend = 0 lnbr = 1 ! ! LBEGIN and LVLEND point to the beginning and ! the end of the current level respectively. ! do while ( lvlend < lnbr ) lbegin = lvlend + 1 lvlend = lnbr do i = lbegin, lvlend ! ! For each node in the current level... ! node = perm(i) jstrt = adj_row(node) jstop = adj_row(node+1) - 1 ! ! Find the unnumbered neighbors of NODE. ! ! FNBR and LNBR point to the first and last neighbors ! of the current node in PERM. ! fnbr = lnbr + 1 do j = jstrt, jstop nbr = adj(j) if ( mask(nbr) /= 0 ) then lnbr = lnbr + 1 mask(nbr) = 0 perm(lnbr) = nbr end if end do ! ! If no neighbors, skip to next node in this level. ! if ( lnbr <= fnbr ) then cycle end if ! ! Sort the neighbors of NODE in increasing order by degree. ! Linear insertion is used. ! k = fnbr do while ( k < lnbr ) l = k k = k + 1 nbr = perm(k) do while ( fnbr < l ) lperm = perm(l) if ( deg(lperm) <= deg(nbr) ) then exit end if perm(l+1) = lperm l = l - 1 end do perm(l+1) = nbr end do end do end do ! ! We now have the Cuthill-McKee ordering. Reverse it. ! call i4vec_reverse ( iccsze, perm ) return end subroutine root_find ( root, adj_num, adj_row, adj, mask, level_num, & level_row, level, node_num ) !*****************************************************************************80 ! !! ROOT_FIND finds a pseudo-peripheral node. ! ! Discussion: ! ! The diameter of a graph is the maximum distance (number of edges) ! between any two nodes of the graph. ! ! The eccentricity of a node is the maximum distance between that ! node and any other node of the graph. ! ! A peripheral node is a node whose eccentricity equals the ! diameter of the graph. ! ! A pseudo-peripheral node is an approximation to a peripheral node; ! it may be a peripheral node, but all we know is that we tried our ! best. ! ! The routine is given a graph, and seeks pseudo-peripheral nodes, ! using a modified version of the scheme of Gibbs, Poole and ! Stockmeyer. It determines such a node for the section subgraph ! specified by MASK and ROOT. ! ! The routine also determines the level structure associated with ! the given pseudo-peripheral node; that is, how far each node ! is from the pseudo-peripheral node. The level structure is ! returned as a list of nodes LS, and pointers to the beginning ! of the list of nodes that are at a distance of 0, 1, 2, ..., ! NODE_NUM-1 from the pseudo-peripheral node. ! ! Licensing: ! ! This code is distributed under the MIT license. ! ! Modified: ! ! 28 October 2003 ! ! Author: ! ! Original FORTRAN77 version by Alan George, Joseph Liu. ! FORTRAN90 version by John Burkardt ! ! Reference: ! ! Alan George, Joseph Liu, ! Computer Solution of Large Sparse Positive Definite Systems, ! Prentice Hall, 1981. ! ! Norman Gibbs, William Poole, Paul Stockmeyer, ! An Algorithm for Reducing the Bandwidth and Profile of a Sparse Matrix, ! SIAM Journal on Numerical Analysis, ! Volume 13, pages 236-250, 1976. ! ! Norman Gibbs, ! Algorithm 509: A Hybrid Profile Reduction Algorithm, ! ACM Transactions on Mathematical Software, ! Volume 2, pages 378-387, 1976. ! ! Parameters: ! ! Input/output, integer ROOT. On input, ROOT is a node in the ! the component of the graph for which a pseudo-peripheral node is ! sought. On output, ROOT is the pseudo-peripheral node obtained. ! ! Input, integer ADJ_NUM, the number of adjacency entries. ! ! Input, integer ADJ_ROW(NODE_NUM+1). Information about ! row I is stored in entries ADJ_ROW(I) through ADJ_ROW(I+1)-1 of ADJ. ! ! Input, integer ADJ(ADJ_NUM), the adjacency structure. ! For each row, it contains the column indices of the nonzero entries. ! ! Input, integer MASK(NODE_NUM), specifies a section subgraph. ! Nodes for which MASK is zero are ignored by FNROOT. ! ! Output, integer LEVEL_NUM, is the number of levels in the ! level structure rooted at the node ROOT. ! ! Output, integer LEVEL_ROW(NODE_NUM+1), LEVEL(NODE_NUM), the ! level structure array pair containing the level structure found. ! ! Input, integer NODE_NUM, the number of nodes. ! implicit none integer adj_num integer node_num integer adj(adj_num) integer adj_row(node_num+1) integer iccsze integer j integer jstrt integer k integer kstop integer kstrt integer level(node_num) integer level_num integer level_num2 integer level_row(node_num+1) integer mask(node_num) integer mindeg integer nabor integer ndeg integer node integer root ! ! Determine the level structure rooted at ROOT. ! call level_set ( root, adj_num, adj_row, adj, mask, level_num, & level_row, level, node_num ) ! ! Count the number of nodes in this level structure. ! iccsze = level_row(level_num+1) - 1 ! ! Extreme case: ! A complete graph has a level set of only a single level. ! Every node is equally good (or bad). ! if ( level_num == 1 ) then return end if ! ! Extreme case: ! A "line graph" 0--0--0--0--0 has every node in its only level. ! By chance, we've stumbled on the ideal root. ! if ( level_num == iccsze ) then return end if ! ! Pick any node from the last level that has minimum degree ! as the starting point to generate a new level set. ! do mindeg = iccsze jstrt = level_row(level_num) root = level(jstrt) if ( jstrt < iccsze ) then do j = jstrt, iccsze node = level(j) ndeg = 0 kstrt = adj_row(node) kstop = adj_row(node+1) - 1 do k = kstrt, kstop nabor = adj(k) if ( 0 < mask(nabor) ) then ndeg = ndeg + 1 end if end do if ( ndeg < mindeg ) then root = node mindeg = ndeg end if end do end if ! ! Generate the rooted level structure associated with this node. ! call level_set ( root, adj_num, adj_row, adj, mask, level_num2, & level_row, level, node_num ) ! ! If the number of levels did not increase, accept the new ROOT. ! if ( level_num2 <= level_num ) then exit end if level_num = level_num2 ! ! In the unlikely case that ROOT is one endpoint of a line graph, ! we can exit now. ! if ( iccsze <= level_num ) then exit end if end do return end subroutine s_to_i4 ( s, ival, ierror, length ) !*****************************************************************************80 ! !! S_TO_I4 reads an I4 from a string. ! ! Licensing: ! ! This code is distributed under the MIT license. ! ! Modified: ! ! 28 June 2000 ! ! Author: ! ! John Burkardt ! ! Parameters: ! ! Input, character ( len = * ) S, a string to be examined. ! ! Output, integer IVAL, the integer value read from the string. ! If the string is blank, then IVAL will be returned 0. ! ! Output, integer IERROR, an error flag. ! 0, no error. ! 1, an error occurred. ! ! Output, integer LENGTH, the number of characters of S ! used to make IVAL. ! implicit none character c integer i integer ierror integer isgn integer istate integer ival integer length character ( len = * ) s ierror = 0 istate = 0 isgn = 1 ival = 0 do i = 1, len_trim ( s ) c = s(i:i) ! ! Haven't read anything. ! if ( istate == 0 ) then if ( c == ' ' ) then else if ( c == '-' ) then istate = 1 isgn = -1 else if ( c == '+' ) then istate = 1 isgn = + 1 else if ( lle ( '0', c ) .and. lle ( c, '9' ) ) then istate = 2 ival = ichar ( c ) - ichar ( '0' ) else ierror = 1 return end if ! ! Have read the sign, expecting digits. ! else if ( istate == 1 ) then if ( c == ' ' ) then else if ( lle ( '0', c ) .and. lle ( c, '9' ) ) then istate = 2 ival = ichar ( c ) - ichar ( '0' ) else ierror = 1 return end if ! ! Have read at least one digit, expecting more. ! else if ( istate == 2 ) then if ( lle ( '0', c ) .and. lle ( c, '9' ) ) then ival = 10 * ival + ichar ( c ) - ichar ( '0' ) else ival = isgn * ival length = i - 1 return end if end if end do ! ! If we read all the characters in the string, see if we're OK. ! if ( istate == 2 ) then ival = isgn * ival length = len_trim ( s ) else ierror = 1 length = 0 end if return end subroutine s_to_i4vec ( s, n, ivec, ierror ) !*****************************************************************************80 ! !! S_TO_I4VEC reads an I4VEC from a string. ! ! Licensing: ! ! This code is distributed under the MIT license. ! ! Modified: ! ! 08 October 2003 ! ! Author: ! ! John Burkardt ! ! Parameters: ! ! Input, character ( len = * ) S, the string to be read. ! ! Input, integer N, the number of values expected. ! ! Output, integer IVEC(N), the values read from the string. ! ! Output, integer IERROR, error flag. ! 0, no errors occurred. ! -K, could not read data for entries -K through N. ! implicit none integer n integer i integer ierror integer ilo integer ivec(n) integer length character ( len = * ) s i = 0 ierror = 0 ilo = 1 do while ( i < n ) i = i + 1 call s_to_i4 ( s(ilo:), ivec(i), ierror, length ) if ( ierror /= 0 ) then ierror = -i exit end if ilo = ilo + length end do return end subroutine s_to_r8 ( s, dval, ierror, length ) !*****************************************************************************80 ! !! S_TO_R8 reads an R8 from a string. ! ! Discussion: ! ! The 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 number. ! ! Legal input is: ! ! 1 blanks, ! 2 '+' or '-' sign, ! 2.5 blanks ! 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 DVAL ! ! '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 MIT license. ! ! Modified: ! ! 07 September 2004 ! ! 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 = rk ) DVAL, the value read from the string. ! ! Output, integer 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 ! to form the number, including any terminating ! characters such as a trailing comma or blanks. ! implicit none integer, parameter :: rk = kind ( 1.0D+00 ) character c logical ch_eqi real ( kind = rk ) dval integer ierror integer ihave integer isgn integer iterm integer jbot integer jsgn integer jtop integer length integer nchar integer ndig real ( kind = rk ) rbot real ( kind = rk ) rexp real ( kind = rk ) rtop character ( len = * ) s nchar = len_trim ( s ) ierror = 0 dval = 0.0D+00 length = -1 isgn = 1 rtop = 0 rbot = 1 jsgn = 1 jtop = 0 jbot = 1 ihave = 1 iterm = 0 do length = length + 1 if ( nchar < length+1 ) then exit end if c = s(length+1:length+1) ! ! Blank character. ! if ( c == ' ' ) 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 ! ! Scientific notation 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. lle ( '0', c ) .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 = rk ) else if ( ihave == 5 ) then rtop = 10.0D+00 * rtop + real ( ndig, kind = rk ) 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 ) 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 write ( *, '(a)' ) ' ' write ( *, '(a)' ) 'S_TO_R8 - Serious error!' write ( *, '(a)' ) ' Illegal or nonnumeric input:' write ( *, '(a)' ) ' ' // trim ( s ) 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 = 10.0D+00 ** ( real ( jsgn * jtop, kind = rk ) & / real ( jbot, kind = rk ) ) end if end if dval = real ( isgn, kind = rk ) * 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 MIT license. ! ! Modified: ! ! 07 September 2004 ! ! Author: ! ! John Burkardt ! ! Parameters: ! ! Input, character ( len = * ) S, the string to be read. ! ! Input, integer N, the number of values expected. ! ! Output, real ( kind = rk ) RVEC(N), the values read from the string. ! ! Output, integer IERROR, error flag. ! 0, no errors occurred. ! -K, could not read data for entries -K through N. ! implicit none integer, parameter :: rk = kind ( 1.0D+00 ) integer n integer i integer ierror integer ilo integer lchar real ( kind = rk ) 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, lchar ) if ( ierror /= 0 ) then ierror = -i exit end if ilo = ilo + lchar 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 MIT license. ! ! Modified: ! ! 14 April 1999 ! ! Author: ! ! John Burkardt ! ! Parameters: ! ! Input, character ( len = * ) S, the string to be examined. ! ! Output, integer NWORD, the number of "words" in the string. ! Words are presumed to be separated by one or more blanks. ! implicit none logical blank integer i integer lens integer 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 sort_heap_external ( n, indx, i, j, isgn ) !*****************************************************************************80 ! !! SORT_HEAP_EXTERNAL externally sorts a list of items into ascending order. ! ! Discussion: ! ! The actual list of data is not passed to the routine. Hence this ! routine may be used to sort integers, reals, numbers, names, ! dates, shoe sizes, and so on. After each call, the routine asks ! the user to compare or interchange two items, until a special ! return value signals that the sorting is completed. ! ! Licensing: ! ! This code is distributed under the MIT license. ! ! Modified: ! ! 05 February 2004 ! ! Author: ! ! Original FORTRAN77 version by Albert Nijenhuis, Herbert Wilf. ! FORTRAN90 version by John Burkardt. ! ! Reference: ! ! Albert Nijenhuis, Herbert Wilf, ! Combinatorial Algorithms, ! Academic Press, 1978, second edition, ! ISBN 0-12-519260-6. ! ! Parameters: ! ! Input, integer N, the number of items to be sorted. ! ! Input/output, integer INDX, the main communication signal. ! ! The user must set INDX to 0 before the first call. ! Thereafter, the user should not change the value of INDX until ! the sorting is done. ! ! On return, if INDX is ! ! greater than 0, ! * interchange items I and J; ! * call again. ! ! less than 0, ! * compare items I and J; ! * set ISGN = -1 if I < J, ISGN = +1 if J < I; ! * call again. ! ! equal to 0, the sorting is done. ! ! Output, integer I, J, the indices of two items. ! On return with INDX positive, elements I and J should be interchanged. ! On return with INDX negative, elements I and J should be compared, and ! the result reported in ISGN on the next call. ! ! Input, integer ISGN, results of comparison of elements ! I and J. (Used only when the previous call returned INDX less than 0). ! ISGN <= 0 means I is less than or equal to J; ! 0 <= ISGN means I is greater than or equal to J. ! implicit none integer i integer, save :: i_save = 0 integer indx integer isgn integer j integer, save :: j_save = 0 integer, save :: k = 0 integer, save :: k1 = 0 integer n integer, save :: n1 = 0 ! ! INDX = 0: This is the first call. ! if ( indx == 0 ) then i_save = 0 j_save = 0 k = n / 2 k1 = k n1 = n ! ! INDX < 0: The user is returning the results of a comparison. ! else if ( indx < 0 ) then if ( indx == -2 ) then if ( isgn < 0 ) then i_save = i_save + 1 end if j_save = k1 k1 = i_save indx = -1 i = i_save j = j_save return end if if ( 0 < isgn ) then indx = 2 i = i_save j = j_save return end if if ( k <= 1 ) then if ( n1 == 1 ) then i_save = 0 j_save = 0 indx = 0 else i_save = n1 n1 = n1 - 1 j_save = 1 indx = 1 end if i = i_save j = j_save return end if k = k - 1 k1 = k ! ! 0 < INDX, the user was asked to make an interchange. ! else if ( indx == 1 ) then k1 = k end if do i_save = 2 * k1 if ( i_save == n1 ) then j_save = k1 k1 = i_save indx = -1 i = i_save j = j_save return else if ( i_save <= n1 ) then j_save = i_save + 1 indx = -2 i = i_save j = j_save return end if if ( k <= 1 ) then exit end if k = k - 1 k1 = k end do if ( n1 == 1 ) then i_save = 0 j_save = 0 indx = 0 i = i_save j = j_save else i_save = n1 n1 = n1 - 1 j_save = 1 indx = 1 i = i_save j = j_save end if 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 MIT license. ! ! Modified: ! ! 18 May 2013 ! ! Author: ! ! John Burkardt ! implicit none character ( len = 8 ) ampm integer d integer h integer m integer mm character ( len = 9 ), parameter, dimension(12) :: month = (/ & 'January ', 'February ', 'March ', 'April ', & 'May ', 'June ', 'July ', 'August ', & 'September', 'October ', 'November ', 'December ' /) integer n integer s integer values(8) integer 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