aboutsummaryrefslogtreecommitdiff
path: root/src/EHFinder_ReadData.F90
blob: 7215d0a47c9cae8d92292a76a85952fcfa49edc4 (plain)
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
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
! Read metric, lapse, shift and conformal factor from files.
! $Header$

#include "cctk.h"
#include "cctk_Parameters.h"
#include "cctk_Arguments.h"
#include "cctk_Functions.h"

!This routine reads in the metric (from ADMBase).
subroutine EHFinder_Read_Metric(CCTK_ARGUMENTS)

  use EHFinder_mod

  implicit none

  DECLARE_CCTK_PARAMETERS
  DECLARE_CCTK_ARGUMENTS
  DECLARE_CCTK_FUNCTIONS

  character(len=512) :: in_files, in_vars
  character(len=64), dimension(6)  :: var_names
  character(len=10) :: iteration_string
  CCTK_INT :: i, nc, ntot, res

! Figure out which iteration number to read, based on the parameters
! last_iteration_number (the last iteration in the numerical evolution
! producing the metric) and saved_iteration_every (how often was the metric
! saved) and the current iteration and save it in a string variable.
  i = min ( last_iteration_number - saved_iteration_every * cctk_iteration, &
            last_iteration_number )
  write(iteration_string,'(i10)') i

! Trim the string variable.
  iteration_string = adjustl(iteration_string)
  nc = len_trim(iteration_string)

! Generate a string with the filenames. Note at present the requirement
! therefore is that all iterations of one variable is written in the
! same file.

  in_files = 'gxx gxy gxz gyy gyz gzz'
  if ( CCTK_EQUALS(file_type,'sep_time_files') ) then
    call create_filenames ( in_files, i )
  end if

! Fill in the string array, used to requesting the metric components at
! the specified iteration number.
  var_names(1) = 'admbase::gxx{cctk_iteration='//iteration_string(1:nc)//'}'
  var_names(2) = 'admbase::gxy{cctk_iteration='//iteration_string(1:nc)//'}'
  var_names(3) = 'admbase::gxz{cctk_iteration='//iteration_string(1:nc)//'}'
  var_names(4) = 'admbase::gyy{cctk_iteration='//iteration_string(1:nc)//'}'
  var_names(5) = 'admbase::gyz{cctk_iteration='//iteration_string(1:nc)//'}'
  var_names(6) = 'admbase::gzz{cctk_iteration='//iteration_string(1:nc)//'}'

! merge all the variable names into a single string.
  in_vars = ' '
  ntot = 0
  do i = 1, 6
    nc = len_trim(var_names(i))
    in_vars(ntot+1:ntot+1+nc+1) = var_names(i)(1:nc+1)
    ntot = ntot + nc + 1
  end do
   
! Call the routine that actualle does the file access. Note failures are
! just silently ignored.
  call IOUtil_RecoverVarsFromDatafiles ( res, cctkGH, in_files, in_vars )

end subroutine EHFinder_Read_Metric


! This routine reads in the lapse (from ADMBase).
subroutine EHFinder_Read_Lapse(CCTK_ARGUMENTS)

  use EHFinder_mod

  implicit none

  DECLARE_CCTK_PARAMETERS
  DECLARE_CCTK_ARGUMENTS
  DECLARE_CCTK_FUNCTIONS

  character(len=128) :: in_files, in_vars
  character(len=10) :: iteration_string
  CCTK_INT :: i, nc, res

! Figure out which iteration number to read, based on the parameters
! last_iteration_number (the last iteration in the numerical evolution
! producing the metric) and saved_iteration_every (how often was the lapse
! saved) and the current iteration and save it in a string variable.
  i = min ( last_iteration_number - saved_iteration_every * cctk_iteration, &
            last_iteration_number )
  write(iteration_string,'(i10)') i 

! Trim the string variable.
  iteration_string = adjustl(iteration_string)
  nc = len_trim(iteration_string)

! Generate a string with the filename. Note at present the requirement
! therefore is that all iterations of one variable is written in the
! same file.
  in_files = 'alp'
  if ( CCTK_EQUALS(file_type,'sep_time_files') ) then
    call create_filenames ( in_files, i )
  end if

! Fill in the string used to requesting the lapse at the specified
! iteration number.
  in_vars = 'admbase::alp{cctk_iteration='//iteration_string(1:nc)//'}'

! Call the routine that actualle does the file access. Note failures are
! just silently ignored.
  call IOUtil_RecoverVarsFromDatafiles ( res, cctkGH, in_files, in_vars )

end subroutine EHFinder_Read_Lapse


! This routine reads in all the shift components (from ADMBAse).
subroutine EHFinder_Read_Shift(CCTK_ARGUMENTS)

  use EHFinder_mod

  implicit none

  DECLARE_CCTK_PARAMETERS
  DECLARE_CCTK_ARGUMENTS
  DECLARE_CCTK_FUNCTIONS

  character(len=256) :: in_files, in_vars
  character(len=64), dimension(3)  :: var_names
  character(len=10) :: iteration_string
  CCTK_INT :: i, nc, ntot, res

! Figure out which iteration number to read, based on the parameters
! last_iteration_number (the last iteration in the numerical evolution
! producing the metric) and saved_iteration_every (how often was the metric
! saved) and the current iteration and save it in a string variable.
  i = min ( last_iteration_number - saved_iteration_every * cctk_iteration, &
            last_iteration_number )
  write(iteration_string,'(i10)') i 

! Trim the string variable.
  iteration_string = adjustl(iteration_string)
  nc = len_trim(iteration_string)

! Generate a string with the filenames. Note at present the requirement
! therefore is that all iterations of one variable is written in the
! same file.
  in_files = 'betax betay betaz'
  if ( CCTK_EQUALS(file_type,'sep_time_files') ) then
    call create_filenames ( in_files, i )
  end if

! Fill in the string array, used to requesting the shift components at
! the specified iteration number.
  var_names(1) = 'admbase::betax{cctk_iteration='//iteration_string(1:nc)//'}'
  var_names(2) = 'admbase::betay{cctk_iteration='//iteration_string(1:nc)//'}'
  var_names(3) = 'admbase::betaz{cctk_iteration='//iteration_string(1:nc)//'}'

! merge all the variable names into a single string.
  in_vars = ' '
  ntot = 0
  do i = 1, 3
    nc = len_trim(var_names(i))
    in_vars(ntot+1:ntot+1+nc+1) = var_names(i)(1:nc+1)
    ntot = ntot + nc + 1
  end do
 
! Call the routine that actualle does the file access. Note failures are
! just silently ignored.
  call IOUtil_RecoverVarsFromDatafiles ( res, cctkGH, in_files, in_vars )

end subroutine EHFinder_Read_Shift


! This routine reads in the conformal factor (from StaticConformal).
subroutine EHFinder_Read_Conformal(CCTK_ARGUMENTS)

  use EHFinder_mod

  implicit none

  DECLARE_CCTK_PARAMETERS
  DECLARE_CCTK_ARGUMENTS
  DECLARE_CCTK_FUNCTIONS

  character(len=128) :: in_files, in_vars
  character(len=10) :: iteration_string
  CCTK_INT :: i, nc, res

! Figure out which iteration number to read, based on the parameters
! last_iteration_number (the last iteration in the numerical evolution
! producing the metric) and saved_iteration_every (how often was the metric
! saved) and the current iteration and save it in a string variable. Note
! that for this variable the default option is to only to read it in once at
! the first timestep.
  if ( read_conformal_factor_once .gt. 0 ) then
    write(iteration_string,'(i10)') 0
  else
    i = min ( last_iteration_number - saved_iteration_every * cctk_iteration, &
              last_iteration_number )
    write(iteration_string,'(i10)') i 
  end if

! Trim the string variable.
  iteration_string = adjustl(iteration_string)
  nc = len_trim(iteration_string)

! Generate a string with the filenames. Note at present the requirement
! therefore is that all iterations of one variable is written in the
! same file.
  in_files = 'psi'
  if ( CCTK_EQUALS(file_type,'sep_time_files') ) then
    call create_filenames ( in_files, i )
  end if

! Fill in the string array, used to requesting the conformal factor at
! the specified iteration number.
  in_vars = 'staticconformal::psi{cctk_iteration='//iteration_string(1:nc)//'}'

! Call the routine that actualle does the file access. Note failures are
! just silently ignored.
  call IOUtil_RecoverVarsFromDatafiles ( res, cctkGH, in_files, in_vars )

end subroutine EHFinder_Read_Conformal


! This routine reads in the excision mask (from SpaceMask).
subroutine EHFinder_Read_Mask(CCTK_ARGUMENTS)

  use EHFinder_mod

  implicit none

  DECLARE_CCTK_PARAMETERS
  DECLARE_CCTK_ARGUMENTS
  DECLARE_CCTK_FUNCTIONS

  character(len=128) :: in_files, in_vars
  character(len=10) :: iteration_string
  CCTK_INT :: i, nc, res

! Figure out which iteration number to read, based on the parameters
! last_iteration_number (the last iteration in the numerical evolution
! producing the metric) and saved_iteration_every (how often was the metric
! saved) and the current iteration and save it in a string variable.
  i = min ( last_iteration_number - saved_iteration_every * cctk_iteration, &
            last_iteration_number )
  write(iteration_string,'(i10)') i 

! Trim the string variable.
  iteration_string = adjustl(iteration_string)
  nc = len_trim(iteration_string)

! Generate a string with the filename. Note at present the requirement
! therefore is that all iterations of one variable is written in the
! same file.
  in_files = 'emask'
  if ( CCTK_EQUALS(file_type,'sep_time_files') ) then
    call create_filenames ( in_files, i )
  end if

! Fill in the string array, used to requesting the old style mask at
! the specified iteration number.
  in_vars = 'spacemask::emask{cctk_iteration='//iteration_string(1:nc)//'}'

! Call the routine that actualle does the file access. Note failures are
! just silently ignored.
  call IOUtil_RecoverVarsFromDatafiles ( res, cctkGH, in_files, in_vars )

end subroutine EHFinder_Read_Mask

subroutine create_filenames ( in_files, i )

  implicit none

  character(len=*), intent(INOUT) :: in_files
  CCTK_INT, intent(IN) :: i

  character(len=len(in_files)) :: new_files
  character(len=len(in_files)) :: tmp1, tmp2, istr
  character(len=len_trim(in_files)) :: old_files
  character(len=7), parameter :: template1 = '_000000'
  character(len=7) :: template2
  logical :: last
  CCTK_INT :: j, nc, lnew_file, ltmp, ltmp2

  last = .false.
  old_files = in_files
  lnew_file = 0
  ltmp = len(old_files)
  ltmp2 = 0
  tmp2 = ''

  write(istr,'(i10)') i
  istr = adjustl(istr)
  nc = len_trim(istr)
  template2 = template1
  template2(8-nc:7) = istr(1:nc)

  j = scan (old_files(1:ltmp), ' ')
  name_loop: do
    if ( j == 0 ) last = .true.
!    print*,'last = ', last
    tmp1 = ''
    if ( .not. last  ) then
      tmp1 = old_files(1:j-1)
      old_files(1:ltmp-j) = old_files(j+1:ltmp)
    else
      tmp1 = old_files(1:ltmp)
      old_files = ''
    end if
!    print*,'tmp1 = ',tmp1(1:j-1)
!    print*,'old_files = ',old_files(1:ltmp)
    ltmp = ltmp - j
!    print*,'ltmp = ',ltmp
    tmp2 = adjustl(trim(tmp2)//' '//trim(tmp1)//template2)
!    print*,'tmp2 = ',trim(tmp2)
!    print* 
    if ( last ) exit name_loop
    j = scan (old_files(1:ltmp), ' ') 
  end do name_loop 

  in_files = tmp2
!  pause
end subroutine create_filenames