Neko  0.8.1
A portable framework for high-order spectral element flow simulations
field_series.f90
Go to the documentation of this file.
1 ! Copyright (c) 2021-2023, The Neko Authors
2 ! All rights reserved.
3 !
4 ! Redistribution and use in source and binary forms, with or without
5 ! modification, are permitted provided that the following conditions
6 ! are met:
7 !
8 ! * Redistributions of source code must retain the above copyright
9 ! notice, this list of conditions and the following disclaimer.
10 !
11 ! * Redistributions in binary form must reproduce the above
12 ! copyright notice, this list of conditions and the following
13 ! disclaimer in the documentation and/or other materials provided
14 ! with the distribution.
15 !
16 ! * Neither the name of the authors nor the names of its
17 ! contributors may be used to endorse or promote products derived
18 ! from this software without specific prior written permission.
19 !
20 ! THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS
21 ! "AS IS" AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT
22 ! LIMITED TO, THE IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS
23 ! FOR A PARTICULAR PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL THE
24 ! COPYRIGHT OWNER OR CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT,
25 ! INCIDENTAL, SPECIAL, EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING,
26 ! BUT NOT LIMITED TO, PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES;
27 ! LOSS OF USE, DATA, OR PROFITS; OR BUSINESS INTERRUPTION) HOWEVER
28 ! CAUSED AND ON ANY THEORY OF LIABILITY, WHETHER IN CONTRACT, STRICT
29 ! LIABILITY, OR TORT (INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN
30 ! ANY WAY OUT OF THE USE OF THIS SOFTWARE, EVEN IF ADVISED OF THE
31 ! POSSIBILITY OF SUCH DAMAGE.
32 !
35  use num_types
36  use field
37  implicit none
38  private
39 
40  type, public :: field_series_t
41  type(field_t), pointer :: f => null()
42  type(field_t), allocatable :: lf(:)
43  integer, private :: len = 0
44  contains
45  procedure, pass(this) :: init => field_series_init
46  procedure, pass(this) :: free => field_series_free
47  procedure, pass(this) :: update => field_series_update
48  procedure, pass(this) :: set => field_series_set
49  procedure, pass(this) :: size => field_series_size
50  end type field_series_t
51 
52 contains
53 
55  subroutine field_series_init(this, f, len)
56  class(field_series_t), intent(inout) :: this
57  type(field_t), intent(inout), target :: f
58  integer :: len
59  character(len=80) :: name
60  character(len=5) :: id_str
61  integer :: i
62 
63  call this%free()
64 
65  this%f => f
66  this%len = len
67 
68  allocate(this%lf(len))
69 
70  do i = 1, this%len
71  write(id_str, '(I0)') i
72  name = trim(f%name)//'_lag'//id_str
73  call this%lf(i)%init(this%f%dof, name)
74  end do
75 
76  end subroutine field_series_init
77 
79  subroutine field_series_free(this)
80  class(field_series_t), intent(inout) :: this
81  integer :: i
82 
83  if (associated(this%f)) then
84  nullify(this%f)
85  end if
86 
87  do i = 1, this%len
88  call this%lf(i)%free()
89  end do
90 
91  end subroutine field_series_free
92 
94  function field_series_size(this) result(len)
95  class(field_series_t), intent(in) :: this
96  integer :: len
97  len = this%len
98  end function field_series_size
99 
101  subroutine field_series_update(this)
102  class(field_series_t), intent(inout) :: this
103  integer :: i
104 
105  do i = this%len, 2, -1
106  this%lf(i) = this%lf(i-1)
107  end do
108 
109  this%lf(1) = this%f
110 
111  end subroutine field_series_update
112 
114  subroutine field_series_set(this, g)
115  class(field_series_t), intent(inout) :: this
116  type(field_t), intent(in) :: g
117  integer :: i
118 
119  do i = 1, this%len
120  this%lf(i) = g
121  end do
122 
123  end subroutine field_series_set
124 
125 end module field_series
Stores a series fields.
subroutine field_series_set(this, g)
Set all fields in a series to g.
subroutine field_series_update(this)
Update a field series (evict oldest entry)
integer function field_series_size(this)
Return the size of the field series.
subroutine field_series_free(this)
Deallocates a field series.
subroutine field_series_init(this, f, len)
Initialize a field series of length len for a field f.
Defines a field.
Definition: field.f90:34