Neko 1.99.1
A portable framework for high-order spectral element flow simulations
Loading...
Searching...
No Matches
data_streamer.F90
Go to the documentation of this file.
1! Copyright (c) 2020-2025, 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 comm, only : neko_comm
40 use mpi_f08, only : mpi_comm
41 use, intrinsic :: iso_c_binding
42 implicit none
43 private
44
54 type, public :: data_streamer_t
56 integer :: if_asynch
58 integer, allocatable :: lglel(:)
59 contains
61 procedure, pass(this) :: init => data_streamer_init
63 procedure, pass(this) :: free => data_streamer_free
65 procedure, pass(this) :: stream => data_streamer_stream
67 procedure, pass(this) :: recieve => data_streamer_recieve
68
69 end type data_streamer_t
70
71contains
72
79 subroutine data_streamer_init(this, coef)
80 class(data_streamer_t), intent(inout) :: this
81 type(coef_t), intent(inout) :: coef
82 integer :: nelb, nelv, nelgv, npts, gdim
83
84 !Assign the set up parameters
85 nelv = coef%msh%nelv
86 npts = coef%Xh%lx*coef%Xh%ly*coef%Xh%lz
87 nelgv = coef%msh%glb_nelv
88 nelb = coef%msh%offset_el
89 gdim = coef%msh%gdim
90
91#ifdef HAVE_ADIOS2
92 call fortran_adios2_initialize(npts, nelv, nelb, nelgv, gdim, neko_comm)
93#else
94 call neko_warning('Is not being built with ADIOS2 support.')
95 call neko_warning('Not able to use stream/compression functionality')
96#endif
97
98
99 end subroutine data_streamer_init
100
103 subroutine data_streamer_free(this)
104 class(data_streamer_t), intent(inout) :: this
105
106#ifdef HAVE_ADIOS2
107 call fortran_adios2_finalize()
108#else
109 call neko_warning('Is not being built with ADIOS2 support.')
110 call neko_warning('Not able to use stream/compression functionality')
111#endif
112
113 end subroutine data_streamer_free
114
117 subroutine data_streamer_stream(this, fld)
118 class(data_streamer_t), intent(inout) :: this
119 real(kind=rp), dimension(:,:,:,:), intent(inout) :: fld
120
121#ifdef HAVE_ADIOS2
122 call fortran_adios2_stream(fld)
123#else
124 call neko_warning('Is not being built with ADIOS2 support.')
125 call neko_warning('Not able to use stream/compression functionality')
126#endif
127
128 end subroutine data_streamer_stream
129
132 subroutine data_streamer_recieve(this, fld)
133 class(data_streamer_t), intent(inout) :: this
134 real(kind=rp), dimension(:,:,:,:), intent(inout) :: fld
135
136#ifdef HAVE_ADIOS2
137 call fortran_adios2_recieve(fld)
138#else
139 call neko_warning('Is not being built with ADIOS2 support.')
140 call neko_warning('Not able to use stream/compression functionality')
141#endif
142
143 end subroutine data_streamer_recieve
144
145
146#ifdef HAVE_ADIOS2
147
157 subroutine fortran_adios2_initialize(npts, nelv, nelb, nelgv, gdim, comm)
158 use, intrinsic :: iso_c_binding
159 implicit none
160 integer, intent(in) :: npts, nelv, nelb, nelgv, gdim
161 type(mpi_comm) :: comm
162
163 interface
164
169 subroutine c_adios2_initialize(npts, nelv, nelb, nelgv, gdim, &
170 comm) bind(C,name="adios2_initialize_")
171 use, intrinsic :: iso_c_binding
172 import c_rp
173 implicit none
174 integer(kind=C_INT) :: npts
175 integer(kind=C_INT) :: nelv
176 integer(kind=C_INT) :: nelb
177 integer(kind=C_INT) :: nelgv
178 integer(kind=C_INT) :: gdim
179 type(*) :: comm
180 end subroutine c_adios2_initialize
181 end interface
182
183 call c_adios2_initialize(npts, nelv, nelb, nelgv, gdim, comm)
184 end subroutine fortran_adios2_initialize
185
188 subroutine fortran_adios2_finalize()
189 use, intrinsic :: iso_c_binding
190 implicit none
191
192 interface
193
194 subroutine c_adios2_finalize() bind(C,name="adios2_finalize_")
195 use, intrinsic :: iso_c_binding
196 implicit none
197 end subroutine c_adios2_finalize
198 end interface
199
200 call c_adios2_finalize()
201 end subroutine fortran_adios2_finalize
202
208 subroutine fortran_adios2_stream(fld)
209 use, intrinsic :: iso_c_binding
210 implicit none
211 real(kind=rp), dimension(:,:,:,:), intent(inout) :: fld
212
213 interface
214
215 subroutine c_adios2_stream(fld) &
216 bind(C,name="adios2_stream_")
217 use, intrinsic :: iso_c_binding
218 import c_rp
219 implicit none
220 real(kind=c_rp), intent(INOUT) :: fld(*)
221 end subroutine c_adios2_stream
222 end interface
223
224 call c_adios2_stream(fld)
225 end subroutine fortran_adios2_stream
226
232 subroutine fortran_adios2_recieve(fld)
233 use, intrinsic :: iso_c_binding
234 implicit none
235 real(kind=rp), dimension(:,:,:,:), intent(inout) :: fld
236
237 interface
238
239 subroutine c_adios2_recieve(fld) &
240 bind(C,name="adios2_recieve_")
241 use, intrinsic :: iso_c_binding
242 import c_rp
243 implicit none
244 real(kind=c_rp), intent(INOUT) :: fld(*)
245 end subroutine c_adios2_recieve
246 end interface
247
248 call c_adios2_recieve(fld)
249 end subroutine fortran_adios2_recieve
250
251#endif
252
253end module data_streamer
Coefficients.
Definition coef.f90:34
Definition comm.F90:1
type(mpi_comm), public neko_comm
MPI communicator.
Definition comm.F90:42
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.
Defines a field.
Definition field.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:284
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.