vecmem 1.22.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 // Copy the data layout to the device, if needed.
345 if (data.layout().ptr() != data.host_layout().ptr()) {
346 assert(data.layout().capacity() > 0u);
347 [[maybe_unused]] bool did_copy =
348 copy_view_impl(data.host_layout(), data.layout(), type::unknown);
349 assert(did_copy);
350 }
351
352 // Initialize the "size variable(s)" correctly on the buffer.
353 if (data.size().ptr() != nullptr) {
354 assert(data.size().capacity() > 0u);
355 do_memset(data.size().capacity() * sizeof(char), data.size().ptr(), 0);
356 }
357 VECMEM_DEBUG_MSG(3,
358 "Prepared an SoA container of capacity %u "
359 "for use on a device (layout: {%u, %p}, size: {%u, %p})",
360 data.capacity(), data.layout().size(),
361 static_cast<void*>(data.layout().ptr()),
362 data.size().size(), static_cast<void*>(data.size().ptr()));
363
364 // Return a new event.
365 return create_event();
366}
367
368template <typename... VARTYPES>
369copy::event_type copy::memset(edm::view<edm::schema<VARTYPES...>> data,
370 int value) const {
371
372 // For buffers, we can do this in one go.
373 if (data.payload().ptr() != nullptr) {
374 memset(data.payload(), value);
375 } else {
376 // Do the operation using the helper function, recursively.
377 memset_impl<0>(data, value);
378 }
379
380 // Return a new event.
381 return create_event();
382}
383
384template <typename... VARTYPES>
385edm::buffer<edm::details::remove_cv_t<edm::schema<VARTYPES...>>> copy::to(
386 const edm::view<edm::schema<VARTYPES...>>& data, memory_resource& resource,
387 [[maybe_unused]] memory_resource* host_access_resource,
388 type::copy_type cptype) const {
389
390 // Create the result buffer object.
392 ((data.size().ptr() != nullptr) ? data::buffer_type::resizable
394 edm::buffer<edm::details::remove_cv_t<edm::schema<VARTYPES...>>> result;
396 edm::schema<VARTYPES...>>::value) {
397 // Set up a buffer that has (a) jagged vector(s).
399 btype};
400 // The idea here is that since we are going to copy data into the newly
401 // created buffer, the only thing that setup(...) needs to do is to set
402 // up "the layout" in non-host-accessible memory. (Zeroing the size
403 // variable(s) is not needed.) When doing a device-to-host copy, copying
404 // "the layout" is not actually needed. Worse yet, the copy object is
405 // not capable of memsetting the size variable(s). Since the device
406 // specific copy object cannot memset host memory.
407 //
408 // Long story short, the entire setup is skipped in the absence of a
409 // host-accessible memory resource.
410 if (host_access_resource != nullptr) {
411 setup(result)->wait();
412 }
413 } else {
414 // Set up a buffer absent of jagged vectors.
415 result = {data.capacity(), resource, btype};
416 // In the absence of jagged vectors there is no "layout" to copy. And
417 // since the size of the buffer will be set during the copy correctly,
418 // there is nothing else that setup(...) would need to do. So no need
419 // to call it here.
420 }
421
422 // Perform the copy.
423 operator()(data, result, cptype)->wait();
424
425 // Return the buffer.
426 return result;
427}
428
429template <typename... VARTYPES>
431 const edm::view<edm::details::add_const_t<edm::schema<VARTYPES...>>>&
432 from_view,
433 edm::view<edm::schema<VARTYPES...>> to_view, type::copy_type cptype) const {
434
435 // First, handle the simple case, when both views have a contiguous memory
436 // layout.
437 if ((from_view.payload().ptr() != nullptr) &&
438 (to_view.payload().ptr() != nullptr) &&
439 (from_view.payload().capacity() == to_view.payload().capacity())) {
440
441 // If the "common size" is zero, we're done.
442 if (from_view.payload().capacity() == 0) {
444 }
445
446 // Let the user know what's happening.
447 VECMEM_DEBUG_MSG(2, "Performing simple SoA copy of %u bytes",
448 from_view.payload().size());
449
450 // Copy the payload with a single copy operation.
451 copy_view_impl(from_view.payload(), to_view.payload(), cptype);
452
453 // If the target view is resizable, set its size.
454 if (to_view.size().ptr() != nullptr) {
455 // If the source is also resizable, the situation should be simple.
456 if (from_view.size().ptr() != nullptr) {
457 // Check that the sizes are the same.
458 if (from_view.size().capacity() != to_view.size().capacity()) {
459 std::ostringstream msg;
460 msg << "from_view.size().capacity() ("
461 << from_view.size().capacity()
462 << ") != to_view.size().capacity() ("
463 << to_view.size().capacity() << ")";
464 throw std::length_error(msg.str());
465 }
466 // Perform a dumb copy.
467 copy_view_impl(from_view.size(), to_view.size(), cptype);
468 } else {
469 // If not, then copy the size(s) recursively.
470 copy_sizes_impl<0>(from_view, to_view, cptype);
471 }
472 }
473
474 // Create a synchronization event.
475 return create_event();
476 }
477
478 // If not, then do an un-optimized copy, variable-by-variable.
479 copy_payload_impl<0>(from_view, to_view, cptype);
480
481 // Return a new event.
482 return create_event();
483}
484
485template <typename... VARTYPES, template <typename> class INTERFACE>
487 const edm::view<edm::details::add_const_t<edm::schema<VARTYPES...>>>&
488 from_view,
489 edm::host<edm::schema<VARTYPES...>, INTERFACE>& to_vec,
490 type::copy_type cptype) const {
491
492 // Resize the output object to the correct size(s).
493 resize_impl<0>(from_view, to_vec, cptype);
494
495 // Perform the memory copy.
496 return operator()(from_view, vecmem::get_data(to_vec), cptype);
497}
498
499template <typename... VARTYPES>
501 const edm::view<edm::schema<VARTYPES...>>& data) const {
502
503 // Start by taking the capacity of the container.
504 typename edm::view<edm::schema<VARTYPES...>>::size_type size =
505 data.capacity();
506
507 // If there are jagged vectors in the container and/or the container is not
508 // resizable, we're done. It is done in two separate if statements to avoid
509 // MSVC trying to be too smart, and giving a warning...
510 if constexpr (std::disjunction_v<
512 return size;
513 } else {
514 // All the rest is put into an else block, to avoid MSVC trying to be
515 // too smart, and giving a warning about unreachable code...
516 if (data.size().ptr() == nullptr) {
517 return size;
518 }
519
520 // A small security check.
521 assert(data.size().size() ==
522 sizeof(typename edm::view<edm::schema<VARTYPES...>>::size_type));
523
524 // Get the exact size of the container.
525 do_copy(sizeof(typename edm::view<edm::schema<VARTYPES...>>::size_type),
526 data.size().ptr(), &size, type::unknown);
527 // We have to wait for this to finish, since the "size" variable is
528 // not going to be available outside of this function. And
529 // asynchronous SYCL memory copies can happen from variables on the
530 // stack as well...
531 create_event()->wait();
532
533 // Return what we got.
534 return size;
535 }
536}
537
538template <typename... VARTYPES>
539std::vector<data::vector_view<int>::size_type> copy::get_sizes(
540 const edm::view<edm::schema<VARTYPES...>>& data) const {
541
542 // Make sure that there's a jagged vector in here.
543 static_assert(
545 "Function can only be used on containers with jagged vectors");
546
547 // Get the sizes using the helper function.
548 return get_sizes_impl<0>(data);
549}
550
551template <typename TYPE>
552bool copy::copy_view_impl(
553 const data::vector_view<std::add_const_t<TYPE>>& from_view,
554 data::vector_view<TYPE> to_view, type::copy_type cptype) const {
555
556 // Get the size of the source view.
557 const typename data::vector_view<std::add_const_t<TYPE>>::size_type size =
559 // If it's zero, don't do an actual copy.
560 if (size == 0u) {
561 return false;
562 }
563
564 // Make sure that the copy can happen.
565 if (to_view.capacity() < size) {
566 std::ostringstream msg;
567 msg << "Target capacity (" << to_view.capacity() << ") < source size ("
568 << size << ")";
569 throw std::length_error(msg.str());
570 }
571
572 // Make sure that if the target view is resizable, that it would be set up
573 // for the correct size.
574 if (to_view.size_ptr() != nullptr) {
575 // Select what type of copy this should be. Keeping in mind that we copy
576 // from a variable on the host stack. So the question is just whether
577 // the target is the host, or a device.
578 type::copy_type size_cptype = type::unknown;
579 switch (cptype) {
582 size_cptype = type::host_to_device;
583 break;
586 size_cptype = type::host_to_host;
587 break;
588 default:
589 break;
590 }
591 // Perform the copy.
592 do_copy(sizeof(typename data::vector_view<TYPE>::size_type), &size,
593 to_view.size_ptr(), size_cptype);
594 // We have to wait for this to finish, since the "size" variable is not
595 // going to be available outside of this function. And asynchronous SYCL
596 // memory copies can happen from variables on the stack as well...
597 create_event()->wait();
598 }
599
600 // Copy the payload.
601 do_copy(size * sizeof(TYPE), from_view.ptr(), to_view.ptr(), cptype);
602 return true;
603}
604
605template <typename TYPE>
606bool copy::copy_view_impl(
607 const data::jagged_vector_view<std::add_const_t<TYPE>>& from_view,
608 data::jagged_vector_view<TYPE> to_view, type::copy_type cptype) const {
609
610 // Sanity checks.
611 if (from_view.size() > to_view.size()) {
612 std::ostringstream msg;
613 msg << "from_view.size() (" << from_view.size()
614 << ") > to_view.size() (" << to_view.size() << ")";
615 throw std::length_error(msg.str());
616 }
617
618 // Get the "outer size" of the container.
619 const typename data::jagged_vector_view<std::add_const_t<TYPE>>::size_type
620 size = from_view.size();
621 if (size == 0u) {
622 return false;
623 }
624
625 // Calculate the contiguous-ness of the memory allocations.
626 const bool from_is_contiguous = is_contiguous(from_view.host_ptr(), size);
627 const bool to_is_contiguous = is_contiguous(to_view.host_ptr(), size);
628 VECMEM_DEBUG_MSG(3, "from_is_contiguous = %d, to_is_contiguous = %d",
629 from_is_contiguous, to_is_contiguous);
630
631 // Get the sizes of the source jagged vector.
632 const auto sizes = get_sizes(from_view);
633
634 // Before even attempting the copy, make sure that the target view either
635 // has the correct sizes, or can be resized correctly. Captire the event
636 // marking the end of the operation. We need to wait for the copy to finish
637 // before returning from this function.
638 auto set_sizes_event = set_sizes(sizes, to_view);
639
640 // Check whether the source and target capacities match up. We can only
641 // perform the "optimised copy" if they do.
642 std::vector<typename data::vector_view<std::add_const_t<TYPE>>::size_type>
643 capacities(size);
644 bool capacities_match = true;
645 for (std::size_t i = 0; i < size; ++i) {
646 if (from_view.host_ptr()[i].capacity() !=
647 to_view.host_ptr()[i].capacity()) {
648 capacities_match = false;
649 break;
650 }
651 capacities[i] = from_view.host_ptr()[i].capacity();
652 }
653
654 // Perform the copy as best as we can.
655 if (from_is_contiguous && to_is_contiguous && capacities_match) {
656 // Perform the copy in one go.
657 copy_views_contiguous_impl(capacities, from_view.host_ptr(),
658 to_view.host_ptr(), cptype);
659 } else {
660 // Do the copy as best as we can. Note that since they are not
661 // contiguous anyway, we use the sizes of the vectors here, not their
662 // capcities.
663 copy_views_impl(sizes, from_view.host_ptr(), to_view.host_ptr(),
664 cptype);
665 }
666
667 // Before returning, make sure that the "size set" operation has finished.
668 // Because the "sizes" variable is just about to go out of scope.
669 set_sizes_event->wait();
670 return true;
671}
672
673template <typename TYPE>
674void copy::copy_views_impl(
675 const std::vector<typename data::vector_view<TYPE>::size_type>& sizes,
676 const data::vector_view<std::add_const_t<TYPE>>* from_view,
677 data::vector_view<TYPE>* to_view, type::copy_type cptype) const {
678
679 // Some security checks.
680 assert(from_view != nullptr);
681 assert(to_view != nullptr);
682
683 // Helper variable(s) used in the copy.
684 const std::size_t size = sizes.size();
685 [[maybe_unused]] std::size_t copy_ops = 0;
686
687 // Perform the copy in multiple steps.
688 for (std::size_t i = 0; i < size; ++i) {
689
690 // Skip empty "inner vectors".
691 if (sizes[i] == 0) {
692 continue;
693 }
694
695 // Some sanity checks.
696 assert(from_view[i].ptr() != nullptr);
697 assert(to_view[i].ptr() != nullptr);
698 assert(sizes[i] <= from_view[i].capacity());
699 assert(sizes[i] <= to_view[i].capacity());
700
701 // Perform the copy.
702 do_copy(sizes[i] * sizeof(TYPE), from_view[i].ptr(), to_view[i].ptr(),
703 cptype);
704 ++copy_ops;
705 }
706
707 // Let the user know what happened.
708 VECMEM_DEBUG_MSG(2,
709 "Copied the payload of a jagged vector of type "
710 "\"%s\" with %lu copy operation(s)",
711 typeid(TYPE).name(), copy_ops);
712}
713
714template <typename TYPE>
715void copy::copy_views_contiguous_impl(
716 const std::vector<typename data::vector_view<TYPE>::size_type>& sizes,
717 const data::vector_view<std::add_const_t<TYPE>>* from_view,
718 data::vector_view<TYPE>* to_view, type::copy_type cptype) const {
719
720 // Some security checks.
721 assert(from_view != nullptr);
722 assert(to_view != nullptr);
723 assert(is_contiguous(from_view, sizes.size()));
724 assert(is_contiguous(to_view, sizes.size()));
725
726 // Helper variable(s) used in the copy.
727 const std::size_t size = sizes.size();
728 const std::size_t total_size =
729 std::accumulate(sizes.begin(), sizes.end(),
730 static_cast<std::size_t>(0)) *
731 sizeof(TYPE);
732
733 // Find the first non-empty element.
734 for (std::size_t i = 0; i < size; ++i) {
735
736 // Jump over empty elements.
737 if (sizes[i] == 0) {
738 continue;
739 }
740
741 // Some sanity checks.
742 assert(from_view[i].ptr() != nullptr);
743 assert(to_view[i].ptr() != nullptr);
744
745 // Perform the copy.
746 do_copy(total_size, from_view[i].ptr(), to_view[i].ptr(), cptype);
747 break;
748 }
749
750 // Let the user know what happened.
751 VECMEM_DEBUG_MSG(2,
752 "Copied the payload of a jagged vector of type "
753 "\"%s\" with 1 copy operation(s)",
754 typeid(TYPE).name());
755}
756
757template <typename TYPE>
758std::vector<typename data::vector_view<TYPE>::size_type> copy::get_sizes_impl(
759 const data::vector_view<TYPE>* data, std::size_t size) const {
760
761 // Create the result vector.
762 std::vector<typename data::vector_view<TYPE>::size_type> result(size, 0);
763
764 // Try to get the "resizable sizes" first.
765 for (std::size_t i = 0; i < size; ++i) {
766 // Find the first "inner vector" that has a non-zero capacity, and is
767 // resizable.
768 if ((data[i].capacity() != 0) && (data[i].size_ptr() != nullptr)) {
769 // Copy the sizes of the inner vectors into the result vector.
771 (size - i),
772 data[i].size_ptr(), result.data() + i, type::unknown);
773 // Wait for the copy operation to finish. With some backends
774 // (khm... SYCL... khm...) copies can be asynchronous even into
775 // non-pinned host memory.
776 create_event()->wait();
777 // At this point the result vector should have been set up
778 // correctly.
779 return result;
780 }
781 }
782
783 // If we're still here, then the buffer is not resizable. So let's just
784 // collect the capacity of each of the inner vectors.
785 for (std::size_t i = 0; i < size; ++i) {
786 result[i] = data[i].capacity();
787 }
788 return result;
789}
790
791template <typename TYPE>
792bool copy::is_contiguous(const data::vector_view<TYPE>* data,
793 std::size_t size) {
794
795 // We should never call this function for an empty jagged vector.
796 assert(size > 0);
797
798 // Check whether all memory blocks are contiguous.
799 auto ptr = data[0].ptr();
800 for (std::size_t i = 1; i < size; ++i) {
801 if ((ptr + data[i - 1].capacity()) != data[i].ptr()) {
802 return false;
803 }
804 ptr = data[i].ptr();
805 }
806 return true;
807}
808
809template <std::size_t INDEX, typename... VARTYPES>
810void copy::memset_impl(edm::view<edm::schema<VARTYPES...>> data,
811 int value) const {
812
813 // Scalars do not have their own dedicated @c memset functions.
814 if constexpr (edm::type::details::is_scalar<typename std::tuple_element<
815 INDEX, std::tuple<VARTYPES...>>::type>::value) {
816 do_memset(sizeof(typename std::tuple_element<
817 INDEX, std::tuple<VARTYPES...>>::type::type),
818 data.template get<INDEX>(), value);
819 } else {
820 // But vectors and jagged vectors do.
821 memset(data.template get<INDEX>(), value);
822 }
823 // Recurse into the next variable.
824 if constexpr (sizeof...(VARTYPES) > (INDEX + 1)) {
825 memset_impl<INDEX + 1>(data, value);
826 }
827}
828
829template <std::size_t INDEX, typename... VARTYPES,
830 template <typename> class INTERFACE>
831void copy::resize_impl(
832 const edm::view<edm::details::add_const_t<edm::schema<VARTYPES...>>>&
833 from_view,
834 edm::host<edm::schema<VARTYPES...>, INTERFACE>& to_vec,
835 [[maybe_unused]] type::copy_type cptype) const {
836
837 // The target is a host container, so the copy type can't be anything
838 // targeting a device.
839 assert((cptype == type::device_to_host) || (cptype == type::host_to_host) ||
840 (cptype == type::unknown));
841
842 // First, handle containers with no jagged vectors in them.
843 if constexpr (std::disjunction_v<
844 edm::type::details::is_jagged_vector<VARTYPES>...> ==
845 false) {
846 // The single size of the source container.
847 typename edm::view<
848 edm::details::add_const_t<edm::schema<VARTYPES...>>>::size_type
849 size = from_view.capacity();
850 // If the container is resizable, take its size.
851 if (from_view.size().ptr() != nullptr) {
852 // A small security check.
853 assert(from_view.size().size() ==
854 sizeof(typename edm::view<edm::details::add_const_t<
855 edm::schema<VARTYPES...>>>::size_type));
856 // Get the exact size of the container.
857 do_copy(sizeof(typename edm::view<edm::details::add_const_t<
858 edm::schema<VARTYPES...>>>::size_type),
859 from_view.size().ptr(), &size, cptype);
860 // We have to wait for this to finish, since the "size" variable is
861 // not going to be available outside of this function. And
862 // asynchronous SYCL memory copies can happen from variables on the
863 // stack as well...
864 create_event()->wait();
865 }
866 // Resize the target container.
867 VECMEM_DEBUG_MSG(4, "Resizing a (non-jagged) container to size %u",
868 size);
869 to_vec.resize(size);
870 } else {
871 // Resize vector and jagged vector variables one by one. Note that
872 // evaluation order matters here. Since all jagged vectors are vectors,
873 // but not all vectors are jagged vectors. ;-)
874 if constexpr (edm::type::details::is_jagged_vector<
875 typename std::tuple_element<
876 INDEX, std::tuple<VARTYPES...>>::type>::value) {
877 // Get the sizes of this jagged vector.
878 auto sizes = get_sizes(from_view.template get<INDEX>());
879 // Set the "outer size" of the jagged vector.
880 VECMEM_DEBUG_MSG(
881 4, "Resizing jagged vector variable at index %lu to size %lu",
882 INDEX, sizes.size());
883 details::resize_jagged_vector(to_vec.template get<INDEX>(),
884 sizes.size());
885 // Set the "inner sizes" of the jagged vector.
886 for (std::size_t i = 0; i < sizes.size(); ++i) {
887 to_vec.template get<INDEX>()[i].resize(sizes[i]);
888 }
889 } else if constexpr (edm::type::details::is_vector<
890 typename std::tuple_element<
891 INDEX,
892 std::tuple<VARTYPES...>>::type>::value) {
893 // Get the size of this vector.
894 auto size = get_size(from_view.template get<INDEX>());
895 // Resize the target vector.
896 VECMEM_DEBUG_MSG(4,
897 "Resizing vector variable at index %lu to size %u",
898 INDEX, size);
899 to_vec.template get<INDEX>().resize(size);
900 }
901 // Call this function recursively.
902 if constexpr (sizeof...(VARTYPES) > (INDEX + 1)) {
903 resize_impl<INDEX + 1>(from_view, to_vec, cptype);
904 }
905 }
906}
907
908template <std::size_t INDEX, typename... VARTYPES>
909void copy::copy_sizes_impl(
910 [[maybe_unused]] const edm::view<
911 edm::details::add_const_t<edm::schema<VARTYPES...>>>& from_view,
912 [[maybe_unused]] edm::view<edm::schema<VARTYPES...>> to_view,
913 [[maybe_unused]] type::copy_type cptype) const {
914
915 // This should only be called for a resizable target container, with
916 // a non-resizable source container.
917 assert(to_view.size().ptr() != nullptr);
918 assert(from_view.size().ptr() == nullptr);
919
920 // First, handle containers with no jagged vectors in them.
921 if constexpr (std::disjunction_v<
922 edm::type::details::is_jagged_vector<VARTYPES>...> ==
923 false) {
924 // Size of the source container.
925 typename edm::view<
926 edm::details::add_const_t<edm::schema<VARTYPES...>>>::size_type
927 size = from_view.capacity();
928 // Choose the copy type.
929 type::copy_type size_cptype = type::unknown;
930 switch (cptype) {
933 size_cptype = type::host_to_device;
934 break;
937 size_cptype = type::host_to_host;
938 break;
939 default:
940 break;
941 }
942 // Set the size of the target container.
943 do_copy(sizeof(typename edm::view<edm::details::add_const_t<
944 edm::schema<VARTYPES...>>>::size_type),
945 &size, to_view.size().ptr(), size_cptype);
946 // We have to wait for this to finish, since the "size" variable is not
947 // going to be available outside of this function. And asynchronous SYCL
948 // memory copies can happen from variables on the stack as well...
949 create_event()->wait();
950 } else {
951 // For the jagged vector case we recursively copy the sizes of every
952 // jagged vector variable. The rest of the variables are not resizable
953 // in such a container, so they are ignored here.
954 if constexpr (edm::type::details::is_jagged_vector<
955 typename std::tuple_element<
956 INDEX, std::tuple<VARTYPES...>>::type>::value) {
957 // Copy the sizes for this variable.
958 const auto sizes = get_sizes(from_view.template get<INDEX>());
959 set_sizes(sizes, to_view.template get<INDEX>())->wait();
960 }
961 // Call this function recursively.
962 if constexpr (sizeof...(VARTYPES) > (INDEX + 1)) {
963 copy_sizes_impl<INDEX + 1>(from_view, to_view, cptype);
964 }
965 }
966}
967
968template <std::size_t INDEX, typename... VARTYPES>
969void copy::copy_payload_impl(
970 const edm::view<edm::details::add_const_t<edm::schema<VARTYPES...>>>&
971 from_view,
972 edm::view<edm::schema<VARTYPES...>> to_view, type::copy_type cptype) const {
973
974 // Scalars do not have their own dedicated @c copy functions.
975 if constexpr (edm::type::details::is_scalar<typename std::tuple_element<
976 INDEX, std::tuple<VARTYPES...>>::type>::value) {
977 do_copy(sizeof(typename std::tuple_element<
978 INDEX, std::tuple<VARTYPES...>>::type::type),
979 from_view.template get<INDEX>(), to_view.template get<INDEX>(),
980 cptype);
981 } else {
982 // But vectors and jagged vectors do.
983 copy_view_impl(from_view.template get<INDEX>(),
984 to_view.template get<INDEX>(), cptype);
985 }
986 // Recurse into the next variable.
987 if constexpr (sizeof...(VARTYPES) > (INDEX + 1)) {
988 copy_payload_impl<INDEX + 1>(from_view, to_view, cptype);
989 }
990}
991
992template <std::size_t INDEX, typename... VARTYPES>
993std::vector<data::vector_view<int>::size_type> copy::get_sizes_impl(
994 const edm::view<edm::schema<VARTYPES...>>& view) const {
995
996 // Make sure that there's a jagged vector in here.
997 static_assert(
998 edm::details::has_jagged_vector<edm::schema<VARTYPES...>>::value,
999 "Function can only be used on containers with jagged vectors");
1000
1001 // Check if the variable at INDEX is a jagged vector.
1002 if constexpr (edm::type::details::is_jagged_vector<
1003 typename std::tuple_element<
1004 INDEX, std::tuple<VARTYPES...>>::type>::value) {
1005 // If it is, get its (inner) sizes.
1006 return get_sizes(view.template get<INDEX>());
1007 } else if constexpr (sizeof...(VARTYPES) > (INDEX + 1)) {
1008 // If it's not, and this was not the last variable, recurse into the
1009 // next variable.
1010 return get_sizes_impl<INDEX + 1>(view);
1011 }
1012
1013 // We should never get here as long as the code was written correctly.
1014#if defined(__GNUC__)
1015 __builtin_unreachable();
1016#elif defined(_MSC_VER)
1017 __assume(0);
1018#endif
1019}
1020
1021} // 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