vecmem 1.18.0
Loading...
Searching...
No Matches
copy.ipp
1/*
2 * VecMem project, part of the ACTS project (R&D line)
3 *
4 * (c) 2021-2025 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/host_memory_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
147 // Check if anything needs to be done.
148 if (data.size() == 0) {
150 }
151
152 // "Set up" the inner vector descriptors, using the host-accessible data.
153 // But only if the jagged vector buffer is resizable.
154 if (data.host_ptr()[0].size_ptr() != nullptr) {
155 do_memset(
156 sizeof(typename data::vector_buffer<TYPE>::size_type) * data.size(),
157 data.host_ptr()[0].size_ptr(), 0);
158 }
159
160 // Check if anything else needs to be done.
161 if (data.ptr() == data.host_ptr()) {
162 return create_event();
163 }
164
165 // Copy the description of the inner vectors of the buffer.
166 do_copy(
167 data.size() *
168 sizeof(
170 data.host_ptr(), data.ptr(), type::host_to_device);
171 VECMEM_DEBUG_MSG(2,
172 "Prepared a jagged device vector buffer of size %lu "
173 "for use on a device",
174 data.size());
175
176 // Return a new event.
177 return create_event();
178}
179
180template <typename TYPE>
181copy::event_type copy::memset(data::jagged_vector_view<TYPE> data,
182 int value) const {
183
184 // Do different things for jagged vectors that are contiguous in memory.
185 if (is_contiguous(data.host_ptr(), data.capacity())) {
186 // For contiguous jagged vectors we can just do a single memset
187 // operation. First, let's calculate the "total size" of the jagged
188 // vector.
189 const std::size_t total_size = std::accumulate(
190 data.host_ptr(), data.host_ptr() + data.size(),
191 static_cast<std::size_t>(0u),
192 [](std::size_t sum, const data::vector_view<TYPE>& iv) {
193 return sum + iv.capacity();
194 });
195 // Find the first "inner vector" that has a non-zero capacity.
196 for (std::size_t i = 0; i < data.size(); ++i) {
197 data::vector_view<TYPE>& iv = data.host_ptr()[i];
198 if ((iv.capacity() != 0u) && (iv.ptr() != nullptr)) {
199 // Call memset with its help.
200 do_memset(total_size * sizeof(TYPE), iv.ptr(), value);
201 return create_event();
202 }
203 }
204 // If we are still here, apparently we didn't need to do anything.
206 } else {
207 // For non-contiguous jagged vectors call memset one-by-one on the
208 // inner vectors. Note that we don't use vecmem::copy::memset here,
209 // as that would require us to wait for each memset individually.
210 for (std::size_t i = 0; i < data.size(); ++i) {
211 data::vector_view<TYPE>& iv = data.host_ptr()[i];
212 do_memset(iv.capacity() * sizeof(TYPE), iv.ptr(), value);
213 }
214 }
215
216 // Return a new event.
217 return create_event();
218}
219
220template <typename TYPE>
222 const data::jagged_vector_view<TYPE>& data, memory_resource& resource,
223 memory_resource* host_access_resource, type::copy_type cptype) const {
224
225 // Create the result buffer object.
227 (((data.capacity() > 0u) && (data.host_ptr()[0].size_ptr() != nullptr))
232 assert(result.size() == data.size());
233
234 // Copy the description of the "inner vectors" if necessary.
235 if (host_access_resource != nullptr) {
236 setup(result)->wait();
237 }
238
239 // Copy the payload of the inner vectors. Explicitly waiting for the copy to
240 // finish before returning the buffer.
241 operator()(data, result, cptype)->wait();
242
243 // Return the newly created object.
244 return result;
245}
246
247template <typename TYPE>
249 const data::jagged_vector_view<std::add_const_t<TYPE>>& from_view,
250 data::jagged_vector_view<TYPE> to_view, type::copy_type cptype) const {
251
252 // Perform the copy. Depending on whether an actual copy has been set up /
253 // performed, return either "an actual", or just a dummy event.
254 if (copy_view_impl(from_view, to_view, cptype)) {
255 return create_event();
256 } else {
258 }
259}
260
261template <typename TYPE, typename ALLOC1, typename ALLOC2>
263 const data::jagged_vector_view<std::add_const_t<TYPE>>& from_view,
264 std::vector<std::vector<TYPE, ALLOC2>, ALLOC1>& to_vec,
265 type::copy_type cptype) const {
266
267 // Resize the output object to the correct size.
268 details::resize_jagged_vector(to_vec, from_view.size());
269 const auto sizes = get_sizes(from_view);
270 assert(sizes.size() == to_vec.size());
271 for (typename data::jagged_vector_view<std::add_const_t<TYPE>>::size_type
272 i = 0;
273 i < from_view.size(); ++i) {
274 to_vec[i].resize(sizes[i]);
275 }
276
277 // Perform the memory copy.
278 return operator()(from_view, vecmem::get_data(to_vec), cptype);
279}
280
281template <typename TYPE>
282std::vector<typename data::vector_view<TYPE>::size_type> copy::get_sizes(
283 const data::jagged_vector_view<TYPE>& data) const {
284
285 // Perform the operation using the private function.
286 return get_sizes_impl(data.host_ptr(), data.size());
287}
288
289template <typename TYPE>
291 const std::vector<typename data::vector_view<TYPE>::size_type>& sizes,
293
294 // Finish early if possible.
295 if ((sizes.size() == 0) && (data.size() == 0)) {
297 }
298 // Make sure that the sizes match up.
299 if (sizes.size() != data.size()) {
300 std::ostringstream msg;
301 msg << "sizes.size() (" << sizes.size() << ") != data.size() ("
302 << data.size() << ")";
303 throw std::length_error(msg.str());
304 }
305 // Make sure that the target jagged vector is either resizable, or it has
306 // the correct sizes/capacities already.
307 bool perform_copy = true;
309 i < data.size(); ++i) {
310 // Perform a capacity check in all cases. Whether or not the target is
311 // resizable.
312 if (data.host_ptr()[i].capacity() < sizes[i]) {
313 std::ostringstream msg;
314 msg << "data.host_ptr()[" << i << "].capacity() ("
315 << data.host_ptr()[i].capacity() << ") < sizes[" << i << "] ("
316 << sizes[i] << ")";
317 throw std::length_error(msg.str());
318 }
319 // Make sure that the inner vectors are consistently either all
320 // resizable, or none of them are.
321 if (data.host_ptr()[i].size_ptr() == nullptr) {
322 perform_copy = false;
323 } else if (perform_copy == false) {
324 throw std::invalid_argument(
325 "Inconsistent target jagged vector view received for resizing");
326 }
327 }
328 // If no copy is necessary, we're done.
329 if (perform_copy == false) {
331 }
332 // Perform the copy with some internal knowledge of how resizable jagged
333 // vector buffers work.
334 do_copy(sizeof(typename data::vector_view<TYPE>::size_type) * sizes.size(),
335 sizes.data(), data.host_ptr()->size_ptr(), type::unknown);
336
337 // Return a new event.
338 return create_event();
339}
340
341template <typename SCHEMA>
342copy::event_type copy::setup(edm::view<SCHEMA> data) const {
343
344 // For empty containers nothing needs to be done.
345 if (data.capacity() == 0) {
347 }
348
349 // Copy the data layout to the device, if needed.
350 if (data.layout().ptr() != data.host_layout().ptr()) {
351 assert(data.layout().capacity() > 0u);
352 [[maybe_unused]] bool did_copy =
353 copy_view_impl(data.host_layout(), data.layout(), type::unknown);
354 assert(did_copy);
355 }
356
357 // Initialize the "size variable(s)" correctly on the buffer.
358 if (data.size().ptr() != nullptr) {
359 assert(data.size().capacity() > 0u);
360 do_memset(data.size().capacity() * sizeof(char), data.size().ptr(), 0);
361 }
362 VECMEM_DEBUG_MSG(3,
363 "Prepared an SoA container of capacity %u "
364 "for use on a device (layout: {%u, %p}, size: {%u, %p})",
365 data.capacity(), data.layout().size(),
366 static_cast<void*>(data.layout().ptr()),
367 data.size().size(), static_cast<void*>(data.size().ptr()));
368
369 // Return a new event.
370 return create_event();
371}
372
373template <typename... VARTYPES>
374copy::event_type copy::memset(edm::view<edm::schema<VARTYPES...>> data,
375 int value) const {
376
377 // For buffers, we can do this in one go.
378 if (data.payload().ptr() != nullptr) {
379 memset(data.payload(), value);
380 } else {
381 // Do the operation using the helper function, recursively.
382 memset_impl<0>(data, value);
383 }
384
385 // Return a new event.
386 return create_event();
387}
388
389template <typename... VARTYPES>
390edm::buffer<edm::details::remove_cv_t<edm::schema<VARTYPES...>>> copy::to(
391 const edm::view<edm::schema<VARTYPES...>>& data, memory_resource& resource,
392 [[maybe_unused]] memory_resource* host_access_resource,
393 type::copy_type cptype) const {
394
395 // Create the result buffer object.
397 ((data.size().ptr() != nullptr) ? data::buffer_type::resizable
399 edm::buffer<edm::details::remove_cv_t<edm::schema<VARTYPES...>>> result;
401 edm::schema<VARTYPES...>>::value) {
402 // Set up a buffer that has (a) jagged vector(s).
404 btype};
405 // The idea here is that since we are going to copy data into the newly
406 // created buffer, the only thing that setup(...) needs to do is to set
407 // up "the layout" in non-host-accessible memory. (Zeroing the size
408 // variable(s) is not needed.) When doing a device-to-host copy, copying
409 // "the layout" is not actually needed. Worse yet, the copy object is
410 // not capable of memsetting the size variable(s). Since the device
411 // specific copy object cannot memset host memory.
412 //
413 // Long story short, the entire setup is skipped in the absence of a
414 // host-accessible memory resource.
415 if (host_access_resource != nullptr) {
416 setup(result)->wait();
417 }
418 } else {
419 // Set up a buffer absent of jagged vectors.
420 result = {data.capacity(), resource, btype};
421 // In the absence of jagged vectors there is no "layout" to copy. And
422 // since the size of the buffer will be set during the copy correctly,
423 // there is nothing else that setup(...) would need to do. So no need
424 // to call it here.
425 }
426
427 // Perform the copy.
428 operator()(data, result, cptype)->wait();
429
430 // Return the buffer.
431 return result;
432}
433
434template <typename... VARTYPES>
436 const edm::view<edm::details::add_const_t<edm::schema<VARTYPES...>>>&
437 from_view,
438 edm::view<edm::schema<VARTYPES...>> to_view, type::copy_type cptype) const {
439
440 // First, handle the simple case, when both views have a contiguous memory
441 // layout.
442 if ((from_view.payload().ptr() != nullptr) &&
443 (to_view.payload().ptr() != nullptr) &&
444 (from_view.payload().capacity() == to_view.payload().capacity())) {
445
446 // If the "common size" is zero, we're done.
447 if (from_view.payload().capacity() == 0) {
449 }
450
451 // Let the user know what's happening.
452 VECMEM_DEBUG_MSG(2, "Performing simple SoA copy of %u bytes",
453 from_view.payload().size());
454
455 // Copy the payload with a single copy operation.
456 copy_view_impl(from_view.payload(), to_view.payload(), cptype);
457
458 // If the target view is resizable, set its size.
459 if (to_view.size().ptr() != nullptr) {
460 // If the source is also resizable, the situation should be simple.
461 if (from_view.size().ptr() != nullptr) {
462 // Check that the sizes are the same.
463 if (from_view.size().capacity() != to_view.size().capacity()) {
464 std::ostringstream msg;
465 msg << "from_view.size().capacity() ("
466 << from_view.size().capacity()
467 << ") != to_view.size().capacity() ("
468 << to_view.size().capacity() << ")";
469 throw std::length_error(msg.str());
470 }
471 // Perform a dumb copy.
472 copy_view_impl(from_view.size(), to_view.size(), cptype);
473 } else {
474 // If not, then copy the size(s) recursively.
475 copy_sizes_impl<0>(from_view, to_view, cptype);
476 }
477 }
478
479 // Create a synchronization event.
480 return create_event();
481 }
482
483 // If not, then do an un-optimized copy, variable-by-variable.
484 copy_payload_impl<0>(from_view, to_view, cptype);
485
486 // Return a new event.
487 return create_event();
488}
489
490template <typename... VARTYPES, template <typename> class INTERFACE>
492 const edm::view<edm::details::add_const_t<edm::schema<VARTYPES...>>>&
493 from_view,
494 edm::host<edm::schema<VARTYPES...>, INTERFACE>& to_vec,
495 type::copy_type cptype) const {
496
497 // Resize the output object to the correct size(s).
498 resize_impl<0>(from_view, to_vec, cptype);
499
500 // Perform the memory copy.
501 return operator()(from_view, vecmem::get_data(to_vec), cptype);
502}
503
504template <typename... VARTYPES>
506 const edm::view<edm::schema<VARTYPES...>>& data) const {
507
508 // Start by taking the capacity of the container.
509 typename edm::view<edm::schema<VARTYPES...>>::size_type size =
510 data.capacity();
511
512 // If there are jagged vectors in the container and/or the container is not
513 // resizable, we're done. It is done in two separate if statements to avoid
514 // MSVC trying to be too smart, and giving a warning...
515 if constexpr (std::disjunction_v<
517 return size;
518 } else {
519 // All the rest is put into an else block, to avoid MSVC trying to be
520 // too smart, and giving a warning about unreachable code...
521 if (data.size().ptr() == nullptr) {
522 return size;
523 }
524
525 // A small security check.
526 assert(data.size().size() ==
527 sizeof(typename edm::view<edm::schema<VARTYPES...>>::size_type));
528
529 // Get the exact size of the container.
530 do_copy(sizeof(typename edm::view<edm::schema<VARTYPES...>>::size_type),
531 data.size().ptr(), &size, type::unknown);
532 // We have to wait for this to finish, since the "size" variable is
533 // not going to be available outside of this function. And
534 // asynchronous SYCL memory copies can happen from variables on the
535 // stack as well...
536 create_event()->wait();
537
538 // Return what we got.
539 return size;
540 }
541}
542
543template <typename... VARTYPES>
544std::vector<data::vector_view<int>::size_type> copy::get_sizes(
545 const edm::view<edm::schema<VARTYPES...>>& data) const {
546
547 // Make sure that there's a jagged vector in here.
548 static_assert(
550 "Function can only be used on containers with jagged vectors");
551
552 // Get the sizes using the helper function.
553 return get_sizes_impl<0>(data);
554}
555
556template <typename TYPE>
557bool copy::copy_view_impl(
558 const data::vector_view<std::add_const_t<TYPE>>& from_view,
559 data::vector_view<TYPE> to_view, type::copy_type cptype) const {
560
561 // Get the size of the source view.
562 const typename data::vector_view<std::add_const_t<TYPE>>::size_type size =
564 // If it's zero, don't do an actual copy.
565 if (size == 0u) {
566 return false;
567 }
568
569 // Make sure that the copy can happen.
570 if (to_view.capacity() < size) {
571 std::ostringstream msg;
572 msg << "Target capacity (" << to_view.capacity() << ") < source size ("
573 << size << ")";
574 throw std::length_error(msg.str());
575 }
576
577 // Make sure that if the target view is resizable, that it would be set up
578 // for the correct size.
579 if (to_view.size_ptr() != nullptr) {
580 // Select what type of copy this should be. Keeping in mind that we copy
581 // from a variable on the host stack. So the question is just whether
582 // the target is the host, or a device.
583 type::copy_type size_cptype = type::unknown;
584 switch (cptype) {
587 size_cptype = type::host_to_device;
588 break;
591 size_cptype = type::host_to_host;
592 break;
593 default:
594 break;
595 }
596 // Perform the copy.
597 do_copy(sizeof(typename data::vector_view<TYPE>::size_type), &size,
598 to_view.size_ptr(), size_cptype);
599 // We have to wait for this to finish, since the "size" variable is not
600 // going to be available outside of this function. And asynchronous SYCL
601 // memory copies can happen from variables on the stack as well...
602 create_event()->wait();
603 }
604
605 // Copy the payload.
606 do_copy(size * sizeof(TYPE), from_view.ptr(), to_view.ptr(), cptype);
607 return true;
608}
609
610template <typename TYPE>
611bool copy::copy_view_impl(
612 const data::jagged_vector_view<std::add_const_t<TYPE>>& from_view,
613 data::jagged_vector_view<TYPE> to_view, type::copy_type cptype) const {
614
615 // Sanity checks.
616 if (from_view.size() > to_view.size()) {
617 std::ostringstream msg;
618 msg << "from_view.size() (" << from_view.size()
619 << ") > to_view.size() (" << to_view.size() << ")";
620 throw std::length_error(msg.str());
621 }
622
623 // Get the "outer size" of the container.
624 const typename data::jagged_vector_view<std::add_const_t<TYPE>>::size_type
625 size = from_view.size();
626 if (size == 0u) {
627 return false;
628 }
629
630 // Calculate the contiguous-ness of the memory allocations.
631 const bool from_is_contiguous = is_contiguous(from_view.host_ptr(), size);
632 const bool to_is_contiguous = is_contiguous(to_view.host_ptr(), size);
633 VECMEM_DEBUG_MSG(3, "from_is_contiguous = %d, to_is_contiguous = %d",
634 from_is_contiguous, to_is_contiguous);
635
636 // Get the sizes of the source jagged vector.
637 const auto sizes = get_sizes(from_view);
638
639 // Before even attempting the copy, make sure that the target view either
640 // has the correct sizes, or can be resized correctly. Captire the event
641 // marking the end of the operation. We need to wait for the copy to finish
642 // before returning from this function.
643 auto set_sizes_event = set_sizes(sizes, to_view);
644
645 // Check whether the source and target capacities match up. We can only
646 // perform the "optimised copy" if they do.
647 std::vector<typename data::vector_view<std::add_const_t<TYPE>>::size_type>
648 capacities(size);
649 bool capacities_match = true;
650 for (std::size_t i = 0; i < size; ++i) {
651 if (from_view.host_ptr()[i].capacity() !=
652 to_view.host_ptr()[i].capacity()) {
653 capacities_match = false;
654 break;
655 }
656 capacities[i] = from_view.host_ptr()[i].capacity();
657 }
658
659 // Perform the copy as best as we can.
660 if (from_is_contiguous && to_is_contiguous && capacities_match) {
661 // Perform the copy in one go.
662 copy_views_contiguous_impl(capacities, from_view.host_ptr(),
663 to_view.host_ptr(), cptype);
664 } else {
665 // Do the copy as best as we can. Note that since they are not
666 // contiguous anyway, we use the sizes of the vectors here, not their
667 // capcities.
668 copy_views_impl(sizes, from_view.host_ptr(), to_view.host_ptr(),
669 cptype);
670 }
671
672 // Before returning, make sure that the "size set" operation has finished.
673 // Because the "sizes" variable is just about to go out of scope.
674 set_sizes_event->wait();
675 return true;
676}
677
678template <typename TYPE>
679void copy::copy_views_impl(
680 const std::vector<typename data::vector_view<TYPE>::size_type>& sizes,
681 const data::vector_view<std::add_const_t<TYPE>>* from_view,
682 data::vector_view<TYPE>* to_view, type::copy_type cptype) const {
683
684 // Some security checks.
685 assert(from_view != nullptr);
686 assert(to_view != nullptr);
687
688 // Helper variable(s) used in the copy.
689 const std::size_t size = sizes.size();
690 [[maybe_unused]] std::size_t copy_ops = 0;
691
692 // Perform the copy in multiple steps.
693 for (std::size_t i = 0; i < size; ++i) {
694
695 // Skip empty "inner vectors".
696 if (sizes[i] == 0) {
697 continue;
698 }
699
700 // Some sanity checks.
701 assert(from_view[i].ptr() != nullptr);
702 assert(to_view[i].ptr() != nullptr);
703 assert(sizes[i] <= from_view[i].capacity());
704 assert(sizes[i] <= to_view[i].capacity());
705
706 // Perform the copy.
707 do_copy(sizes[i] * sizeof(TYPE), from_view[i].ptr(), to_view[i].ptr(),
708 cptype);
709 ++copy_ops;
710 }
711
712 // Let the user know what happened.
713 VECMEM_DEBUG_MSG(2,
714 "Copied the payload of a jagged vector of type "
715 "\"%s\" with %lu copy operation(s)",
716 typeid(TYPE).name(), copy_ops);
717}
718
719template <typename TYPE>
720void copy::copy_views_contiguous_impl(
721 const std::vector<typename data::vector_view<TYPE>::size_type>& sizes,
722 const data::vector_view<std::add_const_t<TYPE>>* from_view,
723 data::vector_view<TYPE>* to_view, type::copy_type cptype) const {
724
725 // Some security checks.
726 assert(from_view != nullptr);
727 assert(to_view != nullptr);
728 assert(is_contiguous(from_view, sizes.size()));
729 assert(is_contiguous(to_view, sizes.size()));
730
731 // Helper variable(s) used in the copy.
732 const std::size_t size = sizes.size();
733 const std::size_t total_size =
734 std::accumulate(sizes.begin(), sizes.end(),
735 static_cast<std::size_t>(0)) *
736 sizeof(TYPE);
737
738 // Find the first non-empty element.
739 for (std::size_t i = 0; i < size; ++i) {
740
741 // Jump over empty elements.
742 if (sizes[i] == 0) {
743 continue;
744 }
745
746 // Some sanity checks.
747 assert(from_view[i].ptr() != nullptr);
748 assert(to_view[i].ptr() != nullptr);
749
750 // Perform the copy.
751 do_copy(total_size, from_view[i].ptr(), to_view[i].ptr(), cptype);
752 break;
753 }
754
755 // Let the user know what happened.
756 VECMEM_DEBUG_MSG(2,
757 "Copied the payload of a jagged vector of type "
758 "\"%s\" with 1 copy operation(s)",
759 typeid(TYPE).name());
760}
761
762template <typename TYPE>
763std::vector<typename data::vector_view<TYPE>::size_type> copy::get_sizes_impl(
764 const data::vector_view<TYPE>* data, std::size_t size) const {
765
766 // Create the result vector.
767 std::vector<typename data::vector_view<TYPE>::size_type> result(size, 0);
768
769 // Try to get the "resizable sizes" first.
770 for (std::size_t i = 0; i < size; ++i) {
771 // Find the first "inner vector" that has a non-zero capacity, and is
772 // resizable.
773 if ((data[i].capacity() != 0) && (data[i].size_ptr() != nullptr)) {
774 // Copy the sizes of the inner vectors into the result vector.
776 (size - i),
777 data[i].size_ptr(), result.data() + i, type::unknown);
778 // Wait for the copy operation to finish. With some backends
779 // (khm... SYCL... khm...) copies can be asynchronous even into
780 // non-pinned host memory.
781 create_event()->wait();
782 // At this point the result vector should have been set up
783 // correctly.
784 return result;
785 }
786 }
787
788 // If we're still here, then the buffer is not resizable. So let's just
789 // collect the capacity of each of the inner vectors.
790 for (std::size_t i = 0; i < size; ++i) {
791 result[i] = data[i].capacity();
792 }
793 return result;
794}
795
796template <typename TYPE>
797bool copy::is_contiguous(const data::vector_view<TYPE>* data,
798 std::size_t size) {
799
800 // We should never call this function for an empty jagged vector.
801 assert(size > 0);
802
803 // Check whether all memory blocks are contiguous.
804 auto ptr = data[0].ptr();
805 for (std::size_t i = 1; i < size; ++i) {
806 if ((ptr + data[i - 1].capacity()) != data[i].ptr()) {
807 return false;
808 }
809 ptr = data[i].ptr();
810 }
811 return true;
812}
813
814template <std::size_t INDEX, typename... VARTYPES>
815void copy::memset_impl(edm::view<edm::schema<VARTYPES...>> data,
816 int value) const {
817
818 // Scalars do not have their own dedicated @c memset functions.
819 if constexpr (edm::type::details::is_scalar<typename std::tuple_element<
820 INDEX, std::tuple<VARTYPES...>>::type>::value) {
821 do_memset(sizeof(typename std::tuple_element<
822 INDEX, std::tuple<VARTYPES...>>::type::type),
823 data.template get<INDEX>(), value);
824 } else {
825 // But vectors and jagged vectors do.
826 memset(data.template get<INDEX>(), value);
827 }
828 // Recurse into the next variable.
829 if constexpr (sizeof...(VARTYPES) > (INDEX + 1)) {
830 memset_impl<INDEX + 1>(data, value);
831 }
832}
833
834template <std::size_t INDEX, typename... VARTYPES,
835 template <typename> class INTERFACE>
836void copy::resize_impl(
837 const edm::view<edm::details::add_const_t<edm::schema<VARTYPES...>>>&
838 from_view,
839 edm::host<edm::schema<VARTYPES...>, INTERFACE>& to_vec,
840 [[maybe_unused]] type::copy_type cptype) const {
841
842 // The target is a host container, so the copy type can't be anything
843 // targeting a device.
844 assert((cptype == type::device_to_host) || (cptype == type::host_to_host) ||
845 (cptype == type::unknown));
846
847 // First, handle containers with no jagged vectors in them.
848 if constexpr (std::disjunction_v<
849 edm::type::details::is_jagged_vector<VARTYPES>...> ==
850 false) {
851 // The single size of the source container.
852 typename edm::view<
853 edm::details::add_const_t<edm::schema<VARTYPES...>>>::size_type
854 size = from_view.capacity();
855 // If the container is resizable, take its size.
856 if (from_view.size().ptr() != nullptr) {
857 // A small security check.
858 assert(from_view.size().size() ==
859 sizeof(typename edm::view<edm::details::add_const_t<
860 edm::schema<VARTYPES...>>>::size_type));
861 // Get the exact size of the container.
862 do_copy(sizeof(typename edm::view<edm::details::add_const_t<
863 edm::schema<VARTYPES...>>>::size_type),
864 from_view.size().ptr(), &size, cptype);
865 // We have to wait for this to finish, since the "size" variable is
866 // not going to be available outside of this function. And
867 // asynchronous SYCL memory copies can happen from variables on the
868 // stack as well...
869 create_event()->wait();
870 }
871 // Resize the target container.
872 VECMEM_DEBUG_MSG(4, "Resizing a (non-jagged) container to size %u",
873 size);
874 to_vec.resize(size);
875 } else {
876 // Resize vector and jagged vector variables one by one. Note that
877 // evaluation order matters here. Since all jagged vectors are vectors,
878 // but not all vectors are jagged vectors. ;-)
879 if constexpr (edm::type::details::is_jagged_vector<
880 typename std::tuple_element<
881 INDEX, std::tuple<VARTYPES...>>::type>::value) {
882 // Get the sizes of this jagged vector.
883 auto sizes = get_sizes(from_view.template get<INDEX>());
884 // Set the "outer size" of the jagged vector.
885 VECMEM_DEBUG_MSG(
886 4, "Resizing jagged vector variable at index %lu to size %lu",
887 INDEX, sizes.size());
888 details::resize_jagged_vector(to_vec.template get<INDEX>(),
889 sizes.size());
890 // Set the "inner sizes" of the jagged vector.
891 for (std::size_t i = 0; i < sizes.size(); ++i) {
892 to_vec.template get<INDEX>()[i].resize(sizes[i]);
893 }
894 } else if constexpr (edm::type::details::is_vector<
895 typename std::tuple_element<
896 INDEX,
897 std::tuple<VARTYPES...>>::type>::value) {
898 // Get the size of this vector.
899 auto size = get_size(from_view.template get<INDEX>());
900 // Resize the target vector.
901 VECMEM_DEBUG_MSG(4,
902 "Resizing vector variable at index %lu to size %u",
903 INDEX, size);
904 to_vec.template get<INDEX>().resize(size);
905 }
906 // Call this function recursively.
907 if constexpr (sizeof...(VARTYPES) > (INDEX + 1)) {
908 resize_impl<INDEX + 1>(from_view, to_vec, cptype);
909 }
910 }
911}
912
913template <std::size_t INDEX, typename... VARTYPES>
914void copy::copy_sizes_impl(
915 [[maybe_unused]] const edm::view<
916 edm::details::add_const_t<edm::schema<VARTYPES...>>>& from_view,
917 [[maybe_unused]] edm::view<edm::schema<VARTYPES...>> to_view,
918 [[maybe_unused]] type::copy_type cptype) const {
919
920 // This should only be called for a resizable target container, with
921 // a non-resizable source container.
922 assert(to_view.size().ptr() != nullptr);
923 assert(from_view.size().ptr() == nullptr);
924
925 // First, handle containers with no jagged vectors in them.
926 if constexpr (std::disjunction_v<
927 edm::type::details::is_jagged_vector<VARTYPES>...> ==
928 false) {
929 // Size of the source container.
930 typename edm::view<
931 edm::details::add_const_t<edm::schema<VARTYPES...>>>::size_type
932 size = from_view.capacity();
933 // Choose the copy type.
934 type::copy_type size_cptype = type::unknown;
935 switch (cptype) {
938 size_cptype = type::host_to_device;
939 break;
942 size_cptype = type::host_to_host;
943 break;
944 default:
945 break;
946 }
947 // Set the size of the target container.
948 do_copy(sizeof(typename edm::view<edm::details::add_const_t<
949 edm::schema<VARTYPES...>>>::size_type),
950 &size, to_view.size().ptr(), size_cptype);
951 // We have to wait for this to finish, since the "size" variable is not
952 // going to be available outside of this function. And asynchronous SYCL
953 // memory copies can happen from variables on the stack as well...
954 create_event()->wait();
955 } else {
956 // For the jagged vector case we recursively copy the sizes of every
957 // jagged vector variable. The rest of the variables are not resizable
958 // in such a container, so they are ignored here.
959 if constexpr (edm::type::details::is_jagged_vector<
960 typename std::tuple_element<
961 INDEX, std::tuple<VARTYPES...>>::type>::value) {
962 // Copy the sizes for this variable.
963 const auto sizes = get_sizes(from_view.template get<INDEX>());
964 set_sizes(sizes, to_view.template get<INDEX>())->wait();
965 }
966 // Call this function recursively.
967 if constexpr (sizeof...(VARTYPES) > (INDEX + 1)) {
968 copy_sizes_impl<INDEX + 1>(from_view, to_view, cptype);
969 }
970 }
971}
972
973template <std::size_t INDEX, typename... VARTYPES>
974void copy::copy_payload_impl(
975 const edm::view<edm::details::add_const_t<edm::schema<VARTYPES...>>>&
976 from_view,
977 edm::view<edm::schema<VARTYPES...>> to_view, type::copy_type cptype) const {
978
979 // Scalars do not have their own dedicated @c copy functions.
980 if constexpr (edm::type::details::is_scalar<typename std::tuple_element<
981 INDEX, std::tuple<VARTYPES...>>::type>::value) {
982 do_copy(sizeof(typename std::tuple_element<
983 INDEX, std::tuple<VARTYPES...>>::type::type),
984 from_view.template get<INDEX>(), to_view.template get<INDEX>(),
985 cptype);
986 } else {
987 // But vectors and jagged vectors do.
988 copy_view_impl(from_view.template get<INDEX>(),
989 to_view.template get<INDEX>(), cptype);
990 }
991 // Recurse into the next variable.
992 if constexpr (sizeof...(VARTYPES) > (INDEX + 1)) {
993 copy_payload_impl<INDEX + 1>(from_view, to_view, cptype);
994 }
995}
996
997template <std::size_t INDEX, typename... VARTYPES>
998std::vector<data::vector_view<int>::size_type> copy::get_sizes_impl(
999 const edm::view<edm::schema<VARTYPES...>>& view) const {
1000
1001 // Make sure that there's a jagged vector in here.
1002 static_assert(
1003 edm::details::has_jagged_vector<edm::schema<VARTYPES...>>::value,
1004 "Function can only be used on containers with jagged vectors");
1005
1006 // Check if the variable at INDEX is a jagged vector.
1007 if constexpr (edm::type::details::is_jagged_vector<
1008 typename std::tuple_element<
1009 INDEX, std::tuple<VARTYPES...>>::type>::value) {
1010 // If it is, get its (inner) sizes.
1011 return get_sizes(view.template get<INDEX>());
1012 } else if constexpr (sizeof...(VARTYPES) > (INDEX + 1)) {
1013 // If it's not, and this was not the last variable, recurse into the
1014 // next variable.
1015 return get_sizes_impl<INDEX + 1>(view);
1016 }
1017
1018 // We should never get here as long as the code was written correctly.
1019#if defined(__GNUC__)
1020 __builtin_unreachable();
1021#elif defined(_MSC_VER)
1022 __assume(0);
1023#endif
1024}
1025
1026} // namespace vecmem
An allocator class that wraps a memory resource.
Definition allocator.hpp:37
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:70
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:282
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 pointer ptr() const
Get a pointer to the vector elements.
Definition jagged_vector_view.ipp:99
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
VECMEM_HOST_AND_DEVICE pointer ptr() const
Get a pointer to the vector elements.
Definition vector_view.ipp:91
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
VECMEM_HOST_AND_DEVICE size_type size() const
Get the size of the vector.
Definition vector_view.ipp:72
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:53
@ device_to_host
Copy operation between a device and the host.
Definition copy.hpp:57
@ device_to_device
Copy operation between two devices.
Definition copy.hpp:61
@ unknown
Unknown copy type, determined at runtime.
Definition copy.hpp:63
@ host_to_host
Copy operation on the host.
Definition copy.hpp:59
@ host_to_device
Copy operation between the host and a device.
Definition copy.hpp:55
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