forked from emolch/kiwi
-
Notifications
You must be signed in to change notification settings - Fork 0
/
Copy pathseismogram_io.f90
250 lines (196 loc) · 8.98 KB
/
seismogram_io.f90
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
!
! Copyright 2011 Sebastian Heimann
!
! Licensed under the Apache License, Version 2.0 (the "License");
! you may not use this file except in compliance with the License.
! You may obtain a copy of the License at
!
! http://www.apache.org/licenses/LICENSE-2.0
!
! Unless required by applicable law or agreed to in writing, software
! distributed under the License is distributed on an "AS IS" BASIS,
! WITHOUT WARRANTIES OR CONDITIONS OF ANY KIND, either express or implied.
! See the License for the specific language governing permissions and
! limitations under the License.
!
module seismogram_io
! this module provides methods for writing seismograms
! at the moment ascii table output and sac output are intended,
! but other formats may be added
use better_varying_string
use read_table
use unit
use util
implicit none
private
public :: writeseismogram, readseismogram
interface writeseismogram
module procedure writeseismogram_vs
module procedure writeseismogram_c
end interface
interface readseismogram
module procedure readseismogram_c
end interface
contains
subroutine writeseismogram_vs( filename, fileformat, seismogram, toffset, deltat, &
network, station, location, channel, nerr )
type(varying_string), intent(in) :: filename, fileformat
real, dimension(:), intent(in) :: seismogram
real(kind=8), intent(in) :: toffset
real, intent(in) :: deltat
type(varying_string), intent(in) :: network, station, location, channel
integer, intent(out) :: nerr
call writeseismogram_c( char(filename), char(fileformat), seismogram, toffset, deltat, &
char(network), char(station), char(location), char(channel), nerr )
end subroutine
subroutine writeseismogram_c( filename, fileformat, seismogram, toffset, deltat, &
network, station, location, channel, nerr )
character(len=*), intent(in) :: filename, fileformat
real, dimension(:), intent(in) :: seismogram
real(kind=8), intent(in) :: toffset
real, intent(in) :: deltat
integer, intent(out) :: nerr
character(len=*), intent(in) :: network, station, location, channel
type(varying_string) :: fileformat_
! dump 1-component seismogram to file
integer :: iunit, i, nlen
real(kind=8) :: ddeltat
real, dimension(size(seismogram)) :: seismogram_copy
character(len=len(filename)+1) :: filename_cstr
character(len=len(network)+1) :: network_cstr
character(len=len(station)+1) :: station_cstr
character(len=len(location)+1) :: location_cstr
character(len=len(channel)+1) :: channel_cstr
! look at the filename extension what kind of output is wanted
nerr = 0
fileformat_ = trim(fileformat)
if (fileformat == '*') then
fileformat_ = 'table'
if ( len(filename) >= 4 ) then
if (filename(len(filename)-3:len(filename)) == '.sac') &
fileformat_ = 'sac'
end if
if ( len(filename) >= 6 ) then
if (filename(len(filename)-5:len(filename)) == '.mseed') &
fileformat_ = 'mseed'
end if
end if
if (fileformat_ == 'sac') then
seismogram_copy(:) = seismogram(:)
filename_cstr = filename//char(0)
nlen = size(seismogram_copy)
call wsac1( filename_cstr, seismogram_copy, nlen, &
toffset, deltat, nerr )
! call die("sac output currently disabled, because there is no 64bit libsacio.")
end if
if (fileformat_ == 'mseed') then
seismogram_copy(:) = seismogram(:)
ddeltat = deltat
filename_cstr = trim(filename)//char(0)
network_cstr = trim(network)//char(0)
station_cstr = trim(station)//char(0)
location_cstr = trim(location)//char(0)
channel_cstr = trim(channel)//char(0)
nlen = size(seismogram_copy)
call writemseed( filename_cstr, seismogram_copy, nlen, &
toffset, ddeltat, network_cstr, station_cstr, location_cstr, channel_cstr, nerr )
end if
if (fileformat_ == 'table') then
call claim_unit( iunit )
open( unit=iunit, file=filename, status='unknown', iostat=nerr )
if ( nerr /= 0 ) then
call release_unit( iunit )
return
end if
do i=1,size(seismogram)
write (iunit,*) toffset+(i-1)*deltat, seismogram(i)
end do
close( iunit )
call release_unit( iunit )
end if
end subroutine
subroutine readseismogram_c( filename, fileformat, seismogram, toffset, deltat, nerr )
character(len=*), intent(in) :: filename, fileformat
real, dimension(:), allocatable, intent(inout) :: seismogram
real(kind=8), intent(out) :: toffset
real, intent(out) :: deltat
integer, intent(out) :: nerr
type(varying_string) :: vsfn
integer :: nlen, maxlen
real, dimension(:),allocatable :: scratch
real, dimension(:,:),allocatable :: tablebuf
type(varying_string) :: fileformat_
character(len=len(filename)+1) :: filename_cstr
real(kind=8) :: ddeltat
! this is not a very brilliant way to determine the file type, but
! anyway, look at the filename extension if this is sac data...
fileformat_ = fileformat
if (fileformat == '*') then
fileformat_ = 'table'
if ( len(filename) >= 4 ) then
if (filename(len(filename)-3:len(filename)) == '.sac') &
fileformat_ = 'sac'
end if
if ( len(filename) >= 6 ) then
if (filename(len(filename)-5:len(filename)) == '.mseed') &
fileformat_ = 'mseed'
end if
end if
if (fileformat_ == 'sac') then
! read data into buffer 'scratch'...
! couldn't find a way to get the correct data length, without reading the
! full thing... so this stupid iterative attempt is used.
! - getnhv('NPTS',...) always returned maxlen when there was more data in file
! - nerr is not set to -803 when there is more data, opposing to what
! the manual tells
! (at least with my copy of libsacio)
! what a crap!
nerr = -1
maxlen = 1024
do while (nerr<0)
if (allocated(scratch)) deallocate(scratch)
allocate( scratch(maxlen) )
filename_cstr = filename//char(0)
call rsac1(filename_cstr, scratch, nlen, toffset, deltat, maxlen, nerr)
! call die("sac input currently disabled, because there is no 64bit libsacio")
if (nerr > 0) then
call error( "rsac1 returned an error" )
if (allocated(scratch)) deallocate(scratch)
return
end if
if (maxlen == nlen) nerr = -1
maxlen = maxlen*2
end do
if (allocated(seismogram)) deallocate(seismogram)
allocate( seismogram(nlen) )
seismogram(:) = scratch(1:nlen)
deallocate(scratch)
end if
if (fileformat_ == 'mseed') then
filename_cstr = filename//char(0)
call readmseed1( filename_cstr, nlen, toffset, ddeltat, nerr )
if (nerr .ne. 0) then
call error( "readmseed1 returned an error" )
return
end if
if (allocated(seismogram)) deallocate(seismogram)
allocate( seismogram(nlen) )
call readmseed2( seismogram )
deltat = real(ddeltat)
end if
if (fileformat_ == 'table') then
vsfn = filename
call readtable_file( tablebuf, vsfn, min_cols=2, min_rows=2, nerr=nerr )
if (nerr /= 0) then
return
end if
nlen = size(tablebuf,2)
if (allocated(seismogram)) deallocate(seismogram)
allocate( seismogram(nlen) )
seismogram(:) = tablebuf(2,:)
toffset = tablebuf(1,1)
deltat = (tablebuf(1,nlen)-tablebuf(1,1))/(nlen-1)
deallocate( tablebuf )
end if
end subroutine
end module