Reads in the bottom bathymetry
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
character(len=80), | intent(in) | :: | FileRoot | |||
character(len=1), | intent(in) | :: | BotBTY | |||
real(kind=8), | intent(in) | :: | DepthB | |||
integer, | intent(in) | :: | PRTFile |
SUBROUTINE ReadBTY( FileRoot, BotBTY, DepthB, PRTFile ) !! Reads in the bottom bathymetry INTEGER, INTENT( IN ) :: PRTFile CHARACTER (LEN= 1), INTENT( IN ) :: BotBTY REAL (KIND=8), INTENT( IN ) :: DepthB CHARACTER (LEN=80), INTENT( IN ) :: FileRoot SELECT CASE ( BotBTY ) CASE ( '~', '*' ) WRITE( PRTFile, * ) '__________________________________________________________________________' WRITE( PRTFile, * ) WRITE( PRTFile, * ) 'Using bottom-bathymetry file' OPEN( UNIT = BTYFile, FILE = TRIM( FileRoot ) // '.bty', STATUS = 'OLD', IOSTAT = IOStat, ACTION = 'READ' ) IF ( IOsTAT /= 0 ) THEN WRITE( PRTFile, * ) 'BTYFile = ', TRIM( FileRoot ) // '.bty' CALL ERROUT( 'ReadBTY', 'Unable to open bathymetry file' ) END IF READ( BTYFile, * ) btyType BathyType: SELECT CASE ( btyType( 1 : 1 ) ) CASE ( 'C' ) WRITE( PRTFile, * ) 'Curvilinear Interpolation' CASE ( 'L' ) WRITE( PRTFile, * ) 'Piecewise linear interpolation' CASE DEFAULT CALL ERROUT( 'ReadBTY', 'Unknown option for selecting bathymetry interpolation' ) END SELECT BathyType READ( BTYFile, * ) NbtyPts WRITE( PRTFile, * ) 'Number of bathymetry points = ', NbtyPts NbtyPts = NbtyPts + 2 ! we'll be extending the bathymetry to infinity on both sides ALLOCATE( Bot( NbtyPts ), Stat = IAllocStat ) IF ( IAllocStat /= 0 ) & CALL ERROUT( 'BELLHOP:ReadBTY', 'Insufficient memory for bathymetry data: reduce # bty points' ) WRITE( PRTFile, * ) BathyTypeB: SELECT CASE ( btyType( 2 : 2 ) ) CASE ( 'S', '' ) WRITE( PRTFile, * ) 'Short format (bathymetry only)' WRITE( PRTFile, * ) ' Range (km) Depth (m)' CASE ( 'L' ) WRITE( PRTFile, * ) 'Long format (bathymetry and geoacoustics)' WRITE( PRTFile, "( ' Range (km) Depth (m) alphaR (m/s) betaR rho (g/cm^3) alphaI betaI', / )" ) CASE DEFAULT CALL ERROUT( 'ReadBTY', 'Unknown option for selecting bathymetry interpolation' ) END SELECT BathyTypeB btyPt: DO ii = 2, NbtyPts - 1 SELECT CASE ( btyType( 2 : 2 ) ) CASE ( 'S', '' ) ! short format READ( BTYFile, * ) Bot( ii )%x ! LP: This condition was previously ii == NbtyPts, ! which will never be satisfied due to the loop bounds IF ( ii < Number_to_Echo .OR. ii == NbtyPts - 1 ) THEN ! echo some values WRITE( PRTFile, FMT = "(2G11.3)" ) Bot( ii )%x END IF CASE ( 'L' ) ! long format READ( BTYFile, * ) Bot( ii )%x, Bot( ii )%HS%alphaR, Bot( ii )%HS%betaR, Bot( ii )%HS%rho, & Bot( ii )%HS%alphaI, Bot( ii )%HS%betaI ! LP: Same change as above IF ( ii < Number_to_Echo .OR. ii == NbtyPts - 1 ) THEN ! echo some values WRITE( PRTFile, FMT="( F10.2, F10.2, 3X, 2F10.2, 3X, F6.2, 3X, 2F10.4 )" ) & Bot( ii )%x, Bot( ii )%HS%alphaR, Bot( ii )%HS%betaR, Bot( ii )%HS%rho, & Bot( ii )%HS%alphaI, Bot( ii )%HS%betaI END IF CASE DEFAULT CALL ERROUT( 'ReadBTY', 'Unknown option for selecting bathymetry option' ) END SELECT IF ( Bot( ii )%x( 2 ) > DepthB ) THEN CALL ERROUT( 'BELLHOP:ReadBTY', 'Bathymetry drops below lowest point in the sound speed profile' ) END IF END DO btypt CLOSE( BTYFile ) Bot( : )%x( 1 ) = 1000.0 * Bot( : )%x( 1 ) ! Convert ranges in km to m CASE DEFAULT ! no bathymetry given, use SSP depth for flat bottom ALLOCATE( Bot( 2 ), Stat = IAllocStat ) IF ( IAllocStat /= 0 ) CALL ERROUT( 'BELLHOP', 'Insufficient memory for bathymetry data' ) Bot( 1 )%x = [ -sqrt( huge( Bot( 1 )%x( 1 ) ) ) / 1.0d5, DepthB ] Bot( 2 )%x = [ sqrt( huge( Bot( 1 )%x( 1 ) ) ) / 1.0d5, DepthB ] END SELECT CALL ComputeBdryTangentNormal( Bot, 'Bot' ) IF ( .NOT. monotonic( Bot%x( 1 ), NBtyPts ) ) THEN CALL ERROUT( 'BELLHOP:ReadBTY', 'Bathymetry ranges are not monotonically increasing' ) END IF END SUBROUTINE ReadBTY