uniform_velocity.F90 coverage: 57.14 %func 14.93 %block
1) module Uniform_Velocity_module
2)
3) use PFLOTRAN_Constants_module
4)
5) implicit none
6)
7) private
8)
9) #include "petsc/finclude/petscsys.h"
10)
11) PetscInt, parameter :: NULL = 0
12) PetscInt, parameter :: STEP = 1
13) PetscInt, parameter :: LINEAR = 2
14)
15) type, public :: uniform_velocity_dataset_type
16) PetscInt :: rank
17) PetscBool :: is_transient
18) PetscBool :: is_cyclic
19) PetscInt :: interpolation_method
20) PetscReal, pointer :: times(:)
21) PetscReal, pointer :: values(:,:)
22) PetscReal, pointer :: cur_value(:)
23) PetscInt :: cur_time_index
24) PetscInt :: max_time_index
25) PetscReal :: time_shift
26) end type uniform_velocity_dataset_type
27)
28)
29) public :: UniformVelocityDatasetCreate, &
30) UniformVelocityDatasetDestroy, &
31) UniformVelocityDatasetRead, &
32) UniformVelocityDatasetUpdate, &
33) UniformVelocityDatasetVerify
34)
35) contains
36)
37) ! ************************************************************************** !
38)
39) function UniformVelocityDatasetCreate()
40) !
41) ! Creates a velocity data set
42) !
43) ! Author: Glenn Hammond
44) ! Date: 06/02/09
45) !
46)
47) implicit none
48)
49) type(uniform_velocity_dataset_type), pointer :: UniformVelocityDatasetCreate
50)
51) type(uniform_velocity_dataset_type), pointer :: dataset
52)
53) allocate(dataset)
54) nullify(dataset%times)
55) nullify(dataset%values)
56) nullify(dataset%cur_value)
57) dataset%cur_time_index = 0
58) dataset%max_time_index = 0
59) dataset%rank = 0
60) dataset%is_cyclic = PETSC_FALSE
61) dataset%is_transient = PETSC_FALSE
62) dataset%interpolation_method = NULL
63) dataset%time_shift = 0.d0
64)
65) UniformVelocityDatasetCreate => dataset
66)
67) end function UniformVelocityDatasetCreate
68)
69) ! ************************************************************************** !
70)
71) subroutine UniformVelocityDatasetRead(dataset,input,option)
72) !
73) ! Reads a velocity data set from the input file
74) !
75) ! Author: Glenn Hammond
76) ! Date: 06/02/09
77) !
78)
79) use Option_module
80) use Input_Aux_module
81) use String_module
82) use Logging_module
83) use Units_module
84)
85) implicit none
86)
87) type(uniform_velocity_dataset_type) :: dataset
88) type(input_type), pointer :: input
89) type(option_type) :: option
90)
91) character(len=MAXSTRINGLENGTH) :: string
92) character(len=MAXWORDLENGTH) :: word, units, internal_units
93) PetscReal :: units_conversion
94)
95) PetscErrorCode :: ierr
96)
97) ! call PetscLogEventBegin(logging%event_flow_condition_read,ierr)
98)
99) dataset%rank = 3
100) dataset%interpolation_method = STEP
101) dataset%is_cyclic = PETSC_FALSE
102)
103) units = ''
104)
105) ! read the velocity data set
106) input%ierr = 0
107) do
108)
109) call InputReadPflotranString(input,option)
110) call InputReadStringErrorMsg(input,option,'VELOCITY_DATASET')
111)
112) if (InputCheckExit(input,option)) exit
113)
114) call InputReadWord(input,option,word,PETSC_TRUE)
115) call InputErrorMsg(input,option,'keyword','VELOCITY_DATASET')
116)
117) select case(trim(word))
118)
119) case('UNITS') ! read default units for condition arguments
120) call InputReadWord(input,option,units,PETSC_TRUE)
121) call InputErrorMsg(input,option,'UNITS','VELOCITY_DATASET')
122) case('CYCLIC')
123) dataset%is_cyclic = PETSC_TRUE
124) case('INTERPOLATION')
125) call InputReadWord(input,option,word,PETSC_TRUE)
126) call InputErrorMsg(input,option,'INTERPOLATION','VELOCITY_DATASET')
127) call StringToLower(word)
128) select case(word)
129) case('step')
130) dataset%interpolation_method = STEP
131) case('linear')
132) dataset%interpolation_method = LINEAR
133) end select
134) case('VELOCITY')
135) call UniVelocityDatasetReadValues(input,option,word,string, &
136) dataset,units)
137) case default
138) call InputKeywordUnrecognized(word,'VELOCITY_DATASET',option)
139) end select
140)
141) enddo
142)
143) if (len_trim(units) > 1) then
144) internal_units = 'meter/sec'
145) units_conversion = UnitsConvertToInternal(units,internal_units,option)
146) dataset%values = dataset%values * units_conversion
147) word = units(index(units,'/')+1:)
148) internal_units = 'sec'
149) units_conversion = UnitsConvertToInternal(word,internal_units,option)
150) dataset%times = dataset%times * units_conversion
151) endif
152)
153) call UniformVelocityDatasetVerify(option, dataset)
154)
155) ! call PetscLogEventEnd(logging%event_flow_condition_read,ierr)
156)
157) end subroutine UniformVelocityDatasetRead
158)
159) ! ************************************************************************** !
160)
161) subroutine UniVelocityDatasetReadValues(input,option,keyword,string,dataset, &
162) units)
163) !
164) ! Read the value(s) of a velocity data set
165) !
166) ! Author: Glenn Hammond
167) ! Date: 06/02/09
168) !
169)
170) use Input_Aux_module
171) use String_module
172) use Option_module
173) use Logging_module
174)
175) implicit none
176)
177) type(input_type), target :: input
178) type(option_type) :: option
179) character(len=MAXSTRINGLENGTH) :: string
180) character(len=MAXWORDLENGTH) :: keyword
181) type(uniform_velocity_dataset_type) :: dataset
182) character(len=MAXWORDLENGTH) :: units
183)
184) type(input_type), pointer :: input2
185) character(len=MAXSTRINGLENGTH) :: string2, filename
186) character(len=MAXWORDLENGTH) :: word
187) character(len=MAXSTRINGLENGTH) :: error_string
188) PetscBool :: read_from_file
189) PetscBool :: read_multiple_values
190) PetscInt :: irank
191) PetscErrorCode :: ierr
192)
193) ! call PetscLogEventBegin(logging%event_flow_condition_read_values,ierr)
194)
195) read_from_file = PETSC_FALSE
196) read_multiple_values = PETSC_TRUE
197)
198) input%ierr = 0
199) string2 = trim(input%buf)
200) call InputReadWord(input,option,word,PETSC_TRUE)
201) if (input%ierr == 0) then
202) call StringToLower(word)
203) if (StringCompare(word,'file',FOUR_INTEGER)) then
204) call InputReadNChars(input,option,filename,MAXSTRINGLENGTH,PETSC_TRUE)
205) input%err_buf = trim(keyword) // ' FILE'
206) input%err_buf2 = 'VELOCITY_DATASET'
207) call InputErrorMsg(input,option)
208) read_from_file = PETSC_TRUE
209) else if (StringCompare(word,'list',FOUR_INTEGER)) then
210) read_multiple_values = PETSC_TRUE
211) else
212) read_multiple_values = PETSC_FALSE
213) endif
214) endif
215)
216) if (read_from_file .or. read_multiple_values) then
217) if (read_from_file) then
218) input2 => InputCreate(IUNIT_TEMP,filename,option)
219) else
220) input2 => input
221) endif
222) call UniVelocityDatasetReadFromFile(input2,dataset,option)
223) if (read_from_file) call InputDestroy(input2)
224) else
225) input%buf = trim(string2)
226) allocate(dataset%times(1))
227) dataset%times = 0.d0
228) allocate(dataset%values(dataset%rank,1))
229) dataset%values = 0.d0
230) do irank=1,dataset%rank
231) call InputReadDouble(input,option,dataset%values(irank,1))
232) write(input%err_buf,'(a,i2)') trim(keyword) // &
233) ' dataset_values, irank = ', irank
234) input%err_buf2 = 'VELOCITY_DATASET'
235) call InputErrorMsg(input,option)
236) enddo
237) call InputReadWord(input,option,word,PETSC_TRUE)
238) if (InputError(input)) then
239) word = trim(keyword) // ' UNITS'
240) call InputDefaultMsg(input,option,word)
241) else
242) units = trim(word)
243) endif
244) endif
245)
246) ! call PetscLogEventEnd(logging%event_flow_condition_read_values,ierr)
247)
248) end subroutine UniVelocityDatasetReadValues
249)
250) ! ************************************************************************** !
251)
252) subroutine UniVelocityDatasetReadFromFile(input,dataset,option)
253) !
254) ! Read values from a external file
255) !
256) ! Author: Glenn Hammond
257) ! Date: 10/31/07
258) !
259)
260) use Input_Aux_module
261) use String_module
262) use Utility_module
263) use Option_module
264)
265) implicit none
266)
267) type(option_type) :: option
268) type(uniform_velocity_dataset_type) :: dataset
269) character(len=MAXSTRINGLENGTH) :: filename
270)
271) character(len=MAXSTRINGLENGTH) :: string
272) PetscReal, pointer :: temp_times(:), temp_array1(:), temp_array2(:), &
273) temp_array3(:)
274) PetscReal :: temp_time
275) PetscInt :: max_size
276) PetscInt :: temp_max_size
277) PetscInt :: count, i, status
278) type(input_type), pointer :: input
279)
280) max_size = 1000
281) allocate(temp_times(max_size))
282) allocate(temp_array1(max_size))
283) temp_times = 0.d0
284) temp_array1 = 0.d0
285)
286) if (dataset%rank > 1) then
287) allocate(temp_array2(max_size))
288) temp_array2 = 0.d0
289) endif
290)
291) if (dataset%rank > 2) then
292) allocate(temp_array3(max_size))
293) temp_array3 = 0.d0
294) endif
295)
296) count = 0
297) do
298) call InputReadPflotranString(input,option)
299) if (InputError(input)) exit
300) if (InputCheckExit(input,option)) exit
301) count = count + 1
302) call InputReadDouble(input,option,temp_times(count))
303) call InputErrorMsg(input,option,'time','VELOCITY_DATASET FILE')
304) call InputReadDouble(input,option,temp_array1(count))
305) call InputErrorMsg(input,option,'array1','VELOCITY_DATASET FILE')
306) if (dataset%rank > 1) then
307) call InputReadDouble(input,option,temp_array2(count))
308) call InputErrorMsg(input,option,'array2','VELOCITY_DATASET FILE')
309) endif
310) if (dataset%rank > 2) then
311) call InputReadDouble(input,option,temp_array3(count))
312) call InputErrorMsg(input,option,'array3','VELOCITY_DATASET FILE')
313) endif
314) if (count+1 > max_size) then
315) temp_max_size = max_size
316) call reallocateRealArray(temp_times,max_size)
317) ! careful. reallocateRealArray double max_size every time.
318) i = temp_max_size
319) call reallocateRealArray(temp_array1,i)
320) if (dataset%rank > 1) then
321) i = temp_max_size
322) call reallocateRealArray(temp_array2,i)
323) endif
324) if (dataset%rank > 2) then
325) i = temp_max_size
326) call reallocateRealArray(temp_array3,i)
327) endif
328) endif
329) enddo
330)
331) if (associated(dataset%times)) deallocate(dataset%times)
332) nullify(dataset%times)
333) if (associated(dataset%values)) deallocate(dataset%values)
334) nullify(dataset%values)
335)
336) allocate(dataset%times(count))
337) allocate(dataset%values(dataset%rank,count))
338)
339) dataset%times(1:count) = temp_times(1:count)
340) dataset%values(1,1:count) = temp_array1(1:count)
341) if (dataset%rank > 1) dataset%values(2,1:count) = temp_array2(1:count)
342) if (dataset%rank > 2) dataset%values(3,1:count) = temp_array3(1:count)
343)
344) deallocate(temp_times)
345) deallocate(temp_array1)
346) if (dataset%rank > 1) deallocate(temp_array2)
347) if (dataset%rank > 2) deallocate(temp_array3)
348)
349) end subroutine UniVelocityDatasetReadFromFile
350)
351) ! ************************************************************************** !
352)
353) subroutine UniformVelocityDatasetVerify(option, dataset)
354) !
355) ! Verifies the data in a dataset
356) !
357) ! Author: Glenn Hammond
358) ! Date: 02/04/08
359) !
360) use Option_module
361)
362) implicit none
363)
364) type(option_type) :: option
365) type(uniform_velocity_dataset_type) :: dataset
366)
367) PetscInt :: array_size
368) character(len=MAXWORDLENGTH) :: size1, size2
369)
370) dataset%max_time_index = 1
371) if (.not.associated(dataset%times)) then
372) array_size = 1
373) allocate(dataset%times(array_size))
374) dataset%times = 0.d0
375) endif
376) if (.not.associated(dataset%values)) then
377) ! allocate to size 1
378) array_size = 1
379) allocate(dataset%values(1:dataset%rank,array_size))
380) dataset%values(1:dataset%rank,1:array_size) = 0.d0
381) endif
382) dataset%max_time_index = size(dataset%times,1)
383) if (dataset%max_time_index > 1) then
384) if (size(dataset%values,2) /= & ! size of 2nd rank
385) dataset%max_time_index) then
386) write(size1,*) size(dataset%times,1)
387) write(size2,*) size(dataset%values,2)
388) option%io_buffer = 'times/values ('//trim(size1)//'/'//trim(size2) // &
389) ') array size mismatch in velocity dataaset'
390) call printErrMsg(option)
391) endif
392) dataset%is_transient = PETSC_TRUE
393) else
394) dataset%is_transient = PETSC_FALSE
395) endif
396) dataset%cur_time_index = 1
397)
398) if (associated(dataset%cur_value)) deallocate(dataset%cur_value)
399) allocate(dataset%cur_value(dataset%rank))
400) dataset%cur_value(1:dataset%rank) = dataset%values(1:dataset%rank,1)
401)
402) dataset%time_shift = dataset%times(dataset%max_time_index)
403)
404) end subroutine UniformVelocityDatasetVerify
405)
406) ! ************************************************************************** !
407)
408) subroutine UniformVelocityDatasetUpdate(option,time,dataset)
409) !
410) ! Updates a velocity dataset
411) !
412) ! Author: Glenn Hammond
413) ! Date: 06/02/09
414) !
415)
416) use Option_module
417)
418) implicit none
419)
420) type(option_type) :: option
421) PetscReal :: time
422) PetscBool :: is_cyclic
423) PetscInt :: interpolation_method
424) type(uniform_velocity_dataset_type) :: dataset
425)
426) PetscInt :: irank
427) PetscInt :: cur_time_index
428) PetscInt :: next_time_index
429) PetscReal :: time_fraction
430)
431) ! potentially for initial condition
432) if (time < 1.d-40 .and. .not.dataset%is_transient) return
433)
434) ! cycle times if at max_time_index and cyclic
435) if (dataset%cur_time_index == dataset%max_time_index .and. &
436) dataset%is_cyclic .and. dataset%max_time_index > 1) then
437)
438) do cur_time_index = 1, dataset%max_time_index
439) dataset%times(cur_time_index) = dataset%times(cur_time_index) + &
440) dataset%time_shift
441) enddo
442) dataset%cur_time_index = 1
443) endif
444)
445) cur_time_index = dataset%cur_time_index
446) next_time_index = min(dataset%cur_time_index+1, &
447) dataset%max_time_index)
448)
449) ! ensure that condition has started
450) if (time >= dataset%times(cur_time_index) .or. &
451) dabs(time-dataset%times(cur_time_index)) < 1.d-40) then
452)
453) ! find appropriate time interval
454) do
455) if (time < dataset%times(next_time_index) .or. &
456) cur_time_index == next_time_index) &
457) exit
458) cur_time_index = next_time_index
459) ! ensure that time index does not go beyond end of array
460) if (next_time_index < dataset%max_time_index) &
461) next_time_index = next_time_index + 1
462) enddo
463)
464) dataset%cur_time_index = cur_time_index
465)
466) select case(dataset%interpolation_method)
467) case(STEP) ! just use the current value
468) do irank=1,dataset%rank
469) dataset%cur_value(irank) = &
470) dataset%values(irank,dataset%cur_time_index)
471) enddo
472) case(LINEAR) ! interpolate the value between times
473) do irank=1,dataset%rank
474) ! interpolate value based on time
475) dataset%cur_time_index = cur_time_index
476) if (cur_time_index < dataset%max_time_index) then
477) time_fraction = (time-dataset%times(cur_time_index)) / &
478) (dataset%times(next_time_index) - &
479) dataset%times(cur_time_index))
480) dataset%cur_value(irank) = dataset%values(irank,cur_time_index) + &
481) time_fraction * &
482) (dataset%values(irank,next_time_index) - &
483) dataset%values(irank,cur_time_index))
484) else
485) dataset%cur_value(irank) = &
486) dataset%values(irank,dataset%max_time_index)
487) endif
488) enddo
489) end select
490) endif
491)
492) end subroutine UniformVelocityDatasetUpdate
493)
494) ! ************************************************************************** !
495)
496) subroutine UniformVelocityDatasetDestroy(dataset)
497) !
498) ! Destroys a velocity dataset
499) !
500) ! Author: Glenn Hammond
501) ! Date: 06/02/09
502) !
503) use Utility_module, only : DeallocateArray
504)
505) implicit none
506)
507) type(uniform_velocity_dataset_type), pointer :: dataset
508)
509) if (.not.associated(dataset)) return
510)
511) call DeallocateArray(dataset%times)
512) call DeallocateArray(dataset%values)
513) call DeallocateArray(dataset%cur_value)
514)
515) deallocate(dataset)
516) nullify(dataset)
517)
518) end subroutine UniformVelocityDatasetDestroy
519)
520) end module Uniform_Velocity_module