Neko  0.8.99
A portable framework for high-order spectral element flow simulations
fld_file_data.f90
Go to the documentation of this file.
1 
8  use num_types, only : rp
9  use math
10  use vector, only : vector_t, vector_ptr_t
11  implicit none
12  private
13 
14  type, public :: fld_file_data_t
15  type(vector_t) :: x
16  type(vector_t) :: y
17  type(vector_t) :: z
18  type(vector_t) :: u
19  type(vector_t) :: v
20  type(vector_t) :: w
21  type(vector_t) :: p
22  type(vector_t) :: t
23  integer, allocatable :: idx(:)
24  type(vector_t), allocatable :: s(:)
25  integer :: gdim
26  integer :: n_scalars = 0
27  real(kind=rp) :: time = 0.0
28  integer :: glb_nelv = 0
29  integer :: nelv = 0
30  integer :: offset_el = 0
31  integer :: lx = 0
32  integer :: ly = 0
33  integer :: lz = 0
34  integer :: t_counter = 0
35  ! meta file information (if any)
36  integer :: meta_nsamples = 0
37  integer :: meta_start_counter = 0
38  character(len=1024) :: fld_series_fname
39  contains
40  procedure, pass(this) :: init => fld_file_data_init
41  procedure, pass(this) :: free => fld_file_data_free
42  procedure, pass(this) :: scale => fld_file_data_scale
43  procedure, pass(this) :: add => fld_file_data_add
44  procedure, pass(this) :: size => fld_file_data_size
45  procedure, pass(this) :: get_list => fld_file_get_list
46  end type fld_file_data_t
47 
48 contains
49 
51  subroutine fld_file_data_init(this, nelv, offset_el)
52  class(fld_file_data_t), intent(inout) :: this
53  integer, intent(in), optional :: nelv, offset_el
54  call this%free()
55  if (present(nelv)) this%nelv = nelv
56  if (present(offset_el)) this%offset_el = offset_el
57 
58  end subroutine fld_file_data_init
60  function fld_file_data_size(this) result(i)
61  class(fld_file_data_t) :: this
62  integer :: i
63  i = 0
64  if(this%u%n .gt. 0) i = i + 1
65  if(this%v%n .gt. 0) i = i + 1
66  if(this%w%n .gt. 0) i = i + 1
67  if(this%p%n .gt. 0) i = i + 1
68  if(this%t%n .gt. 0) i = i + 1
69  i = i + this%n_scalars
70 
71  end function fld_file_data_size
72 
74  subroutine fld_file_get_list(this, ptr_list, n)
75  class(fld_file_data_t), target, intent(in) :: this
76  integer, intent(in) :: n
77  integer :: i, j
78  type(vector_ptr_t), intent(inout) :: ptr_list(n)
79  i = 1
80  if(this%u%n .gt. 0) then
81  ptr_list(i)%ptr => this%u
82  i = i + 1
83  end if
84  if(this%v%n .gt. 0) then
85  ptr_list(i)%ptr => this%v
86  i = i + 1
87  end if
88  if(this%w%n .gt. 0) then
89  ptr_list(i)%ptr => this%w
90  i = i + 1
91  end if
92  if(this%p%n .gt. 0) then
93  ptr_list(i)%ptr => this%p
94  i = i + 1
95  end if
96  if(this%t%n .gt. 0) then
97  ptr_list(i)%ptr => this%t
98  i = i + 1
99  end if
100  do j = 1, this%n_scalars
101  ptr_list(i)%ptr => this%s(j)
102  i = i +1
103  end do
104 
105  end subroutine fld_file_get_list
106 
107 
108 
110  subroutine fld_file_data_scale(this, c)
111  class(fld_file_data_t), intent(inout) :: this
112  real(kind=rp), intent(in) :: c
113  integer :: i
114 
115  if(this%u%n .gt. 0) call cmult(this%u%x,c,this%u%n)
116  if(this%v%n .gt. 0) call cmult(this%v%x,c,this%v%n)
117  if(this%w%n .gt. 0) call cmult(this%w%x,c,this%w%n)
118  if(this%p%n .gt. 0) call cmult(this%p%x,c,this%p%n)
119  if(this%t%n .gt. 0) call cmult(this%t%x,c,this%t%n)
120 
121  do i = 1, this%n_scalars
122  if(this%s(i)%n .gt. 0) call cmult(this%s(i)%x,c,this%s(i)%n)
123  end do
124 
125  end subroutine fld_file_data_scale
126 
128  subroutine fld_file_data_add(this, fld_data_add)
129  class(fld_file_data_t), intent(inout) :: this
130  class(fld_file_data_t), intent(in) :: fld_data_add
131  integer :: i
132 
133  if(this%u%n .gt. 0) call add2(this%u%x,fld_data_add%u%x,this%u%n)
134  if(this%v%n .gt. 0) call add2(this%v%x,fld_data_add%v%x,this%v%n)
135  if(this%w%n .gt. 0) call add2(this%w%x,fld_data_add%w%x,this%w%n)
136  if(this%p%n .gt. 0) call add2(this%p%x,fld_data_add%p%x,this%p%n)
137  if(this%t%n .gt. 0) call add2(this%t%x,fld_data_add%t%x,this%t%n)
138 
139  do i = 1, this%n_scalars
140  if(this%s(i)%n .gt. 0) call add2(this%s(i)%x,fld_data_add%s(i)%x,this%s(i)%n)
141  end do
142  end subroutine fld_file_data_add
143 
145  subroutine fld_file_data_free(this)
146  class(fld_file_data_t), intent(inout) :: this
147  integer :: i
148  call this%x%free()
149  call this%y%free()
150  call this%z%free()
151  call this%u%free()
152  call this%v%free()
153  call this%w%free()
154  call this%p%free()
155  call this%t%free()
156  if (allocated(this%s)) then
157  do i = 1, this%n_scalars
158  call this%s(i)%free()
159  end do
160  end if
161  this%n_scalars = 0
162  this%time = 0.0
163  this%glb_nelv = 0
164  this%nelv = 0
165  this%offset_el = 0
166  this%lx = 0
167  this%ly = 0
168  this%lz = 0
169  this%t_counter = 0
170  this%meta_nsamples = 0
171  this%meta_start_counter = 0
172  end subroutine fld_file_data_free
173 
174 end module fld_file_data
Simple module to handle fld file series. Provides an interface to the different fields sotred in a fl...
integer function fld_file_data_size(this)
Get number of fields in this fld file.
subroutine fld_file_data_add(this, fld_data_add)
Add the values in another fld file to this.
subroutine fld_file_get_list(this, ptr_list, n)
Get a list with pointers to the fields in the fld file.
subroutine fld_file_data_free(this)
Deallocate fld file data type.
subroutine fld_file_data_init(this, nelv, offset_el)
Initialise a fld_file_data object with nelv elements with a offset_nel.
subroutine fld_file_data_scale(this, c)
Scale the values stored in this fld_file_data.
Definition: math.f90:60
subroutine, public cmult(a, c, n)
Multiplication by constant c .
Definition: math.f90:277
subroutine, public add2(a, b, n)
Vector addition .
Definition: math.f90:547
integer, parameter, public rp
Global precision used in computations.
Definition: num_types.f90:12
Defines a vector.
Definition: vector.f90:34