Reads in the top altimetry
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
character(len=80), | intent(in) | :: | FileRoot | |||
character(len=1), | intent(in) | :: | TopATI | |||
real(kind=8), | intent(in) | :: | DepthT | |||
integer, | intent(in) | :: | PRTFile |
SUBROUTINE ReadATI( FileRoot, TopATI, DepthT, PRTFile ) !! Reads in the top altimetry INTEGER, INTENT( IN ) :: PRTFile CHARACTER (LEN= 1), INTENT( IN ) :: TopATI REAL (KIND=8), INTENT( IN ) :: DepthT ! REAL (KIND=8), ALLOCATABLE :: phi( : ) ! LP: Removed as seems to have been moved to ComputeBdryTangentNormal CHARACTER (LEN=80), INTENT( IN ) :: FileRoot SELECT CASE ( TopATI ) CASE ( '~', '*' ) WRITE( PRTFile, * ) '__________________________________________________________________________' WRITE( PRTFile, * ) WRITE( PRTFile, * ) 'Using top-altimetry file' OPEN( UNIT = ATIFile, FILE = TRIM( FileRoot ) // '.ati', STATUS = 'OLD', IOSTAT = IOStat, ACTION = 'READ' ) IF ( IOsTAT /= 0 ) THEN WRITE( PRTFile, * ) 'ATIFile = ', TRIM( FileRoot ) // '.ati' CALL ERROUT( 'ReadATI', 'Unable to open altimetry file' ) END IF READ( ATIFile, * ) atiType AltiType: SELECT CASE ( atiType( 1 : 1 ) ) CASE ( 'C' ) WRITE( PRTFile, * ) 'Curvilinear Interpolation' CASE ( 'L' ) WRITE( PRTFile, * ) 'Piecewise linear interpolation' CASE DEFAULT CALL ERROUT( 'ReadATI', 'Unknown option for selecting altimetry interpolation' ) END SELECT AltiType READ( ATIFile, * ) NatiPts WRITE( PRTFile, * ) 'Number of altimetry points = ', NatiPts NatiPts = NatiPts + 2 ! we'll be extending the altimetry to infinity to the left and right ALLOCATE( Top( NatiPts ), Stat = IAllocStat ) ! ALLOCATE( phi( NatiPts ), Stat = IAllocStat ) ! LP: Removed as seems to have been moved to ComputeBdryTangentNormal IF ( IAllocStat /= 0 ) & CALL ERROUT( 'BELLHOP:ReadATI', 'Insufficient memory for altimetry data: reduce # ati points' ) WRITE( PRTFile, * ) WRITE( PRTFile, * ) ' Range (km) Depth (m)' atiPt: DO ii = 2, NatiPts - 1 SELECT CASE ( atiType( 2 : 2 ) ) CASE ( 'S', '' ) READ( ATIFile, * ) Top( ii )%x ! LP: This condition was previously ii == NatiPts, ! which will never be satisfied due to the loop bounds IF ( ii < Number_to_Echo .OR. ii == NatiPts - 1 ) THEN ! echo some values WRITE( PRTFile, FMT = "(2G11.3)" ) Top( ii )%x END IF CASE ( 'L' ) READ( ATIFile, * ) Top( ii )%x, Top( ii )%HS%alphaR, Top( ii )%HS%betaR, Top( ii )%HS%rho, & Top( ii )%HS%alphaI, Top( ii )%HS%betaI ! LP: Same change as above IF ( ii < Number_to_Echo .OR. ii == NatiPts - 1 ) THEN ! echo some values WRITE( PRTFile, FMT = "(7G11.3)" ) Top( ii )%x, Top( ii )%HS%alphaR, Top( ii )%HS%betaR, Top( ii )%HS%rho, & Top( ii )%HS%alphaI, Top( ii )%HS%betaI END IF CASE DEFAULT CALL ERROUT( 'ReadATI', 'Unknown option for selecting altimetry option' ) END SELECT IF ( Top( ii )%x( 2 ) < DepthT ) THEN CALL ERROUT( 'BELLHOP:ReadATI', 'Altimetry rises above highest point in the sound speed profile' ) END IF END DO atiPt CLOSE( ATIFile ) Top( : )%x( 1 ) = 1000.0 * Top( : )%x( 1 ) ! Convert ranges in km to m CASE DEFAULT ! no altimetry given, use SSP depth for flat top ALLOCATE( Top( 2 ), Stat = IAllocStat ) IF ( IAllocStat /= 0 ) CALL ERROUT( 'BELLHOP', 'Insufficient memory for altimetry data' ) Top( 1 )%x = [ -sqrt( huge( Top( 1 )%x( 1 ) ) ) / 1.0d5, DepthT ] Top( 2 )%x = [ sqrt( huge( Top( 1 )%x( 1 ) ) ) / 1.0d5, DepthT ] END SELECT CALL ComputeBdryTangentNormal( Top, 'Top' ) IF ( .NOT. monotonic( Top%x( 1 ), NAtiPts ) ) THEN CALL ERROUT( 'BELLHOP:ReadATI', 'Altimetry ranges are not monotonically increasing' ) END IF END SUBROUTINE ReadATI