Neko  0.9.0
A portable framework for high-order spectral element flow simulations
data_streamer.F90
Go to the documentation of this file.
1 ! Copyright (c) 2020-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, only: rp, c_rp
36  use field, only: field_t
37  use coefs, only: coef_t
38  use utils, only: neko_warning
39  use device
40  use comm
41  use neko_mpi_types
42  use neko_config
43  use, intrinsic :: iso_c_binding
44  implicit none
45  private
46 
56  type, public :: data_streamer_t
58  integer :: if_asynch
60  integer, allocatable :: lglel(:)
61  contains
63  procedure, pass(this) :: init => data_streamer_init
65  procedure, pass(this) :: free => data_streamer_free
67  procedure, pass(this) :: stream => data_streamer_stream
69  procedure, pass(this) :: recieve => data_streamer_recieve
70 
71  end type data_streamer_t
72 
73 contains
74 
81  subroutine data_streamer_init(this, coef)
82  class(data_streamer_t), intent(inout) :: this
83  type(coef_t), intent(inout) :: coef
84  integer :: nelb, nelv, nelgv, npts, gdim
85 
86  !Assign the set up parameters
87  nelv = coef%msh%nelv
88  npts = coef%Xh%lx*coef%Xh%ly*coef%Xh%lz
89  nelgv = coef%msh%glb_nelv
90  nelb = coef%msh%offset_el
91  gdim = coef%msh%gdim
92 
93 #ifdef HAVE_ADIOS2
94  call fortran_adios2_initialize(npts, nelv, nelb, nelgv, gdim, neko_comm)
95 #else
96  call neko_warning('Is not being built with ADIOS2 support.')
97  call neko_warning('Not able to use stream/compression functionality')
98 #endif
99 
100 
101  end subroutine data_streamer_init
102 
105  subroutine data_streamer_free(this)
106  class(data_streamer_t), intent(inout) :: this
107 
108 #ifdef HAVE_ADIOS2
109  call fortran_adios2_finalize()
110 #else
111  call neko_warning('Is not being built with ADIOS2 support.')
112  call neko_warning('Not able to use stream/compression functionality')
113 #endif
114 
115  end subroutine data_streamer_free
116 
119  subroutine data_streamer_stream(this, fld)
120  class(data_streamer_t), intent(inout) :: this
121  real(kind=rp), dimension(:,:,:,:), intent(inout) :: fld
122 
123 #ifdef HAVE_ADIOS2
124  call fortran_adios2_stream(fld)
125 #else
126  call neko_warning('Is not being built with ADIOS2 support.')
127  call neko_warning('Not able to use stream/compression functionality')
128 #endif
129 
130  end subroutine data_streamer_stream
131 
134  subroutine data_streamer_recieve(this, fld)
135  class(data_streamer_t), intent(inout) :: this
136  real(kind=rp), dimension(:,:,:,:), intent(inout) :: fld
137 
138 #ifdef HAVE_ADIOS2
139  call fortran_adios2_recieve(fld)
140 #else
141  call neko_warning('Is not being built with ADIOS2 support.')
142  call neko_warning('Not able to use stream/compression functionality')
143 #endif
144 
145  end subroutine data_streamer_recieve
146 
147 
148 #ifdef HAVE_ADIOS2
149 
159  subroutine fortran_adios2_initialize(npts, nelv, nelb, nelgv, gdim, comm)
160  use, intrinsic :: iso_c_binding
161  implicit none
162  integer, intent(in) :: npts, nelv, nelb, nelgv, gdim
163  type(mpi_comm) :: comm
164 
165  interface
166 
171  subroutine c_adios2_initialize(npts, nelv, nelb, nelgv, gdim, &
172  comm) bind(C,name="adios2_initialize_")
173  use, intrinsic :: iso_c_binding
174  import c_rp
175  implicit none
176  integer(kind=C_INT) :: npts
177  integer(kind=C_INT) :: nelv
178  integer(kind=C_INT) :: nelb
179  integer(kind=C_INT) :: nelgv
180  integer(kind=C_INT) :: gdim
181  type(*) :: comm
182  end subroutine c_adios2_initialize
183  end interface
184 
185  call c_adios2_initialize(npts, nelv, nelb, nelgv, gdim, comm)
186  end subroutine fortran_adios2_initialize
187 
190  subroutine fortran_adios2_finalize()
191  use, intrinsic :: iso_c_binding
192  implicit none
193 
194  interface
195 
196  subroutine c_adios2_finalize() bind(C,name="adios2_finalize_")
197  use, intrinsic :: iso_c_binding
198  implicit none
199  end subroutine c_adios2_finalize
200  end interface
201 
202  call c_adios2_finalize()
203  end subroutine fortran_adios2_finalize
204 
210  subroutine fortran_adios2_stream(fld)
211  use, intrinsic :: iso_c_binding
212  implicit none
213  real(kind=rp), dimension(:,:,:,:), intent(inout) :: fld
214 
215  interface
216 
217  subroutine c_adios2_stream(fld) &
218  bind(C,name="adios2_stream_")
219  use, intrinsic :: iso_c_binding
220  import c_rp
221  implicit none
222  real(kind=c_rp), intent(INOUT) :: fld(*)
223  end subroutine c_adios2_stream
224  end interface
225 
226  call c_adios2_stream(fld)
227  end subroutine fortran_adios2_stream
228 
234  subroutine fortran_adios2_recieve(fld)
235  use, intrinsic :: iso_c_binding
236  implicit none
237  real(kind=rp), dimension(:,:,:,:), intent(inout) :: fld
238 
239  interface
240 
241  subroutine c_adios2_recieve(fld) &
242  bind(C,name="adios2_recieve_")
243  use, intrinsic :: iso_c_binding
244  import c_rp
245  implicit none
246  real(kind=c_rp), intent(INOUT) :: fld(*)
247  end subroutine c_adios2_recieve
248  end interface
249 
250  call c_adios2_recieve(fld)
251  end subroutine fortran_adios2_recieve
252 
253 #endif
254 
255 end module data_streamer
Coefficients.
Definition: coef.f90:34
Definition: comm.F90:1
type(mpi_comm) neko_comm
MPI communicator.
Definition: comm.F90:16
Implements type data_streamer_t.
subroutine data_streamer_recieve(this, fld)
reciever
subroutine data_streamer_init(this, coef)
Constructor Wraps the adios2 set-up.
subroutine data_streamer_stream(this, fld)
streamer
subroutine data_streamer_free(this)
Destructor wraps the adios2 finalize routine. Closes insitu writer.
Device abstraction, common interface for various accelerators.
Definition: device.F90:34
Defines a field.
Definition: field.f90:34
Build configurations.
Definition: neko_config.f90:34
MPI derived types.
Definition: mpi_types.f90:34
integer, parameter, public c_rp
Definition: num_types.f90:13
integer, parameter, public rp
Global precision used in computations.
Definition: num_types.f90:12
Utilities.
Definition: utils.f90:35
subroutine, public neko_warning(warning_msg)
Reports a warning to standard output.
Definition: utils.f90:245
Coefficients defined on a given (mesh, ) tuple. Arrays use indices (i,j,k,e): element e,...
Definition: coef.f90:55
Provides access to data streaming by interfacing with c++ ADIOS2 subroutines.