vecmem 1.23.0
Loading...
Searching...
No Matches
copy.ipp
1/*
2 * VecMem project, part of the ACTS project (R&D line)
3 *
4 * (c) 2021-2026 CERN for the benefit of the ACTS project
5 *
6 * Mozilla Public License Version 2.0
7 */
8#pragma once
9
10// VecMem include(s).
11#include "vecmem/containers/details/resize_jagged_vector.hpp"
12#include "vecmem/containers/jagged_vector.hpp"
13#include "vecmem/edm/details/schema_traits.hpp"
14#include "vecmem/memory/get_default_resource.hpp"
15#include "vecmem/utils/debug.hpp"
16#include "vecmem/utils/type_traits.hpp"
17
18// System include(s).
19#include <algorithm>
20#include <cassert>
21#include <numeric>
22#include <sstream>
23#include <stdexcept>
24#include <type_traits>
25
26namespace vecmem {
27
28template <typename TYPE>
29copy::event_type copy::setup(data::vector_view<TYPE> data) const {
30
31 // Check if anything needs to be done.
32 if ((data.size_ptr() == nullptr) || (data.capacity() == 0)) {
34 }
35
36 // Initialize the "size variable" correctly on the buffer.
38 data.size_ptr(), 0);
39 VECMEM_DEBUG_MSG(2,
40 "Prepared a device vector buffer of capacity %u "
41 "for use on a device (ptr: %p)",
42 data.capacity(), static_cast<void*>(data.size_ptr()));
43
44 // Return a new event.
45 return create_event();
46}
47
48template <typename TYPE>
49copy::event_type copy::memset(data::vector_view<TYPE> data, int value) const {
50
51 // Check if anything needs to be done.
52 if (data.capacity() == 0) {
54 }
55
56 // Call memset with the correct arguments.
57 do_memset(data.capacity() * sizeof(TYPE), data.ptr(), value);
58 VECMEM_DEBUG_MSG(2, "Set %u vector elements to %i at ptr: %p",
59 data.capacity(), value, static_cast<void*>(data.ptr()));
60
61 // Return a new event.
62 return create_event();
63}
64
65template <typename TYPE>
67 const vecmem::data::vector_view<TYPE>& data, memory_resource& resource,
68 type::copy_type cptype) const {
69
70 // Set up the result buffer. No need to call setup(...) on it, as the buffer
71 // is not resizable.
73 resource);
74
75 // Copy the payload of the vector. Explicitly waiting for the copy to finish
76 // before returning the buffer.
77 operator()(data, result, cptype)->wait();
78
79 // Return the buffer.
80 return result;
81}
82
83template <typename TYPE>
85 const data::vector_view<std::add_const_t<TYPE>>& from_view,
86 data::vector_view<TYPE> to_view, type::copy_type cptype) const {
87
88 // Perform the copy. Depending on whether an actual copy has been set up /
89 // performed, return either "an actual", or just a dummy event.
90 if (copy_view_impl(from_view, to_view, cptype)) {
91 return create_event();
92 } else {
94 }
95}
96
97template <typename TYPE, typename ALLOC>
99 const data::vector_view<std::add_const_t<TYPE>>& from_view,
100 std::vector<TYPE, ALLOC>& to_vec, type::copy_type cptype) const {
101
102 // Figure out the size of the buffer.
103 const typename data::vector_view<std::add_const_t<TYPE>>::size_type size =
104 get_size(from_view);
105
106 // If the size is zero, don't do an actual copy.
107 if (size == 0u) {
109 }
110
111 // Make the target vector the correct size.
112 to_vec.resize(size);
113 // Perform the memory copy.
114 do_copy(size * sizeof(TYPE), from_view.ptr(), to_vec.data(), cptype);
115
116 // Return a new event.
117 return create_event();
118}
119
120template <typename TYPE>
122 const data::vector_view<TYPE>& data) const {
123
124 // Handle the simple case, when the view/buffer is not resizable.
125 if (data.size_ptr() == nullptr) {
126 return data.capacity();
127 }
128
129 // If it *is* resizable, don't assume that the size is host-accessible.
130 // Explicitly copy it for access.
133 data.size_ptr(), &result, type::unknown);
134
135 // Wait for the copy operation to finish. With some backends
136 // (khm... SYCL... khm...) copies can be asynchronous even into
137 // non-pinned host memory.
138 create_event()->wait();
139
140 // Return what we got.
141 return result;
142}
143
144template <typename TYPE>
146 const data::vector_view<TYPE>& data, memory_resource& pinnedHostMr) const {
147
148 // Helper result value type.
149 using value_type = typename data::vector_view<TYPE>::size_type;
150
151 // Handle the case when the view/buffer is not resizable.
152 if (data.size_ptr() == nullptr) {
153 // Create the result value, in "default" memory.
154 auto value = make_unique_alloc<value_type>(*(get_default_resource()));
155 *value = data.capacity();
156 return {std::move(value), vecmem::copy::create_event()};
157 }
158
159 // Create the result value, in pinned host memory.
160 auto value = make_unique_alloc<value_type>(pinnedHostMr);
161
162 // Set up the copy into a pinned host memory buffer.
163 do_copy(sizeof(value_type), data.size_ptr(), value.get(), type::unknown);
164
165 // Return the appropriate "future" value.
166 return {std::move(value), create_event()};
167}
168
169template <typename TYPE>
170copy::event_type copy::setup(data::jagged_vector_view<TYPE> data) const {
171
172 // Check if anything needs to be done.
173 if (data.size() == 0) {
175 }
176
177 // "Set up" the inner vector descriptors, using the host-accessible data.
178 // But only if the jagged vector buffer is resizable.
179 if (data.host_ptr()[0].size_ptr() != nullptr) {
180 do_memset(
181 sizeof(typename data::vector_buffer<TYPE>::size_type) * data.size(),
182 data.host_ptr()[0].size_ptr(), 0);
183 }
184
185 // Check if anything else needs to be done.
186 if (data.ptr() == data.host_ptr()) {
187 return create_event();
188 }
189
190 // Copy the description of the inner vectors of the buffer.
191 do_copy(
192 data.size() *
193 sizeof(
195 data.host_ptr(), data.ptr(), type::host_to_device);
196 VECMEM_DEBUG_MSG(2,
197 "Prepared a jagged device vector buffer of size %lu "
198 "for use on a device",
199 data.size());
200
201 // Return a new event.
202 return create_event();
203}
204
205template <typename TYPE>
206copy::event_type copy::memset(data::jagged_vector_view<TYPE> data,
207 int value) const {
208
209 // Do different things for jagged vectors that are contiguous in memory.
210 if (is_contiguous(data.host_ptr(), data.capacity())) {
211 // For contiguous jagged vectors we can just do a single memset
212 // operation. First, let's calculate the "total size" of the jagged
213 // vector.
214 const std::size_t total_size = std::accumulate(
215 data.host_ptr(), data.host_ptr() + data.size(),
216 static_cast<std::size_t>(0u),
217 [](std::size_t sum, const data::vector_view<TYPE>& iv) {
218 return sum + iv.capacity();
219 });
220 // Find the first "inner vector" that has a non-zero capacity.
221 for (std::size_t i = 0; i < data.size(); ++i) {
222 data::vector_view<TYPE>& iv = data.host_ptr()[i];
223 if ((iv.capacity() != 0u) && (iv.ptr() != nullptr)) {
224 // Call memset with its help.
225 do_memset(total_size * sizeof(TYPE), iv.ptr(), value);
226 return create_event();
227 }
228 }
229 // If we are still here, apparently we didn't need to do anything.
231 } else {
232 // For non-contiguous jagged vectors call memset one-by-one on the
233 // inner vectors. Note that we don't use vecmem::copy::memset here,
234 // as that would require us to wait for each memset individually.
235 for (std::size_t i = 0; i < data.size(); ++i) {
236 data::vector_view<TYPE>& iv = data.host_ptr()[i];
237 do_memset(iv.capacity() * sizeof(TYPE), iv.ptr(), value);
238 }
239 }
240
241 // Return a new event.
242 return create_event();
243}
244
245template <typename TYPE>
247 const data::jagged_vector_view<TYPE>& data, memory_resource& resource,
248 memory_resource* host_access_resource, type::copy_type cptype) const {
249
250 // Create the result buffer object.
252 (((data.capacity() > 0u) && (data.host_ptr()[0].size_ptr() != nullptr))
257 assert(result.size() == data.size());
258
259 // Copy the description of the "inner vectors" if necessary.
260 if (host_access_resource != nullptr) {
261 setup(result)->wait();
262 }
263
264 // Copy the payload of the inner vectors. Explicitly waiting for the copy to
265 // finish before returning the buffer.
266 operator()(data, result, cptype)->wait();
267
268 // Return the newly created object.
269 return result;
270}
271
272template <typename TYPE>
274 const data::jagged_vector_view<std::add_const_t<TYPE>>& from_view,
275 data::jagged_vector_view<TYPE> to_view, type::copy_type cptype) const {
276
277 // Perform the copy. Depending on whether an actual copy has been set up /
278 // performed, return either "an actual", or just a dummy event.
279 if (copy_view_impl(from_view, to_view, cptype)) {
280 return create_event();
281 } else {
283 }
284}
285
286template <typename TYPE, typename ALLOC1, typename ALLOC2>
288 const data::jagged_vector_view<std::add_const_t<TYPE>>& from_view,
289 std::vector<std::vector<TYPE, ALLOC2>, ALLOC1>& to_vec,
290 type::copy_type cptype) const {
291
292 // Resize the output object to the correct size.
293 details::resize_jagged_vector(to_vec, from_view.size());
294 const auto sizes = get_sizes(from_view);
295 assert(sizes.size() == to_vec.size());
296 for (typename data::jagged_vector_view<std::add_const_t<TYPE>>::size_type
297 i = 0;
298 i < from_view.size(); ++i) {
299 to_vec[i].resize(sizes[i]);
300 }
301
302 // Perform the memory copy.
303 return operator()(from_view, vecmem::get_data(to_vec), cptype);
304}
305
306template <typename TYPE>
307std::vector<typename data::vector_view<TYPE>::size_type> copy::get_sizes(
308 const data::jagged_vector_view<TYPE>& data) const {
309
310 // Create the result vector.
311 std::vector<typename data::vector_view<TYPE>::size_type> result(data.size(),
312 0);
313
314 // Try to get the "resizable sizes" first.
315 for (std::size_t i = 0; i < data.size(); ++i) {
316 // Find the first "inner vector" that has a non-zero capacity, and is
317 // resizable.
318 if ((data.host_ptr()[i].capacity() != 0) &&
319 (data.host_ptr()[i].size_ptr() != nullptr)) {
320 // Copy the sizes of the inner vectors into the result vector.
322 (data.size() - i),
323 data.host_ptr()[i].size_ptr(), result.data() + i,
325 // Wait for the copy operation to finish. With some backends
326 // (khm... SYCL... khm...) copies can be asynchronous even into
327 // non-pinned host memory.
328 create_event()->wait();
329 // At this point the result vector should have been set up
330 // correctly.
331 return result;
332 }
333 }
334
335 // If we're still here, then the buffer is not resizable. So let's just
336 // collect the capacity of each of the inner vectors.
337 for (std::size_t i = 0; i < data.size(); ++i) {
338 result[i] = data.host_ptr()[i].capacity();
339 }
340 return result;
341}
342
343template <typename TYPE>
345 const std::vector<typename data::vector_view<TYPE>::size_type>& sizes,
347
348 // Finish early if possible.
349 if ((sizes.size() == 0) && (data.size() == 0)) {
351 }
352 // Make sure that the sizes match up.
353 if (sizes.size() != data.size()) {
354 std::ostringstream msg;
355 msg << "sizes.size() (" << sizes.size() << ") != data.size() ("
356 << data.size() << ")";
357 throw std::length_error(msg.str());
358 }
359 // Make sure that the target jagged vector is either resizable, or it has
360 // the correct sizes/capacities already.
361 bool perform_copy = true;
363 i < data.size(); ++i) {
364 // Perform a capacity check in all cases. Whether or not the target is
365 // resizable.
366 if (data.host_ptr()[i].capacity() < sizes[i]) {
367 std::ostringstream msg;
368 msg << "data.host_ptr()[" << i << "].capacity() ("
369 << data.host_ptr()[i].capacity() << ") < sizes[" << i << "] ("
370 << sizes[i] << ")";
371 throw std::length_error(msg.str());
372 }
373 // Make sure that the inner vectors are consistently either all
374 // resizable, or none of them are.
375 if (data.host_ptr()[i].size_ptr() == nullptr) {
376 perform_copy = false;
377 } else if (perform_copy == false) {
378 throw std::invalid_argument(
379 "Inconsistent target jagged vector view received for resizing");
380 }
381 }
382 // If no copy is necessary, we're done.
383 if (perform_copy == false) {
385 }
386 // Perform the copy with some internal knowledge of how resizable jagged
387 // vector buffers work.
388 do_copy(sizeof(typename data::vector_view<TYPE>::size_type) * sizes.size(),
389 sizes.data(), data.host_ptr()->size_ptr(), type::unknown);
390
391 // Return a new event.
392 return create_event();
393}
394
395template <typename TYPE>
396async_sizes<typename data::vector_view<TYPE>::size_type> copy::get_sizes(
397 const data::jagged_vector_view<TYPE>& data,
398 memory_resource& pinnedHostMr) const {
399
400 // Try to get the "resizable sizes" first.
401 for (std::size_t i = 0; i < data.size(); ++i) {
402 // Find the first "inner vector" that has a non-zero capacity, and is
403 // resizable.
404 if ((data.host_ptr()[i].capacity() != 0) &&
405 (data.host_ptr()[i].size_ptr() != nullptr)) {
406
407 // Create the result value in pinned host memory.
408 vector<typename data::vector_view<TYPE>::size_type> result(
409 data.size(), 0, &pinnedHostMr);
410
411 // Copy the sizes of the inner vectors into the result vector.
413 (data.size() - i),
414 data.host_ptr()[i].size_ptr(), result.data() + i,
416 // Return the appropriate "future" value.
417 return {std::move(result), create_event()};
418 }
419 }
420
421 // If we're still here, then the buffer is not resizable. So let's just
422 // collect the capacity of each of the inner vectors.
423 vector<typename data::vector_view<TYPE>::size_type> result(data.size(), 0);
424 for (std::size_t i = 0; i < data.size(); ++i) {
425 result[i] = data.host_ptr()[i].capacity();
426 }
427 return {std::move(result), vecmem::copy::create_event()};
428}
429
430template <typename SCHEMA>
431copy::event_type copy::setup(edm::view<SCHEMA> data) const {
432
433 // Copy the data layout to the device, if needed.
434 if (data.layout().ptr() != data.host_layout().ptr()) {
435 assert(data.layout().capacity() > 0u);
436 [[maybe_unused]] bool did_copy =
437 copy_view_impl(data.host_layout(), data.layout(), type::unknown);
438 assert(did_copy);
439 }
440
441 // Initialize the "size variable(s)" correctly on the buffer.
442 if (data.size().ptr() != nullptr) {
443 assert(data.size().capacity() > 0u);
444 do_memset(data.size().capacity() * sizeof(char), data.size().ptr(), 0);
445 }
446 VECMEM_DEBUG_MSG(3,
447 "Prepared an SoA container of capacity %u "
448 "for use on a device (layout: {%u, %p}, size: {%u, %p})",
449 data.capacity(), data.layout().size(),
450 static_cast<void*>(data.layout().ptr()),
451 data.size().size(), static_cast<void*>(data.size().ptr()));
452
453 // Return a new event.
454 return create_event();
455}
456
457template <typename... VARTYPES>
458copy::event_type copy::memset(edm::view<edm::schema<VARTYPES...>> data,
459 int value) const {
460
461 // For buffers, we can do this in one go.
462 if (data.payload().ptr() != nullptr) {
463 memset(data.payload(), value);
464 } else {
465 // Do the operation using the helper function, recursively.
466 memset_impl<0>(data, value);
467 }
468
469 // Return a new event.
470 return create_event();
471}
472
473template <typename... VARTYPES>
474edm::buffer<edm::details::remove_cv_t<edm::schema<VARTYPES...>>> copy::to(
475 const edm::view<edm::schema<VARTYPES...>>& data, memory_resource& resource,
476 [[maybe_unused]] memory_resource* host_access_resource,
477 type::copy_type cptype) const {
478
479 // Create the result buffer object.
481 ((data.size().ptr() != nullptr) ? data::buffer_type::resizable
483 edm::buffer<edm::details::remove_cv_t<edm::schema<VARTYPES...>>> result;
485 edm::schema<VARTYPES...>>::value) {
486 // Set up a buffer that has (a) jagged vector(s).
488 btype};
489 // The idea here is that since we are going to copy data into the newly
490 // created buffer, the only thing that setup(...) needs to do is to set
491 // up "the layout" in non-host-accessible memory. (Zeroing the size
492 // variable(s) is not needed.) When doing a device-to-host copy, copying
493 // "the layout" is not actually needed. Worse yet, the copy object is
494 // not capable of memsetting the size variable(s). Since the device
495 // specific copy object cannot memset host memory.
496 //
497 // Long story short, the entire setup is skipped in the absence of a
498 // host-accessible memory resource.
499 if (host_access_resource != nullptr) {
500 setup(result)->wait();
501 }
502 } else {
503 // Set up a buffer absent of jagged vectors.
504 result = {data.capacity(), resource, btype};
505 // In the absence of jagged vectors there is no "layout" to copy. And
506 // since the size of the buffer will be set during the copy correctly,
507 // there is nothing else that setup(...) would need to do. So no need
508 // to call it here.
509 }
510
511 // Perform the copy.
512 operator()(data, result, cptype)->wait();
513
514 // Return the buffer.
515 return result;
516}
517
518template <typename... VARTYPES>
520 const edm::view<edm::details::add_const_t<edm::schema<VARTYPES...>>>&
521 from_view,
522 edm::view<edm::schema<VARTYPES...>> to_view, type::copy_type cptype) const {
523
524 // First, handle the simple case, when both views have a contiguous memory
525 // layout.
526 if ((from_view.payload().ptr() != nullptr) &&
527 (to_view.payload().ptr() != nullptr) &&
528 (from_view.payload().capacity() == to_view.payload().capacity())) {
529
530 // If the "common size" is zero, we're done.
531 if (from_view.payload().capacity() == 0) {
533 }
534
535 // Let the user know what's happening.
536 VECMEM_DEBUG_MSG(2, "Performing simple SoA copy of %u bytes",
537 from_view.payload().size());
538
539 // Copy the payload with a single copy operation.
540 copy_view_impl(from_view.payload(), to_view.payload(), cptype);
541
542 // If the target view is resizable, set its size.
543 if (to_view.size().ptr() != nullptr) {
544 // If the source is also resizable, the situation should be simple.
545 if (from_view.size().ptr() != nullptr) {
546 // Check that the sizes are the same.
547 if (from_view.size().capacity() != to_view.size().capacity()) {
548 std::ostringstream msg;
549 msg << "from_view.size().capacity() ("
550 << from_view.size().capacity()
551 << ") != to_view.size().capacity() ("
552 << to_view.size().capacity() << ")";
553 throw std::length_error(msg.str());
554 }
555 // Perform a dumb copy.
556 copy_view_impl(from_view.size(), to_view.size(), cptype);
557 } else {
558 // If not, then copy the size(s) recursively.
559 copy_sizes_impl<0>(from_view, to_view, cptype);
560 }
561 }
562
563 // Create a synchronization event.
564 return create_event();
565 }
566
567 // If not, then do an un-optimized copy, variable-by-variable.
568 copy_payload_impl<0>(from_view, to_view, cptype);
569
570 // Return a new event.
571 return create_event();
572}
573
574template <typename... VARTYPES, template <typename> class INTERFACE>
576 const edm::view<edm::details::add_const_t<edm::schema<VARTYPES...>>>&
577 from_view,
578 edm::host<edm::schema<VARTYPES...>, INTERFACE>& to_vec,
579 type::copy_type cptype) const {
580
581 // Resize the output object to the correct size(s).
582 resize_impl<0>(from_view, to_vec, cptype);
583
584 // Perform the memory copy.
585 return operator()(from_view, vecmem::get_data(to_vec), cptype);
586}
587
588template <typename... VARTYPES>
590 const edm::view<edm::schema<VARTYPES...>>& data) const {
591
592 // Helper result value type.
593 using value_type = typename edm::view<edm::schema<VARTYPES...>>::size_type;
594
595 // Start by taking the capacity of the container.
596 value_type size = data.capacity();
597
598 // If there are jagged vectors in the container and/or the container is not
599 // resizable, we're done. It is done in two separate if statements to avoid
600 // MSVC trying to be too smart, and giving a warning...
601 if constexpr (std::disjunction_v<
603 return size;
604 } else {
605 // All the rest is put into an else block, to avoid MSVC trying to be
606 // too smart, and giving a warning about unreachable code...
607 if (data.size().ptr() == nullptr) {
608 return size;
609 }
610
611 // A small security check.
612 assert(data.size().size() == sizeof(value_type));
613
614 // Get the exact size of the container.
615 do_copy(sizeof(value_type), data.size().ptr(), &size, type::unknown);
616 // We have to wait for this to finish, since the "size" variable is
617 // not going to be available outside of this function. And
618 // asynchronous SYCL memory copies can happen from variables on the
619 // stack as well...
620 create_event()->wait();
621
622 // Return what we got.
623 return size;
624 }
625}
626
627template <typename... VARTYPES>
628async_size<typename edm::view<edm::schema<VARTYPES...>>::size_type>
630 memory_resource& pinnedHostMr) const {
631
632 // Helper result value type.
633 using value_type = typename edm::view<edm::schema<VARTYPES...>>::size_type;
634
635 // If there are jagged vectors in the container and/or the container is not
636 // resizable, we're done. It is done in two separate if statements to avoid
637 // MSVC trying to be too smart, and giving a warning...
638 if constexpr (std::disjunction_v<
640 // Create the result value, in "default" memory.
641 auto value = make_unique_alloc<value_type>(*(get_default_resource()));
642 *value = data.capacity();
643 return {std::move(value), vecmem::copy::create_event()};
644 } else {
645 // All the rest is put into an else block, to avoid MSVC trying to be
646 // too smart, and giving a warning about unreachable code...
647 if (data.size().ptr() == nullptr) {
648 // Create the result value, in "default" memory.
649 auto value =
650 make_unique_alloc<value_type>(*(get_default_resource()));
651 *value = data.capacity();
652 return {std::move(value), vecmem::copy::create_event()};
653 }
654
655 // A small security check.
656 assert(data.size().size() == sizeof(value_type));
657
658 // Create the result value, in pinned host memory.
659 auto value = make_unique_alloc<value_type>(pinnedHostMr);
660
661 // Get the exact size of the container.
662 do_copy(sizeof(value_type), data.size().ptr(), value.get(),
664
665 // Return the appropriate "future" value.
666 return {std::move(value), create_event()};
667 }
668}
669
670template <typename... VARTYPES>
671std::vector<data::vector_view<int>::size_type> copy::get_sizes(
672 const edm::view<edm::schema<VARTYPES...>>& data) const {
673
674 // Make sure that there's a jagged vector in here.
675 static_assert(
677 "Function can only be used on containers with jagged vectors");
678
679 // Get the sizes using the helper function.
680 return get_sizes_impl<0>(data);
681}
682
683template <typename... VARTYPES>
686 memory_resource& pinnedHostMr) const {
687
688 // Get the sizes using the helper function.
689 return get_sizes_impl<0>(data, pinnedHostMr);
690}
691
692template <typename TYPE>
693bool copy::copy_view_impl(
694 const data::vector_view<std::add_const_t<TYPE>>& from_view,
695 data::vector_view<TYPE> to_view, type::copy_type cptype) const {
696
697 // Get the size of the source view.
698 const typename data::vector_view<std::add_const_t<TYPE>>::size_type size =
699 get_size(from_view);
700 // If it's zero, don't do an actual copy.
701 if (size == 0u) {
702 return false;
703 }
704
705 // Make sure that the copy can happen.
706 if (to_view.capacity() < size) {
707 std::ostringstream msg;
708 msg << "Target capacity (" << to_view.capacity() << ") < source size ("
709 << size << ")";
710 throw std::length_error(msg.str());
711 }
712
713 // Make sure that if the target view is resizable, that it would be set up
714 // for the correct size.
715 if (to_view.size_ptr() != nullptr) {
716 // Select what type of copy this should be. Keeping in mind that we copy
717 // from a variable on the host stack. So the question is just whether
718 // the target is the host, or a device.
719 type::copy_type size_cptype = type::unknown;
720 switch (cptype) {
723 size_cptype = type::host_to_device;
724 break;
727 size_cptype = type::host_to_host;
728 break;
729 default:
730 break;
731 }
732 // Perform the copy.
733 do_copy(sizeof(typename data::vector_view<TYPE>::size_type), &size,
734 to_view.size_ptr(), size_cptype);
735 // We have to wait for this to finish, since the "size" variable is not
736 // going to be available outside of this function. And asynchronous SYCL
737 // memory copies can happen from variables on the stack as well...
738 create_event()->wait();
739 }
740
741 // Copy the payload.
742 do_copy(size * sizeof(TYPE), from_view.ptr(), to_view.ptr(), cptype);
743 return true;
744}
745
746template <typename TYPE>
747bool copy::copy_view_impl(
748 const data::jagged_vector_view<std::add_const_t<TYPE>>& from_view,
749 data::jagged_vector_view<TYPE> to_view, type::copy_type cptype) const {
750
751 // Sanity checks.
752 if (from_view.size() > to_view.size()) {
753 std::ostringstream msg;
754 msg << "from_view.size() (" << from_view.size()
755 << ") > to_view.size() (" << to_view.size() << ")";
756 throw std::length_error(msg.str());
757 }
758
759 // Get the "outer size" of the container.
760 const typename data::jagged_vector_view<std::add_const_t<TYPE>>::size_type
761 size = from_view.size();
762 if (size == 0u) {
763 return false;
764 }
765
766 // Calculate the contiguous-ness of the memory allocations.
767 const bool from_is_contiguous = is_contiguous(from_view.host_ptr(), size);
768 const bool to_is_contiguous = is_contiguous(to_view.host_ptr(), size);
769 VECMEM_DEBUG_MSG(3, "from_is_contiguous = %d, to_is_contiguous = %d",
770 from_is_contiguous, to_is_contiguous);
771
772 // Get the sizes of the source jagged vector.
773 const auto sizes = get_sizes(from_view);
774
775 // Before even attempting the copy, make sure that the target view either
776 // has the correct sizes, or can be resized correctly. Captire the event
777 // marking the end of the operation. We need to wait for the copy to finish
778 // before returning from this function.
779 auto set_sizes_event = set_sizes(sizes, to_view);
780
781 // Check whether the source and target capacities match up. We can only
782 // perform the "optimised copy" if they do.
783 std::vector<typename data::vector_view<std::add_const_t<TYPE>>::size_type>
784 capacities(size);
785 bool capacities_match = true;
786 for (std::size_t i = 0; i < size; ++i) {
787 if (from_view.host_ptr()[i].capacity() !=
788 to_view.host_ptr()[i].capacity()) {
789 capacities_match = false;
790 break;
791 }
792 capacities[i] = from_view.host_ptr()[i].capacity();
793 }
794
795 // Perform the copy as best as we can.
796 if (from_is_contiguous && to_is_contiguous && capacities_match) {
797 // Perform the copy in one go.
798 copy_views_contiguous_impl(capacities, from_view.host_ptr(),
799 to_view.host_ptr(), cptype);
800 } else {
801 // Do the copy as best as we can. Note that since they are not
802 // contiguous anyway, we use the sizes of the vectors here, not their
803 // capcities.
804 copy_views_impl(sizes, from_view.host_ptr(), to_view.host_ptr(),
805 cptype);
806 }
807
808 // Before returning, make sure that the "size set" operation has finished.
809 // Because the "sizes" variable is just about to go out of scope.
810 set_sizes_event->wait();
811 return true;
812}
813
814template <typename TYPE>
815void copy::copy_views_impl(
816 const std::vector<typename data::vector_view<TYPE>::size_type>& sizes,
817 const data::vector_view<std::add_const_t<TYPE>>* from_view,
818 data::vector_view<TYPE>* to_view, type::copy_type cptype) const {
819
820 // Some security checks.
821 assert(from_view != nullptr);
822 assert(to_view != nullptr);
823
824 // Helper variable(s) used in the copy.
825 const std::size_t size = sizes.size();
826 [[maybe_unused]] std::size_t copy_ops = 0;
827
828 // Perform the copy in multiple steps.
829 for (std::size_t i = 0; i < size; ++i) {
830
831 // Skip empty "inner vectors".
832 if (sizes[i] == 0) {
833 continue;
834 }
835
836 // Some sanity checks.
837 assert(from_view[i].ptr() != nullptr);
838 assert(to_view[i].ptr() != nullptr);
839 assert(sizes[i] <= from_view[i].capacity());
840 assert(sizes[i] <= to_view[i].capacity());
841
842 // Perform the copy.
843 do_copy(sizes[i] * sizeof(TYPE), from_view[i].ptr(), to_view[i].ptr(),
844 cptype);
845 ++copy_ops;
846 }
847
848 // Let the user know what happened.
849 VECMEM_DEBUG_MSG(2,
850 "Copied the payload of a jagged vector of type "
851 "\"%s\" with %lu copy operation(s)",
852 typeid(TYPE).name(), copy_ops);
853}
854
855template <typename TYPE>
856void copy::copy_views_contiguous_impl(
857 const std::vector<typename data::vector_view<TYPE>::size_type>& sizes,
858 const data::vector_view<std::add_const_t<TYPE>>* from_view,
859 data::vector_view<TYPE>* to_view, type::copy_type cptype) const {
860
861 // Some security checks.
862 assert(from_view != nullptr);
863 assert(to_view != nullptr);
864 assert(is_contiguous(from_view, sizes.size()));
865 assert(is_contiguous(to_view, sizes.size()));
866
867 // Helper variable(s) used in the copy.
868 const std::size_t size = sizes.size();
869 const std::size_t total_size =
870 std::accumulate(sizes.begin(), sizes.end(),
871 static_cast<std::size_t>(0)) *
872 sizeof(TYPE);
873
874 // Find the first non-empty element.
875 for (std::size_t i = 0; i < size; ++i) {
876
877 // Jump over empty elements.
878 if (sizes[i] == 0) {
879 continue;
880 }
881
882 // Some sanity checks.
883 assert(from_view[i].ptr() != nullptr);
884 assert(to_view[i].ptr() != nullptr);
885
886 // Perform the copy.
887 do_copy(total_size, from_view[i].ptr(), to_view[i].ptr(), cptype);
888 break;
889 }
890
891 // Let the user know what happened.
892 VECMEM_DEBUG_MSG(2,
893 "Copied the payload of a jagged vector of type "
894 "\"%s\" with 1 copy operation(s)",
895 typeid(TYPE).name());
896}
897
898template <typename TYPE>
899bool copy::is_contiguous(const data::vector_view<TYPE>* data,
900 std::size_t size) {
901
902 // We should never call this function for an empty jagged vector.
903 assert(size > 0);
904
905 // Check whether all memory blocks are contiguous.
906 auto ptr = data[0].ptr();
907 for (std::size_t i = 1; i < size; ++i) {
908 if ((ptr + data[i - 1].capacity()) != data[i].ptr()) {
909 return false;
910 }
911 ptr = data[i].ptr();
912 }
913 return true;
914}
915
916template <std::size_t INDEX, typename... VARTYPES>
917void copy::memset_impl(edm::view<edm::schema<VARTYPES...>> data,
918 int value) const {
919
920 // Scalars do not have their own dedicated @c memset functions.
921 if constexpr (edm::type::details::is_scalar<typename std::tuple_element<
922 INDEX, std::tuple<VARTYPES...>>::type>::value) {
923 do_memset(sizeof(typename std::tuple_element<
924 INDEX, std::tuple<VARTYPES...>>::type::type),
925 data.template get<INDEX>(), value);
926 } else {
927 // But vectors and jagged vectors do.
928 memset(data.template get<INDEX>(), value);
929 }
930 // Recurse into the next variable.
931 if constexpr (sizeof...(VARTYPES) > (INDEX + 1)) {
932 memset_impl<INDEX + 1>(data, value);
933 }
934}
935
936template <std::size_t INDEX, typename... VARTYPES,
937 template <typename> class INTERFACE>
938void copy::resize_impl(
939 const edm::view<edm::details::add_const_t<edm::schema<VARTYPES...>>>&
940 from_view,
941 edm::host<edm::schema<VARTYPES...>, INTERFACE>& to_vec,
942 [[maybe_unused]] type::copy_type cptype) const {
943
944 // The target is a host container, so the copy type can't be anything
945 // targeting a device.
946 assert((cptype == type::device_to_host) || (cptype == type::host_to_host) ||
947 (cptype == type::unknown));
948
949 // First, handle containers with no jagged vectors in them.
950 if constexpr (std::disjunction_v<
951 edm::type::details::is_jagged_vector<VARTYPES>...> ==
952 false) {
953 // The single size of the source container.
954 typename edm::view<
955 edm::details::add_const_t<edm::schema<VARTYPES...>>>::size_type
956 size = from_view.capacity();
957 // If the container is resizable, take its size.
958 if (from_view.size().ptr() != nullptr) {
959 // A small security check.
960 assert(from_view.size().size() ==
961 sizeof(typename edm::view<edm::details::add_const_t<
962 edm::schema<VARTYPES...>>>::size_type));
963 // Get the exact size of the container.
964 do_copy(sizeof(typename edm::view<edm::details::add_const_t<
965 edm::schema<VARTYPES...>>>::size_type),
966 from_view.size().ptr(), &size, cptype);
967 // We have to wait for this to finish, since the "size" variable is
968 // not going to be available outside of this function. And
969 // asynchronous SYCL memory copies can happen from variables on the
970 // stack as well...
971 create_event()->wait();
972 }
973 // Resize the target container.
974 VECMEM_DEBUG_MSG(4, "Resizing a (non-jagged) container to size %u",
975 size);
976 to_vec.resize(size);
977 } else {
978 // Resize vector and jagged vector variables one by one. Note that
979 // evaluation order matters here. Since all jagged vectors are vectors,
980 // but not all vectors are jagged vectors. ;-)
981 if constexpr (edm::type::details::is_jagged_vector<
982 typename std::tuple_element<
983 INDEX, std::tuple<VARTYPES...>>::type>::value) {
984 // Get the sizes of this jagged vector.
985 auto sizes = get_sizes(from_view.template get<INDEX>());
986 // Set the "outer size" of the jagged vector.
987 VECMEM_DEBUG_MSG(
988 4, "Resizing jagged vector variable at index %lu to size %lu",
989 INDEX, sizes.size());
990 details::resize_jagged_vector(to_vec.template get<INDEX>(),
991 sizes.size());
992 // Set the "inner sizes" of the jagged vector.
993 for (std::size_t i = 0; i < sizes.size(); ++i) {
994 to_vec.template get<INDEX>()[i].resize(sizes[i]);
995 }
996 } else if constexpr (edm::type::details::is_vector<
997 typename std::tuple_element<
998 INDEX,
999 std::tuple<VARTYPES...>>::type>::value) {
1000 // Get the size of this vector.
1001 auto size = get_size(from_view.template get<INDEX>());
1002 // Resize the target vector.
1003 VECMEM_DEBUG_MSG(4,
1004 "Resizing vector variable at index %lu to size %u",
1005 INDEX, size);
1006 to_vec.template get<INDEX>().resize(size);
1007 }
1008 // Call this function recursively.
1009 if constexpr (sizeof...(VARTYPES) > (INDEX + 1)) {
1010 resize_impl<INDEX + 1>(from_view, to_vec, cptype);
1011 }
1012 }
1013}
1014
1015template <std::size_t INDEX, typename... VARTYPES>
1016void copy::copy_sizes_impl(
1017 [[maybe_unused]] const edm::view<
1018 edm::details::add_const_t<edm::schema<VARTYPES...>>>& from_view,
1019 [[maybe_unused]] edm::view<edm::schema<VARTYPES...>> to_view,
1020 [[maybe_unused]] type::copy_type cptype) const {
1021
1022 // This should only be called for a resizable target container, with
1023 // a non-resizable source container.
1024 assert(to_view.size().ptr() != nullptr);
1025 assert(from_view.size().ptr() == nullptr);
1026
1027 // First, handle containers with no jagged vectors in them.
1028 if constexpr (std::disjunction_v<
1029 edm::type::details::is_jagged_vector<VARTYPES>...> ==
1030 false) {
1031 // Size of the source container.
1032 typename edm::view<
1033 edm::details::add_const_t<edm::schema<VARTYPES...>>>::size_type
1034 size = from_view.capacity();
1035 // Choose the copy type.
1036 type::copy_type size_cptype = type::unknown;
1037 switch (cptype) {
1040 size_cptype = type::host_to_device;
1041 break;
1042 case type::host_to_host:
1044 size_cptype = type::host_to_host;
1045 break;
1046 default:
1047 break;
1048 }
1049 // Set the size of the target container.
1050 do_copy(sizeof(typename edm::view<edm::details::add_const_t<
1051 edm::schema<VARTYPES...>>>::size_type),
1052 &size, to_view.size().ptr(), size_cptype);
1053 // We have to wait for this to finish, since the "size" variable is not
1054 // going to be available outside of this function. And asynchronous SYCL
1055 // memory copies can happen from variables on the stack as well...
1056 create_event()->wait();
1057 } else {
1058 // For the jagged vector case we recursively copy the sizes of every
1059 // jagged vector variable. The rest of the variables are not resizable
1060 // in such a container, so they are ignored here.
1061 if constexpr (edm::type::details::is_jagged_vector<
1062 typename std::tuple_element<
1063 INDEX, std::tuple<VARTYPES...>>::type>::value) {
1064 // Copy the sizes for this variable.
1065 const auto sizes = get_sizes(from_view.template get<INDEX>());
1066 set_sizes(sizes, to_view.template get<INDEX>())->wait();
1067 }
1068 // Call this function recursively.
1069 if constexpr (sizeof...(VARTYPES) > (INDEX + 1)) {
1070 copy_sizes_impl<INDEX + 1>(from_view, to_view, cptype);
1071 }
1072 }
1073}
1074
1075template <std::size_t INDEX, typename... VARTYPES>
1076void copy::copy_payload_impl(
1077 const edm::view<edm::details::add_const_t<edm::schema<VARTYPES...>>>&
1078 from_view,
1079 edm::view<edm::schema<VARTYPES...>> to_view, type::copy_type cptype) const {
1080
1081 // Scalars do not have their own dedicated @c copy functions.
1082 if constexpr (edm::type::details::is_scalar<typename std::tuple_element<
1083 INDEX, std::tuple<VARTYPES...>>::type>::value) {
1084 do_copy(sizeof(typename std::tuple_element<
1085 INDEX, std::tuple<VARTYPES...>>::type::type),
1086 from_view.template get<INDEX>(), to_view.template get<INDEX>(),
1087 cptype);
1088 } else {
1089 // But vectors and jagged vectors do.
1090 copy_view_impl(from_view.template get<INDEX>(),
1091 to_view.template get<INDEX>(), cptype);
1092 }
1093 // Recurse into the next variable.
1094 if constexpr (sizeof...(VARTYPES) > (INDEX + 1)) {
1095 copy_payload_impl<INDEX + 1>(from_view, to_view, cptype);
1096 }
1097}
1098
1099template <std::size_t INDEX, typename... VARTYPES>
1100std::vector<data::vector_view<int>::size_type> copy::get_sizes_impl(
1101 const edm::view<edm::schema<VARTYPES...>>& view) const {
1102
1103 // Make sure that there's a jagged vector in here.
1104 static_assert(
1105 edm::details::has_jagged_vector<edm::schema<VARTYPES...>>::value,
1106 "Function can only be used on containers with jagged vectors");
1107
1108 // Check if the variable at INDEX is a jagged vector.
1109 if constexpr (edm::type::details::is_jagged_vector<
1110 typename std::tuple_element<
1111 INDEX, std::tuple<VARTYPES...>>::type>::value) {
1112 // If it is, get its (inner) sizes.
1113 return get_sizes(view.template get<INDEX>());
1114 } else if constexpr (sizeof...(VARTYPES) > (INDEX + 1)) {
1115 // If it's not, and this was not the last variable, recurse into the
1116 // next variable.
1117 return get_sizes_impl<INDEX + 1>(view);
1118 }
1119
1120 // We should never get here as long as the code was written correctly.
1121#if defined(__GNUC__)
1122 __builtin_unreachable();
1123#elif defined(_MSC_VER)
1124 __assume(0);
1125#endif
1126}
1127
1128template <std::size_t INDEX, typename... VARTYPES>
1129async_sizes<data::vector_view<int>::size_type> copy::get_sizes_impl(
1130 const edm::view<edm::schema<VARTYPES...>>& view,
1131 memory_resource& pinnedHostMr) const {
1132
1133 // Make sure that there's a jagged vector in here.
1134 static_assert(
1135 edm::details::has_jagged_vector<edm::schema<VARTYPES...>>::value,
1136 "Function can only be used on containers with jagged vectors");
1137
1138 // Check if the variable at INDEX is a jagged vector.
1139 if constexpr (edm::type::details::is_jagged_vector<
1140 typename std::tuple_element<
1141 INDEX, std::tuple<VARTYPES...>>::type>::value) {
1142 // If it is, get its (inner) sizes.
1143 return get_sizes(view.template get<INDEX>(), pinnedHostMr);
1144 } else if constexpr (sizeof...(VARTYPES) > (INDEX + 1)) {
1145 // If it's not, and this was not the last variable, recurse into the
1146 // next variable.
1147 return get_sizes_impl<INDEX + 1>(view, pinnedHostMr);
1148 }
1149
1150 // We should never get here as long as the code was written correctly.
1151#if defined(__GNUC__)
1152 __builtin_unreachable();
1153#elif defined(_MSC_VER)
1154 __assume(0);
1155#endif
1156}
1157
1158} // namespace vecmem
An allocator class that wraps a memory resource.
Definition allocator.hpp:37
Return type for asynchronous size retrievals.
Definition async_size.hpp:25
VECMEM_NODISCARD event_type set_sizes(const std::vector< typename data::vector_view< TYPE >::size_type > &sizes, data::jagged_vector_view< TYPE > data) const
Helper function for setting the sizes of a resizable jagged vector.
VECMEM_NODISCARD event_type operator()(const data::vector_view< std::add_const_t< TYPE > > &from, data::vector_view< TYPE > to, type::copy_type cptype=type::unknown) const
Copy a 1-dimensional vector's data between two existing memory blocks.
VECMEM_NODISCARD event_type setup(data::vector_view< TYPE > data) const
Set up the internal state of a vector buffer correctly on a device.
virtual void do_copy(std::size_t size, const void *from, void *to, type::copy_type cptype) const
Perform a "low level" memory copy.
Definition copy.cpp:31
std::unique_ptr< abstract_event > event_type
Event type used by the copy class.
Definition copy.hpp:73
data::vector_view< TYPE >::size_type get_size(const data::vector_view< TYPE > &data) const
Helper function for getting the size of a resizable 1D buffer.
Definition copy.ipp:121
virtual void do_memset(std::size_t size, void *ptr, int value) const
Perform a "low level" memory filling operation.
Definition copy.cpp:44
virtual VECMEM_NODISCARD event_type create_event() const
Create an event for synchronization.
Definition copy.cpp:54
std::vector< typename data::vector_view< TYPE >::size_type > get_sizes(const data::jagged_vector_view< TYPE > &data) const
Helper function for getting the sizes of a resizable jagged vector.
Definition copy.ipp:307
VECMEM_NODISCARD event_type memset(data::vector_view< TYPE > data, int value) const
Set all bytes of the vector to some value.
data::vector_buffer< std::remove_cv_t< TYPE > > to(const data::vector_view< TYPE > &data, memory_resource &resource, type::copy_type cptype=type::unknown) const
Copy a 1-dimensional vector to the specified memory resource.
Definition copy.ipp:66
Object owning all the data of a jagged vector.
Definition jagged_vector_buffer.hpp:30
typename base_type::value_type value_type
Use the base class's value_type.
Definition jagged_vector_buffer.hpp:38
A view for jagged vectors.
Definition jagged_vector_view.hpp:45
VECMEM_HOST_AND_DEVICE pointer host_ptr() const
Access the host accessible array describing the inner vectors.
Definition jagged_vector_view.ipp:105
std::size_t size_type
We cannot use boolean types.
Definition jagged_vector_view.hpp:53
VECMEM_HOST_AND_DEVICE size_type capacity() const
Get the maximum capacity of the "outer" vector.
Definition jagged_vector_view.ipp:92
VECMEM_HOST_AND_DEVICE size_type size() const
Get the "outer" size of the jagged vector.
Definition jagged_vector_view.ipp:86
Object owning the data held by it.
Definition vector_buffer.hpp:29
typename base_type::size_type size_type
Size type definition coming from the base class.
Definition vector_buffer.hpp:35
Class holding data about a 1 dimensional vector/array.
Definition vector_view.hpp:38
unsigned int size_type
We cannot use boolean types.
Definition vector_view.hpp:47
VECMEM_HOST_AND_DEVICE size_pointer size_ptr() const
Get a pointer to the size of the vector.
Definition vector_view.ipp:84
VECMEM_HOST_AND_DEVICE size_type capacity() const
Get the maximum capacity of the vector.
Definition vector_view.ipp:78
Technical base type for buffer<schema<VARTYPES...>>
Definition buffer.hpp:28
Technical base type for view<schema<VARTYPES...>>
Definition view.hpp:29
VECMEM_HOST std::vector< typename vector_view< T >::size_type > get_capacities(const jagged_vector_view< T > &data)
Get the capacities of the inner vectors of a jagged vector.
Definition jagged_vector_view.ipp:111
buffer_type
"Overall type" for a buffer object
Definition buffer_type.hpp:13
@ fixed_size
The buffer has a fixed number of elements.
@ resizable
The buffer is resizable/expandable.
void resize_jagged_vector(std::vector< std::vector< T, ALLOC1 >, ALLOC2 > &vec, std::size_t size)
Resize a generic jagged vector.
Definition resize_jagged_vector.hpp:23
VECMEM_HOST std::vector< vecmem::data::vector_view< int >::size_type > get_capacities(const view< schema< VARTYPES... > > &soa)
Helper function to get the capacities of the jagged vectors in a view.
Definition view.ipp:198
Main namespace for the vecmem classes/functions.
Definition atomic_ref.hpp:16
VECMEM_HOST data::vector_view< T > get_data(array< T, N > &a)
Helper function creating a vecmem::data::vector_view object.
Definition array.ipp:217
copy_type
Types of memory copies to handle.
Definition copy.hpp:56
@ device_to_host
Copy operation between a device and the host.
Definition copy.hpp:60
@ device_to_device
Copy operation between two devices.
Definition copy.hpp:64
@ unknown
Unknown copy type, determined at runtime.
Definition copy.hpp:66
@ host_to_host
Copy operation on the host.
Definition copy.hpp:62
@ host_to_device
Copy operation between the host and a device.
Definition copy.hpp:58
Technical base type for has_jagged_vector<schema<VARTYPES...>>
Definition schema_traits.hpp:238
Meta type describing the "schema" of an SoA container.
Definition schema.hpp:46
Definition schema_traits.hpp:108