function r8_acos ( c ) c*********************************************************************72 c cc R8_ACOS computes the arc cosine function, with argument truncation. c c Discussion: c c If you call your system ACOS routine with an input argument that is c even slightly outside the range [-1.0, 1.0 ], you may get an unpleasant c surprise (I did). c c This routine simply truncates arguments outside the range. c c Licensing: c c This code is distributed under the MIT license. c c Modified: c c 19 October 2012 c c Author: c c John Burkardt c c Parameters: c c Input, double precision C, the argument. c c Output, double precision R8_ACOS, an angle whose cosine is C. c implicit none double precision c double precision c2 double precision r8_acos c2 = c c2 = max ( c2, -1.0D+00 ) c2 = min ( c2, +1.0D+00 ) r8_acos = acos ( c2 ) return end subroutine r8_swap ( x, y ) c*********************************************************************72 c cc R8_SWAP switches two R8's. c c Licensing: c c This code is distributed under the MIT license. c c Modified: c c 30 November 1998 c c Author: c c John Burkardt c c Parameters: c c Input/output, double precision X, Y. On output, the values of X and c Y have been interchanged. c implicit none double precision x double precision y double precision z z = x x = y y = z return end function r8mat_det_4d ( a ) c*********************************************************************72 c cc R8MAT_DET_4D computes the determinant of a 4 by 4 R8MAT. c c Discussion: c c An R8MAT is an array of R8's. c c Licensing: c c This code is distributed under the MIT license. c c Modified: c c 31 August 2008 c c Author: c c John Burkardt c c Parameters: c c Input, double precision A(4,4), the matrix whose determinant is desired. c c Output, double precision R8MAT_DET_4D, the determinant of the matrix. c implicit none double precision a(4,4) double precision r8mat_det_4d r8mat_det_4d = & a(1,1) * ( & a(2,2) * ( a(3,3) * a(4,4) - a(3,4) * a(4,3) ) & - a(2,3) * ( a(3,2) * a(4,4) - a(3,4) * a(4,2) ) & + a(2,4) * ( a(3,2) * a(4,3) - a(3,3) * a(4,2) ) ) & - a(1,2) * ( & a(2,1) * ( a(3,3) * a(4,4) - a(3,4) * a(4,3) ) & - a(2,3) * ( a(3,1) * a(4,4) - a(3,4) * a(4,1) ) & + a(2,4) * ( a(3,1) * a(4,3) - a(3,3) * a(4,1) ) ) & + a(1,3) * ( & a(2,1) * ( a(3,2) * a(4,4) - a(3,4) * a(4,2) ) & - a(2,2) * ( a(3,1) * a(4,4) - a(3,4) * a(4,1) ) & + a(2,4) * ( a(3,1) * a(4,2) - a(3,2) * a(4,1) ) ) & - a(1,4) * ( & a(2,1) * ( a(3,2) * a(4,3) - a(3,3) * a(4,2) ) & - a(2,2) * ( a(3,1) * a(4,3) - a(3,3) * a(4,1) ) & + a(2,3) * ( a(3,1) * a(4,2) - a(3,2) * a(4,1) ) ) return end subroutine r8mat_solve ( n, rhs_num, a, info ) c*********************************************************************72 c cc R8MAT_SOLVE uses Gauss-Jordan elimination to solve an N by N linear system. c c Discussion: c c An R8MAT is an array of R8's. c c Licensing: c c This code is distributed under the MIT license. c c Modified: c c 08 July 2009 c c Author: c c John Burkardt c c Parameters: c c Input, integer N, the order of the matrix. c c Input, integer RHS_NUM, the number of right hand sides. c RHS_NUM must be at least 0. c c Input/output, double precision A(N,N+rhs_num), contains in rows and c columns 1 to N the coefficient matrix, and in columns N+1 through c N+rhs_num, the right hand sides. On output, the coefficient matrix c area has been destroyed, while the right hand sides have c been overwritten with the corresponding solutions. c c Output, integer INFO, singularity flag. c 0, the matrix was not singular, the solutions were computed; c J, factorization failed on step J, and the solutions could not c be computed. c implicit none integer n integer rhs_num double precision a(n,n+rhs_num) double precision apivot double precision factor integer i integer info integer ipivot integer j integer k info = 0 do j = 1, n c c Choose a pivot row. c ipivot = j apivot = a(j,j) do i = j+1, n if ( abs ( apivot ) .lt. abs ( a(i,j) ) ) then apivot = a(i,j) ipivot = i end if end do if ( apivot .eq. 0.0D+00 ) then info = j return end if c c Interchange. c do i = 1, n + rhs_num call r8_swap ( a(ipivot,i), a(j,i) ) end do c c A(J,J) becomes 1. c a(j,j) = 1.0D+00 do k = j + 1, n + rhs_num a(j,k) = a(j,k) / apivot end do c c A(I,J) becomes 0. c do i = 1, n if ( i .ne. j ) then factor = a(i,j) a(i,j) = 0.0D+00 do k = j + 1, n + rhs_num a(i,k) = a(i,k) - factor * a(j,k) end do end if end do end do return end subroutine r8mat_transpose_print ( m, n, a, title ) c*********************************************************************72 c cc R8MAT_TRANSPOSE_PRINT prints an R8MAT, transposed. c c Discussion: c c An R8MAT is an array of R8's. c c Licensing: c c This code is distributed under the MIT license. c c Modified: c c 28 April 2008 c c Author: c c John Burkardt c c Parameters: c c Input, integer M, N, the number of rows and columns. c c Input, double precision A(M,N), an M by N matrix to be printed. c c Input, character*(*) TITLE, a title. c implicit none integer m integer n double precision a(m,n) character*(*) 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 ) c*********************************************************************72 c cc R8MAT_TRANSPOSE_PRINT_SOME prints some of an R8MAT transposed. c c Discussion: c c An R8MAT is an array of R8's. c c Licensing: c c This code is distributed under the MIT license. c c Modified: c c 28 April 2008 c c Author: c c John Burkardt c c Parameters: c c Input, integer M, N, the number of rows and columns. c c Input, double precision A(M,N), an M by N matrix to be printed. c c Input, integer ILO, JLO, the first row and column to print. c c Input, integer IHI, JHI, the last row and column to print. c c Input, character * ( * ) TITLE, a title. c implicit none integer incx parameter ( incx = 5 ) integer m integer n double precision a(m,n) character * ( 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 * ( * ) title write ( *, '(a)' ) ' ' write ( *, '(a)' ) trim ( title ) if ( m .le. 0 .or. n .le. 0 ) then write ( *, '(a)' ) ' ' write ( *, '(a)' ) ' (None)' return 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), '(i8,6x)') i end do write ( *, '('' Row'',5a14)' ) ctemp(1:inc) write ( *, '(a)' ) ' Col' 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 ( *, '(2x,i8,a,5a14)' ) j, ':', ( ctemp(i), i = 1, inc ) end do end do return end subroutine r8vec_angle_3d ( u, v, angle ) c*********************************************************************72 c cc R8VEC_ANGLE_3D computes the angle between two vectors in 3D. c c Licensing: c c This code is distributed under the MIT license. c c Modified: c c 07 July 2009 c c Author: c c John Burkardt c c Parameters: c c Input, double precision U(3), V(3), the vectors. c c Output, double precision ANGLE, the angle between the two vectors. c implicit none double precision angle double precision angle_cos double precision r8_acos double precision r8vec_dot_product double precision r8vec_norm double precision u(3) double precision u_norm double precision uv_dot double precision v(3) double precision v_norm uv_dot = r8vec_dot_product ( 3, u, v ) u_norm = r8vec_norm ( 3, u ) v_norm = r8vec_norm ( 3, v ) angle_cos = uv_dot / u_norm / v_norm angle = r8_acos ( angle_cos ) return end subroutine r8vec_cross_3d ( v1, v2, v3 ) c*********************************************************************72 c cc R8VEC_CROSS_3D computes the cross product of two vectors in 3D. c c Discussion: c c The cross product in 3D can be regarded as the determinant of the c symbolic matrix: c c | i j k | c det | x1 y1 z1 | c | x2 y2 z2 | c c = ( y1 * z2 - z1 * y2 ) * i c + ( z1 * x2 - x1 * z2 ) * j c + ( x1 * y2 - y1 * x2 ) * k c c Licensing: c c This code is distributed under the MIT license. c c Modified: c c 07 August 2005 c c Author: c c John Burkardt c c Parameters: c c Input, double precision V1(3), V2(3), the two vectors. c c Output, double precision V3(3), the cross product vector. c implicit none double precision v1(3) double precision v2(3) double precision v3(3) v3(1) = v1(2) * v2(3) - v1(3) * v2(2) v3(2) = v1(3) * v2(1) - v1(1) * v2(3) v3(3) = v1(1) * v2(2) - v1(2) * v2(1) return end function r8vec_dot_product ( n, v1, v2 ) c*********************************************************************72 c cc R8VEC_DOT_PRODUCT finds the dot product of a pair of R8VEC's. c c Discussion: c c An R8VEC is a vector of R8 values. c c In FORTRAN90, the system routine DOT_PRODUCT should be called c directly. c c Licensing: c c This code is distributed under the MIT license. c c Modified: c c 27 May 2008 c c Author: c c John Burkardt c c Parameters: c c Input, integer N, the dimension of the vectors. c c Input, double precision V1(N), V2(N), the vectors. c c Output, double precision R8VEC_DOT_PRODUCT, the dot product. c implicit none integer n integer i double precision r8vec_dot_product double precision v1(n) double precision v2(n) double precision value value = 0.0D+00 do i = 1, n value = value + v1(i) * v2(i) end do r8vec_dot_product = value return end function r8vec_max ( n, a ) c*********************************************************************72 c cc R8VEC_MAX() returns the maximum value in an R8VEC. c c Discussion: c c An R8VEC is a vector of R8's. c c Licensing: c c This code is distributed under the MIT license. c c Modified: c c 31 May 2009 c c Author: c c John Burkardt c c Parameters: c c Input, integer N, the number of entries in the array. c c Input, double precision A(N), the array. c c Output, double precision AMAX, the value of the largest entry. c implicit none integer n double precision a(n) integer i double precision r8vec_max double precision value value = a(1) do i = 2, n value = max ( value, a(i) ) end do r8vec_max = value return end function r8vec_norm ( n, a ) c*********************************************************************72 c cc R8VEC_NORM() returns the L2 norm of an R8VEC. c c Discussion: c c An R8VEC is a vector of R8 values. c c The vector L2 norm is defined as: c c R8VEC_NORM = sqrt ( sum ( 1 <= I <= N ) A(I)^2 ). c c Licensing: c c This code is distributed under the MIT license. c c Modified: c c 27 May 2008 c c Author: c c John Burkardt c c Parameters: c c Input, integer N, the number of entries in A. c c Input, double precision A(N), the vector whose L2 norm is desired. c c Output, double precision R8VEC_NORM, the L2 norm of A. c implicit none integer n double precision a(n) integer i double precision r8vec_norm double precision value value = 0.0D+00 do i = 1, n value = value + a(i) * a(i) end do value = sqrt ( value ) r8vec_norm = value return end function r8vec_normsq ( n, a ) c*********************************************************************72 c cc R8VEC_NORMSQ() returns the square of the L2 norm of an R8VEC. c c Discussion: c c An R8VEC is a vector of R8's. c c The square of the vector L2 norm is defined as: c c R8VEC_NORMSQ = sum ( 1 <= I <= N ) V(I)^2. c c Licensing: c c This code is distributed under the MIT license. c c Modified: c c 28 October 2010 c c Author: c c John Burkardt c c Parameters: c c Input, integer N, the vector dimension. c c Input, double precision A(N), the vector. c c Output, double precision R8VEC_NORMSQ, the squared L2 norm. c implicit none integer n double precision a(n) integer i double precision r8vec_normsq double precision value value = 0.0D+00 do i = 1, n value = value + a(i) * a(i) end do r8vec_normsq = value return end subroutine r8vec_print ( n, a, title ) c*********************************************************************72 c cc R8VEC_PRINT prints an R8VEC. c c Discussion: c c An R8VEC is a vector of R8's. c c Licensing: c c This code is distributed under the MIT license. c c Modified: c c 12 January 2007 c c Author: c c John Burkardt c c Parameters: c c Input, integer N, the number of components of the vector. c c Input, double precision A(N), the vector to be printed. c c Input, character * ( * ) TITLE, a title. c implicit none integer n double precision a(n) integer i character ( len = * ) title write ( *, '(a)' ) ' ' write ( *, '(a)' ) trim ( title ) write ( *, '(a)' ) ' ' do i = 1, n write ( *, '(2x,i8,a,1x,g16.8)' ) i, ':', a(i) end do return end subroutine r8vec_transpose_print ( n, a, title ) c*********************************************************************72 c cc R8VEC_TRANSPOSE_PRINT prints an R8VEC "transposed". c c Discussion: c c An R8VEC is a vector of R8's. c c Example: c c A = (/ 1.0, 2.1, 3.2, 4.3, 5.4, 6.5, 7.6, 8.7, 9.8, 10.9, 11.0 /) c TITLE = 'My vector: ' c c My vector: 1.0 2.1 3.2 4.3 5.4 c 6.5 7.6 8.7 9.8 10.9 c 11.0 c c Licensing: c c This code is distributed under the MIT license. c c Modified: c c 12 November 2010 c c Author: c c John Burkardt c c Parameters: c c Input, integer N, the number of components of the vector. c c Input, double precision A(N), the vector to be printed. c c Input, character * ( * ) TITLE, a title. c implicit none integer n double precision a(n) integer i integer ihi integer ilo character * ( * ) title integer title_length title_length = len_trim ( title ) do ilo = 1, n, 5 ihi = min ( ilo + 5 - 1, n ) if ( ilo .eq. 1 ) then write ( *, '(a,2x,5g14.6)' ) title(1:title_length), a(ilo:ihi) else write ( *, '(a,2x,5g14.6)' ) & ( ' ', i = 1, title_length ), a(ilo:ihi) end if end do return end subroutine tetrahedron_centroid ( tetra, centroid ) c*********************************************************************72 c cc TETRAHEDRON_CENTROID computes the centroid of a tetrahedron. c c Licensing: c c This code is distributed under the MIT license. c c Modified: c c 30 December 2004 c c Author: c c John Burkardt c c Parameters: c c Input, double precision TETRA(3,4) the tetrahedron vertices. c c Output, double precision CENTROID(3), the coordinates of the centroid. c implicit none double precision centroid(3) integer i integer j double precision t double precision tetra(3,4) do i = 1, 3 t = 0.0D+00 do j = 1, 4 t = t + tetra(i,j) end do centroid(i) = t / 4.0D+00 end do return end subroutine tetrahedron_circumsphere ( tetra, r, pc ) c*********************************************************************72 c cc TETRAHEDRON_CIRCUMSPHERE computes the circumsphere of a tetrahedron. c c Discussion: c c The circumsphere, or circumscribed sphere, of a tetrahedron is the c sphere that passes through the four vertices. The circumsphere is c not necessarily the smallest sphere that contains the tetrahedron. c c Surprisingly, the diameter of the sphere can be found by solving c a 3 by 3 linear system. This is because the vectors P2 - P1, c P3 - P1 and P4 - P1 are secants of the sphere, and each forms a c right triangle with the diameter through P1. Hence, the dot product of c P2 - P1 with that diameter is equal to the square of the length c of P2 - P1, and similarly for P3 - P1 and P4 - P1. This determines c the diameter vector originating at P1, and hence the radius and c center. c c Licensing: c c This code is distributed under the MIT license. c c Modified: c c 10 August 2005 c c Author: c c John Burkardt c c Reference: c c Adrian Bowyer, John Woodwark, c A Programmer's Geometry, c Butterworths, 1983. c c Parameters: c c Input, double precision TETRA(3,4) the tetrahedron vertices. c c Output, double precision R, PC(3), the center of the c circumscribed sphere, and its radius. If the linear system is c singular, then R = -1, PC(1:3) = 0. c implicit none double precision a(3,4) integer i integer info integer j double precision pc(3) double precision r double precision tetra(3,4) c c Set up the linear system. c do j = 1, 3 do i = 1, 3 a(i,j) = tetra(i,i+1) end do end do do j = 1, 3 do i = 1, 3 a(i,j) = a(i,j) - tetra(j,1) end do end do do i = 1, 3 a(i,4) = 0.0D+00 do j = 1, 3 a(i,4) = a(i,4) + a(i,j)**2 end do end do c c Solve the linear system. c call r8mat_solve ( 3, 1, a, info ) c c If the system was singular, return a consolation prize. c if ( info .ne. 0 ) then r = -1.0D+00 do i = 1, 3 pc(i) = 0.0D+00 end do return end if c c Compute the radius and center. c r = 0.0D+00 do i = 1, 3 r = r + a(i,4)**2 end do r = 0.5D+00 * sqrt ( r ) do i = 1, 3 pc(i) = tetra(i,1) + 0.5D+00 * a(i,4) end do return end subroutine tetrahedron_dihedral_angles ( tetra, angle ) c*********************************************************************72 c cc TETRAHEDRON_DIHEDRAL_ANGLES computes dihedral angles of a tetrahedron. c c Licensing: c c This code is distributed under the MIT license. c c Modified: c c 07 July 2009 c c Author: c c John Burkardt c c Parameters: c c Input, double precision TETRA(3,4), the vertices of the tetrahedron. c c Output, double precision ANGLE(6), the dihedral angles along the c axes AB, AC, AD, BC, BD and CD, respectively. c implicit none double precision ab(3) double precision abc_normal(3) double precision abd_normal(3) double precision ac(3) double precision acd_normal(3) double precision ad(3) double precision angle(6) double precision bc(3) double precision bcd_normal(3) double precision bd(3) double precision cd(3) integer i double precision r8_pi parameter ( r8_pi = 3.141592653589793D+00 ) double precision tetra(3,4) call tetrahedron_edges ( tetra, ab, ac, ad, bc, bd, cd ) call r8vec_cross_3d ( ac, ab, abc_normal ) call r8vec_cross_3d ( ab, ad, abd_normal ) call r8vec_cross_3d ( ad, ac, acd_normal ) call r8vec_cross_3d ( bc, bd, bcd_normal ) call r8vec_angle_3d ( abc_normal, abd_normal, angle(1) ) call r8vec_angle_3d ( abc_normal, acd_normal, angle(2) ) call r8vec_angle_3d ( abd_normal, acd_normal, angle(3) ) call r8vec_angle_3d ( abc_normal, bcd_normal, angle(4) ) call r8vec_angle_3d ( abd_normal, bcd_normal, angle(5) ) call r8vec_angle_3d ( acd_normal, bcd_normal, angle(6) ) do i = 1, 6 angle(i) = r8_pi - angle(i) end do return end subroutine tetrahedron_edges ( tetra, ab, ac, ad, bc, bd, cd ) c*********************************************************************72 c cc TETRAHEDRON_EDGES computes the edges of a tetrahedron. c c Discussion: c c The vertices are A, B, C, D. The edge from A to B is denoted by AB. c c Licensing: c c This code is distributed under the MIT license. c c Modified: c c 11 May 2014 c c Author: c c John Burkardt. c c Parameters: c c Input, double precision TETRA(3,4), the vertices of the tetrahedron. c c Output, double precision AB(3), AC(3), AD(3), BC(3), BD(3), CD(3), c vectors that represent the edges of the tetrahedron. c implicit none double precision ab(3) double precision ac(3) double precision ad(3) double precision bc(3) double precision bd(3) double precision cd(3) integer i double precision tetra(3,4) do i = 1, 3 ab(i) = tetra(i,2) - tetra(i,1) ac(i) = tetra(i,3) - tetra(i,1) ad(i) = tetra(i,4) - tetra(i,1) bc(i) = tetra(i,3) - tetra(i,2) bd(i) = tetra(i,4) - tetra(i,2) cd(i) = tetra(i,4) - tetra(i,3) end do return end subroutine tetrahedron_edge_length ( tetra, edge_length ) c*********************************************************************72 c cc TETRAHEDRON_EDGE_LENGTH returns edge lengths of a tetrahedron. c c Licensing: c c This code is distributed under the MIT license. c c Modified: c c 09 August 2005 c c Author: c c John Burkardt c c Parameters: c c Input, double precision TETRA(3,4), the tetrahedron vertices. c c Output, double precision EDGE_LENGTH(6), the length of the edges. c implicit none double precision edge_length(6) integer i integer j1 integer j2 integer k double precision r8vec_norm double precision t(3) double precision tetra(3,4) k = 0 do j1 = 1, 3 do j2 = j1 + 1, 4 k = k + 1 do i = 1, 3 t(i) = tetra(i,j2) - tetra(i,j1) end do edge_length(k) = r8vec_norm ( 3, t ) end do end do return end subroutine tetrahedron_face_angles ( tetra, angles ) c*********************************************************************72 c cc TETRAHEDRON_FACE_ANGLES returns the 12 face angles of a tetrahedron. c c Discussion: c c The tetrahedron has 4 triangular faces. This routine computes the c 3 planar angles associated with each face. c c Licensing: c c This code is distributed under the MIT license. c c Modified: c c 03 July 2009 c c Author: c c John Burkardt c c Parameters: c c Input, double precision TETRA(3,4) the tetrahedron vertices. c c Output, double precision ANGLES(3,4), the face angles. c implicit none double precision angles(3,4) integer i integer j double precision tri(3,3) double precision tetra(3,4) c c Face 123 c do j = 1, 3 do i = 1, 3 tri(i,j) = tetra(i,j) end do end do call triangle_angles_3d ( tri, angles(1:3,1) ) c c Face 124 c do i = 1, 3 tri(i,1) = tetra(i,1) end do do i = 1, 3 tri(i,2) = tetra(i,2) end do do i = 1, 3 tri(i,3) = tetra(i,4) end do call triangle_angles_3d ( tri, angles(1:3,2) ) c c Face 134 c do i = 1, 3 tri(i,1) = tetra(i,1) end do do i = 1, 3 tri(i,2) = tetra(i,3) end do do i = 1, 3 tri(i,3) = tetra(i,4) end do call triangle_angles_3d ( tri, angles(1:3,3) ) c c Face 234 c do i = 1, 3 tri(i,1) = tetra(i,2) end do do i = 1, 3 tri(i,2) = tetra(i,3) end do do i = 1, 3 tri(i,3) = tetra(i,4) end do call triangle_angles_3d ( tri, angles(1:3,4) ) return end subroutine tetrahedron_face_areas ( tetra, areas ) c*********************************************************************72 c cc TETRAHEDRON_FACE_AREAS returns the 4 face areas of a tetrahedron. c c Discussion: c c The tetrahedron has 4 triangular faces. This routine computes the c area of each face. c c Licensing: c c This code is distributed under the MIT license. c c Modified: c c 07 July 2009 c c Author: c c John Burkardt c c Parameters: c c Input, double precision TETRA(3,4) the tetrahedron vertices. c c Output, double precision AREAS(4), the face areas. c implicit none double precision areas(4) integer i integer j double precision tri(3,3) double precision tetra(3,4) c c Face 123 c do j = 1, 3 do i = 1, 3 tri(i,j) = tetra(i,j) end do end do call triangle_area_3d ( tri, areas(1) ) c c Face 124 c do i = 1, 3 tri(i,1) = tetra(i,1) end do do i = 1, 3 tri(i,2) = tetra(i,2) end do do i = 1, 3 tri(i,3) = tetra(i,4) end do call triangle_area_3d ( tri, areas(2) ) c c Face 134 c do i = 1, 3 tri(i,1) = tetra(i,1) end do do i = 1, 3 tri(i,2) = tetra(i,3) end do do i = 1, 3 tri(i,3) = tetra(i,4) end do call triangle_area_3d ( tri, areas(3) ) c c Face 234 c do i = 1, 3 tri(i,1) = tetra(i,2) end do do i = 1, 3 tri(i,2) = tetra(i,3) end do do i = 1, 3 tri(i,3) = tetra(i,4) end do call triangle_area_3d ( tri, areas(4) ) return end subroutine tetrahedron_insphere ( tetra, r, pc ) c*********************************************************************72 c cc TETRAHEDRON_INSPHERE finds the insphere of a tetrahedron. c c Discussion: c c The insphere of a tetrahedron is the inscribed sphere, which touches c each face of the tetrahedron at a single point. c c The points of contact are the centroids of the triangular faces c of the tetrahedron. Therefore, the point of contact for a face c can be computed as the average of the vertices of that face. c c The sphere can then be determined as the unique sphere through c the four given centroids. c c Licensing: c c This code is distributed under the MIT license. c c Modified: c c 08 August 2005 c c Author: c c John Burkardt c c Reference: c c Philip Schneider, David Eberly, c Geometric Tools for Computer Graphics, c Elsevier, 2002, c ISBN: 1558605940, c LC: T385.G6974. c c Parameters: c c Input, double precision TETRA(3,4), the vertices of the tetrahedron. c c Output, double precision R, PC(3), the radius and the center c of the sphere. c implicit none double precision b(4,4) double precision r8mat_det_4d double precision r8vec_norm double precision gamma integer i integer j double precision l123 double precision l124 double precision l134 double precision l234 double precision n123(3) double precision n124(3) double precision n134(3) double precision n234(3) double precision pc(3) double precision r double precision tetra(3,4) double precision v21(3) double precision v31(3) double precision v41(3) double precision v32(3) double precision v42(3) double precision v43(3) call tetrahedron_edges ( tetra, v21, v31, v41, v32, v42, v43 ) call r8vec_cross_3d ( v21, v31, n123 ) call r8vec_cross_3d ( v41, v21, n124 ) call r8vec_cross_3d ( v31, v41, n134 ) call r8vec_cross_3d ( v42, v32, n234 ) l123 = r8vec_norm ( 3, n123 ) l124 = r8vec_norm ( 3, n124 ) l134 = r8vec_norm ( 3, n134 ) l234 = r8vec_norm ( 3, n234 ) do i = 1, 3 pc(i) = ( l234 * tetra(i,1) & + l134 * tetra(i,2) & + l124 * tetra(i,3) & + l123 * tetra(i,4) ) & / ( l234 + l134 + l124 + l123 ) end do do j = 1, 4 do i = 1, 3 b(i,j) = tetra(i,j) end do b(4,j) = 1.0D+00 end do gamma = abs ( r8mat_det_4d ( b ) ) r = gamma / ( l234 + l134 + l124 + l123 ) return end subroutine tetrahedron_quality1 ( tetra, quality ) c*********************************************************************72 c cc TETRAHEDRON_QUALITY1: "quality" of a tetrahedron. c c Discussion: c c The quality of a tetrahedron is 3 times the ratio of the radius of c the inscribed sphere divided by that of the circumscribed sphere. c c An equilateral tetrahredron achieves the maximum possible quality of 1. c c Licensing: c c This code is distributed under the MIT license. c c Modified: c c 09 August 2005 c c Author: c c John Burkardt c c Parameters: c c Input, double precision TETRA(3,4), the tetrahedron vertices. c c Output, double precision QUALITY, the quality of the tetrahedron. c implicit none double precision pc(3) double precision quality double precision r_in double precision r_out double precision tetra(3,4) call tetrahedron_circumsphere ( tetra, r_out, pc ) call tetrahedron_insphere ( tetra, r_in, pc ) quality = 3.0D+00 * r_in / r_out return end subroutine tetrahedron_quality2 ( tetra, quality2 ) c*********************************************************************72 c cc TETRAHEDRON_QUALITY2: "quality" of a tetrahedron. c c Discussion: c c The quality measure #2 of a tetrahedron is: c c QUALITY2 = 2 * sqrt ( 6 ) * RIN / LMAX c c where c c RIN = radius of the inscribed sphere; c LMAX = length of longest side of the tetrahedron. c c An equilateral tetrahredron achieves the maximum possible quality of 1. c c Licensing: c c This code is distributed under the MIT license. c c Modified: c c 16 August 2005 c c Author: c c John Burkardt c c Reference: c c Qiang Du, Desheng Wang, c The Optimal Centroidal Voronoi Tesselations and the Gersho's c Conjecture in the Three-Dimensional Space, c Computers and Mathematics with Applications, c Volume 49, 2005, pages 1355-1373. c c Parameters: c c Input, double precision TETRA(3,4), the tetrahedron vertices. c c Output, double precision QUALITY2, the quality of the tetrahedron. c implicit none double precision edge_length(6) double precision l_max double precision pc(3) double precision quality2 double precision r_in double precision r8vec_max double precision tetra(3,4) call tetrahedron_edge_length ( tetra, edge_length ) l_max = r8vec_max ( 6, edge_length ) call tetrahedron_insphere ( tetra, r_in, pc ) quality2 = 2.0D+00 * sqrt ( 6.0D+00 ) * r_in / l_max return end subroutine tetrahedron_quality3 ( tetra, quality3 ) c*********************************************************************72 c cc TETRAHEDRON_QUALITY3 computes the mean ratio of a tetrahedron. c c Discussion: c c This routine computes QUALITY3, the eigenvalue or mean ratio of c a tetrahedron. c c QUALITY3 = 12 * ( 3 * volume )^(2/3) / (sum of squares of edge lengths). c c This value may be used as a shape quality measure for the tetrahedron. c c For an equilateral tetrahedron, the value of this quality measure c will be 1. For any other tetrahedron, the value will be between c 0 and 1. c c Licensing: c c This code is distributed under the MIT license. c c Modified: c c 17 August 2005 c c Author: c c Original FORTRAN77 version by Barry Joe. c FORTRAN90 version by John Burkardt. c c Reference: c c Barry Joe, c GEOMPACK - a software package for the generation of meshes c using geometric algorithms, c Advances in Engineering Software, c Volume 13, pages 325-331, 1991. c c Parameters: c c Input, double precision TETRA(3,4), the vertices of the tetrahedron. c c Output, double precision QUALITY3, the mean ratio of the tetrahedron. c implicit none double precision ab(3) double precision ac(3) double precision ad(3) double precision bc(3) double precision bd(3) double precision cd(3) double precision denom double precision lab double precision lac double precision lad double precision lbc double precision lbd double precision lcd double precision quality3 double precision r8vec_normsq double precision tetra(3,4) double precision volume c c Compute the vectors representing the sides of the tetrahedron. c call tetrahedron_edges ( tetra, ab, ac, ad, bc, bd, cd ) c c Compute the squares of the lengths of the sides. c lab = r8vec_normsq ( 3, ab ) lac = r8vec_normsq ( 3, ac ) lad = r8vec_normsq ( 3, ad ) lbc = r8vec_normsq ( 3, bc ) lbd = r8vec_normsq ( 3, bd ) lcd = r8vec_normsq ( 3, cd ) c c Compute the volume. c volume = abs ( & ab(1) * ( ac(2) * ad(3) - ac(3) * ad(2) ) & + ab(2) * ( ac(3) * ad(1) - ac(1) * ad(3) ) & + ab(3) * ( ac(1) * ad(2) - ac(2) * ad(1) ) ) / 6.0D+00 denom = lab + lac + lad + lbc + lbd + lcd if ( denom .eq. 0.0D+00 ) then quality3 = 0.0D+00 else quality3 = 12.0D+00 & * ( 3.0D+00 * volume )**( 2.0D+00 / 3.0D+00 ) / denom end if return end subroutine tetrahedron_quality4 ( tetra, quality4 ) c*********************************************************************72 c cc TETRAHEDRON_QUALITY4 computes the minimum solid angle of a tetrahedron. c c Discussion: c c This routine computes a quality measure for a tetrahedron, based c on the sine of half the minimum of the four solid angles. c c The quality measure for an equilateral tetrahedron should be 1, c since the solid angles of such a tetrahedron are each equal to pi. c c Licensing: c c This code is distributed under the MIT license. c c Modified: c c 17 August 2005 c c Author: c c Original FORTRAN77 version by Barry Joe. c This version by John Burkardt. c c Reference: c c Barry Joe, c GEOMPACK - a software package for the generation of meshes c using geometric algorithms, c Advances in Engineering Software, c Volume 13, pages 325-331, 1991. c c Parameters: c c Input, double precision TETRA(3,4), the vertices of the tetrahedron. c c Output, double precision QUALITY4, the value of the quality measure. c implicit none double precision ab(3) double precision ac(3) double precision ad(3) double precision bc(3) double precision bd(3) double precision cd(3) double precision denom double precision l1 double precision l2 double precision l3 double precision lab double precision lac double precision lad double precision lbc double precision lbd double precision lcd double precision quality4 double precision r8vec_norm double precision tetra(3,4) double precision volume c c Compute the vectors that represent the sides. c call tetrahedron_edges ( tetra, ab, ac, ad, bc, bd, cd ) c c Compute the lengths of the sides. c lab = r8vec_norm ( 3, ab ) lac = r8vec_norm ( 3, ac ) lad = r8vec_norm ( 3, ad ) lbc = r8vec_norm ( 3, bc ) lbd = r8vec_norm ( 3, bd ) lcd = r8vec_norm ( 3, cd ) c c Compute the volume c volume = abs ( & ab(1) * ( ac(2) * ad(3) - ac(3) * ad(2) ) & + ab(2) * ( ac(3) * ad(1) - ac(1) * ad(3) ) & + ab(3) * ( ac(1) * ad(2) - ac(2) * ad(1) ) ) / 6.0D+00 quality4 = 1.0D+00 l1 = lab + lac l2 = lab + lad l3 = lac + lad denom = ( l1 + lbc ) * ( l1 - lbc ) & * ( l2 + lbd ) * ( l2 - lbd ) & * ( l3 + lcd ) * ( l3 - lcd ) if ( denom .le. 0.0D+00 ) then quality4 = 0.0D+00 else quality4 = min ( quality4, 12.0D+00 * volume / sqrt ( denom ) ) end if l1 = lab + lbc l2 = lab + lbd l3 = lbc + lbd denom = ( l1 + lac ) * ( l1 - lac ) & * ( l2 + lad ) * ( l2 - lad ) & * ( l3 + lcd ) * ( l3 - lcd ) if ( denom .le. 0.0D+00 ) then quality4 = 0.0D+00 else quality4 = min ( quality4, 12.0D+00 * volume / sqrt ( denom ) ) end if l1 = lac + lbc l2 = lac + lcd l3 = lbc + lcd denom = ( l1 + lab ) * ( l1 - lab ) & * ( l2 + lad ) * ( l2 - lad ) & * ( l3 + lbd ) * ( l3 - lbd ) if ( denom .le. 0.0D+00 ) then quality4 = 0.0D+00 else quality4 = min ( quality4, 12.0D+00 * volume / sqrt ( denom ) ) end if l1 = lad + lbd l2 = lad + lcd l3 = lbd + lcd denom = ( l1 + lab ) * ( l1 - lab ) & * ( l2 + lac ) * ( l2 - lac ) & * ( l3 + lbc ) * ( l3 - lbc ) if ( denom .le. 0.0D+00 ) then quality4 = 0.0D+00 else quality4 = min ( quality4, 12.0D+00 * volume / sqrt ( denom ) ) end if quality4 = quality4 * 1.5D+00 * sqrt ( 6.0D+00 ) return end subroutine tetrahedron_solid_angles ( tetra, angle ) c*********************************************************************72 c cc TETRAHEDRON_SOLID_ANGLES computes solid angles of a tetrahedron. c c Licensing: c c This code is distributed under the MIT license. c c Modified: c c 07 July 2009 c c Author: c c John Burkardt c c Parameters: c c Input, double precision TETRA(3,4), the vertices of the tetrahedron. c c Output, double precision ANGLE(4), the solid angles. c implicit none double precision angle(4) double precision dihedral_angles(6) double precision, parameter :: r8_pi = 3.141592653589793D+00 double precision tetra(3,4) call tetrahedron_dihedral_angles ( tetra, dihedral_angles ) angle(1) = dihedral_angles(1) & + dihedral_angles(2) & + dihedral_angles(3) - r8_pi angle(2) = dihedral_angles(1) & + dihedral_angles(4) & + dihedral_angles(5) - r8_pi angle(3) = dihedral_angles(2) & + dihedral_angles(4) & + dihedral_angles(6) - r8_pi angle(4) = dihedral_angles(3) & + dihedral_angles(5) & + dihedral_angles(6) - r8_pi return end subroutine tetrahedron_volume ( tetra, volume ) c*********************************************************************72 c cc TETRAHEDRON_VOLUME computes the volume of a tetrahedron. c c Licensing: c c This code is distributed under the MIT license. c c Modified: c c 30 December 2004 c c Author: c c John Burkardt c c Parameters: c c Input, double precision TETRA(3,4), the vertices of the tetrahedron. c c Output, double precision VOLUME, the volume of the tetrahedron. c implicit none double precision a(4,4) integer i integer j double precision r8mat_det_4d double precision tetra(3,4) double precision volume do j = 1, 4 do i = 1, 3 a(i,j) = tetra(i,j) end do a(4,j) = 1.0D+00 end do volume = abs ( r8mat_det_4d ( a ) ) / 6.0D+00 return end subroutine triangle_angles_3d ( t, angle ) c*********************************************************************72 c cc TRIANGLE_ANGLES_3D computes the angles of a triangle in 3D. c c Discussion: c c The law of cosines is used: c c C * C = A * A + B * B - 2 * A * B * COS ( GAMMA ) c c where GAMMA is the angle opposite side C. c c Licensing: c c This code is distributed under the MIT license. c c Modified: c c 04 May 2005 c c Author: c c John Burkardt c c Parameters: c c Input, double precision T(3,3), the triangle vertices. c c Output, double precision ANGLE(3), the angles opposite c sides P1-P2, P2-P3 and P3-P1, in radians. c implicit none double precision a double precision angle(3) double precision b double precision c integer i double precision r8_acos double precision r8_pi parameter ( r8_pi = 3.141592653589793D+00 ) double precision t(3,3) c c Compute the length of each side. c a = 0.0D+00 b = 0.0D+00 c = 0.0D+00 do i = 1, 3 a = a + ( t(i,1) - t(i,2) )**2 b = b + ( t(i,2) - t(i,3) )**2 c = c + ( t(i,3) - t(i,1) )**2 end do a = sqrt ( a ) b = sqrt ( b ) c = sqrt ( c ) c c Take care of a ridiculous special case. c if ( a .eq. 0.0D+00 .and. & b .eq. 0.0D+00 .and. & c .eq. 0.0D+00 ) then do i = 1, 3 angle(i) = 2.0D+00 * r8_pi / 3.0D+00 end do return end if if ( c .eq. 0.0D+00 .or. a .eq. 0.0D+00 ) then angle(1) = r8_pi else angle(1) = r8_acos ( ( c * c + a * a - b * b ) & / ( 2.0D+00 * c * a ) ) end if if ( a .eq. 0.0D+00 .or. b .eq. 0.0D+00 ) then angle(2) = r8_pi else angle(2) = r8_acos ( ( a * a + b * b - c * c ) & / ( 2.0D+00 * a * b ) ) end if if ( b .eq. 0.0D+00 .or. c .eq. 0.0D+00 ) then angle(3) = r8_pi else angle(3) = r8_acos ( ( b * b + c * c - a * a ) & / ( 2.0D+00 * b * c ) ) end if return end subroutine triangle_area_3d ( t, area ) c*********************************************************************72 c cc TRIANGLE_AREA_3D computes the area of a triangle in 3D. c c Discussion: c c This routine uses the fact that the norm of the cross product c of two vectors is the area of the parallelogram they form. c c Therefore, the area of the triangle is half of that value. c c Licensing: c c This code is distributed under the MIT license. c c Modified: c c 27 December 2004 c c Author: c c John Burkardt c c Reference: c c Adrian Bowyer, John Woodwark, c A Programmer's Geometry, c Butterworths, 1983. c c Parameters: c c Input, double precision T(3,3), the triangle vertices. c c Output, double precision AREA, the area of the triangle. c implicit none double precision area double precision cross(3) double precision r8vec_norm double precision t(3,3) c c Compute the cross product vector. c cross(1) = ( t(2,2) - t(2,1) ) * ( t(3,3) - t(3,1) ) & - ( t(3,2) - t(3,1) ) * ( t(2,3) - t(2,1) ) cross(2) = ( t(3,2) - t(3,1) ) * ( t(1,3) - t(1,1) ) & - ( t(1,2) - t(1,1) ) * ( t(3,3) - t(3,1) ) cross(3) = ( t(1,2) - t(1,1) ) * ( t(2,3) - t(2,1) ) & - ( t(2,2) - t(2,1) ) * ( t(1,3) - t(1,1) ) area = 0.5D+00 * r8vec_norm ( 3, cross ) return end