9.1
general documentation
cs_array.h
Go to the documentation of this file.
1#ifndef __CS_ARRAY_H__
2#define __CS_ARRAY_H__
3
4/*============================================================================
5 * Array handling utilities.
6 *============================================================================*/
7
8/*
9 This file is part of code_saturne, a general-purpose CFD tool.
10
11 Copyright (C) 1998-2025 EDF S.A.
12
13 This program is free software; you can redistribute it and/or modify it under
14 the terms of the GNU General Public License as published by the Free Software
15 Foundation; either version 2 of the License, or (at your option) any later
16 version.
17
18 This program is distributed in the hope that it will be useful, but WITHOUT
19 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
20 FOR A PARTICULAR PURPOSE. See the GNU General Public License for more
21 details.
22
23 You should have received a copy of the GNU General Public License along with
24 this program; if not, write to the Free Software Foundation, Inc., 51 Franklin
25 Street, Fifth Floor, Boston, MA 02110-1301, USA.
26*/
27
28/*----------------------------------------------------------------------------*/
29
30/*----------------------------------------------------------------------------
31 * Standard C and C++ library headers
32 *----------------------------------------------------------------------------*/
33
34#include <string.h>
35
36/*----------------------------------------------------------------------------
37 * Local headers
38 *----------------------------------------------------------------------------*/
39
40#include "base/cs_defs.h"
41#include "base/cs_base_accel.h"
42#if defined (__CUDACC__)
43#include "base/cs_array_cuda.h"
44#include "base/cs_base_cuda.h"
45#endif
46
47#include "base/cs_dispatch.h"
48#include "base/cs_mdspan.h"
49#include "base/cs_parall.h"
50
51/*----------------------------------------------------------------------------*/
52
53/*=============================================================================
54 * Macro definitions
55 *============================================================================*/
56
57/* Define the way to apply a subset. In case of a copy, the subset is related
58 to the reference (input array) and/or the destination array (output
59 array) */
60
61#define CS_ARRAY_SUBSET_NULL -1
62#define CS_ARRAY_SUBSET_IN 0
63#define CS_ARRAY_SUBSET_OUT 1
64#define CS_ARRAY_SUBSET_INOUT 2
65
66/*============================================================================
67 * Type definitions
68 *============================================================================*/
69
70/*============================================================================
71 * Global variables
72 *============================================================================*/
73
74/*=============================================================================
75 * Public inline function prototypes
76 *============================================================================*/
77
78/*=============================================================================
79 * Public function prototypes
80 *============================================================================*/
81
82#if defined(__cplusplus)
83
84/*----------------------------------------------------------------------------*/
99/*----------------------------------------------------------------------------*/
100
101template <typename T, size_t stride, typename... Arrays>
102void
104 const T *ref_val,
105 Arrays&&... arrays)
106{
107 /* Expand the parameter pack */
108 T* array_ptrs[] = {arrays ... };
109
110#if defined (__CUDACC__)
111 bool is_available_on_device = cs_check_device_ptr(ref_val);
112 for (T* array : array_ptrs)
113 is_available_on_device = is_available_on_device
114 && (cs_check_device_ptr(array)
116
117 if (is_available_on_device) {
118 cudaStream_t stream_ = cs_cuda_get_stream(0);
119 cs_arrays_set_value<T, stride>(stream_,
120 false,
121 n_elts,
122 ref_val,
123 arrays...);
124 return;
125 }
126#endif
127
128 auto set_value = [=](cs_lnum_t i_elt) {
129 for (T* array : array_ptrs)
130 memcpy(array + i_elt*stride, ref_val, stride*sizeof(T));
131 };
132
133 /* Loop through each index and assign values */
134# ifdef _OPENMP
135 int n_threads = cs_parall_n_threads(n_elts, CS_THR_MIN);
136# pragma omp parallel for num_threads(n_threads)
137# endif
138 for (cs_lnum_t i = 0; i < n_elts; i++)
139 set_value(i);
140}
141
142/*----------------------------------------------------------------------------*/
157/*----------------------------------------------------------------------------*/
158
159template <typename T, size_t stride, typename... Arrays>
160void
162 const T ref_val,
163 Arrays&&... arrays)
164{
165 /* Expand the parameter pack */
166 T* array_ptrs[] = {arrays ... };
167
168#if defined (__CUDACC__)
169 bool is_available_on_device = true;
170 for (T* array : array_ptrs)
171 is_available_on_device = is_available_on_device
172 && (cs_check_device_ptr(array)
174
175 if (is_available_on_device) {
176 cudaStream_t stream_ = cs_cuda_get_stream(0);
177 cs_arrays_set_value<T, stride>(stream_,
178 false,
179 n_elts,
180 ref_val,
181 arrays...);
182 return;
183 }
184#endif
185
186 auto set_value = [=](cs_lnum_t i_elt) {
187 for (T* array : array_ptrs)
188 for (size_t k = 0; k < stride; k++) {
189 array[i_elt*stride + k] = ref_val;
190 }
191 };
192
193 /* Loop through each index and assign values */
194# ifdef _OPENMP
195 int n_threads = cs_parall_n_threads(n_elts, CS_THR_MIN);
196# pragma omp parallel for num_threads(n_threads)
197# endif
198 for (cs_lnum_t i = 0; i < n_elts; i++)
199 set_value(i);
200}
201
202/*----------------------------------------------------------------------------*/
221/*----------------------------------------------------------------------------*/
222
223template <typename T, size_t stride, typename... Arrays>
224void
226 const cs_lnum_t n_elts,
227 const T *ref_val,
228 Arrays&&... arrays)
229{
230 /* Expand the parameter pack */
231 T* array_ptrs[] = {arrays ... };
232
233#if defined (__CUDACC__)
234 if (ctx.use_gpu()) {
235 cudaStream_t stream_ = ctx.cuda_stream();
236 cs_arrays_set_value<T, stride>(stream_,
237 true,
238 n_elts,
239 ref_val,
240 arrays...);
241 return;
242 }
243#endif
244
245 auto set_value = [=](cs_lnum_t i_elt) {
246 for (T* array : array_ptrs)
247 memcpy(array + i_elt*stride, ref_val, stride*sizeof(T));
248 };
249
250 /* Loop through each index and assign values */
251# ifdef _OPENMP
252 int n_threads = ctx.n_cpu_threads();
253 if (n_threads < 0)
254 n_threads = cs_parall_n_threads(n_elts, CS_THR_MIN);
255# pragma omp parallel for num_threads(n_threads)
256# endif
257 for (cs_lnum_t i = 0; i < n_elts; i++)
258 set_value(i);
259}
260
261/*----------------------------------------------------------------------------*/
280/*----------------------------------------------------------------------------*/
281
282template <typename T, size_t stride, typename... Arrays>
283void
285 const cs_lnum_t n_elts,
286 const T ref_val,
287 Arrays&&... arrays)
288{
289 /* Expand the parameter pack */
290 T* array_ptrs[] = {arrays ... };
291
292#if defined (__CUDACC__)
293 if (ctx.use_gpu()) {
294 cudaStream_t stream_ = ctx.cuda_stream();
295 cs_arrays_set_value<T, stride>(stream_,
296 true,
297 n_elts,
298 ref_val,
299 arrays...);
300 return;
301 }
302#endif
303
304 auto set_value = [=](cs_lnum_t i_elt) {
305 for (T* array : array_ptrs)
306 for (size_t k = 0; k < stride; k++) {
307 array[i_elt*stride + k] = ref_val;
308 }
309 };
310
311 /* Loop through each index and assign values */
312# ifdef _OPENMP
313 int n_threads = ctx.n_cpu_threads();
314 if (n_threads < 0)
315 n_threads = cs_parall_n_threads(n_elts, CS_THR_MIN);
316# pragma omp parallel for num_threads(n_threads)
317# endif
318 for (cs_lnum_t i = 0; i < n_elts; i++)
319 set_value(i);
320}
321
322/*----------------------------------------------------------------------------*/
339/*----------------------------------------------------------------------------*/
340
341template <typename T, size_t stride, typename... Arrays>
342void
344 const cs_lnum_t n_elts,
345 Arrays&&... arrays)
346{
347 /* Expand the parameter pack */
348 T* array_ptrs[] = {arrays ... };
349
350#if defined (__CUDACC__)
351 if (ctx.use_gpu()) {
352 cudaStream_t stream_ = ctx.cuda_stream();
353 cs_arrays_set_zero<T, stride>(stream_,
354 true,
355 n_elts,
356 arrays...);
357 return;
358 }
359#endif
360
361 constexpr T c = static_cast<T>(0);
362
363 cs_lnum_t n_vals = n_elts*stride;
364
365# ifdef _OPENMP
366 int n_threads = ctx.n_cpu_threads();
367 if (n_threads < 0)
368 n_threads = cs_parall_n_threads(n_vals, CS_THR_MIN);
369# pragma omp parallel num_threads(n_threads)
370# endif
371 {
372 for (T* array : array_ptrs) {
373# ifdef _OPENMP
374# pragma omp for nowait
375# endif
376 for (cs_lnum_t i = 0; i < n_vals; i++)
377 array[i] = c;
378 }
379 };
380}
381
382/*----------------------------------------------------------------------------*/
398/*----------------------------------------------------------------------------*/
399
400template <typename T, size_t stride, typename... Arrays>
401void
403 const cs_lnum_t *elt_ids,
404 const T *ref_val,
405 Arrays&&... arrays)
406{
407 if (n_elts < 1)
408 return;
409
410 if (elt_ids == NULL)
411 cs_arrays_set_value<T, stride>(n_elts, ref_val, arrays...);
412 else {
413 /* Expand the parameter pack */
414 T* array_ptrs[] = {arrays ... };
415
416 auto set_value = [=](cs_lnum_t i_elt) {
417 for (T* array : array_ptrs)
418 memcpy(array + i_elt*stride, ref_val, stride*sizeof(T));
419 };
420
421 /* Loop through each index on the subset and assign values */
422# ifdef _OPENMP
423 int n_threads = cs_parall_n_threads(n_elts, CS_THR_MIN);
424# pragma omp parallel for num_threads(n_threads)
425# endif
426 for (cs_lnum_t i = 0; i < n_elts; i++)
427 set_value(elt_ids[i]);
428 }
429}
430
431/*----------------------------------------------------------------------------*/
447/*----------------------------------------------------------------------------*/
448
449template <typename T, size_t stride, typename... Arrays>
450void
452 const cs_lnum_t *elt_ids,
453 const T ref_val,
454 Arrays&&... arrays)
455{
456 if (n_elts < 1)
457 return;
458
459 if (elt_ids == NULL)
460 cs_arrays_set_value<T, stride>(n_elts, ref_val, arrays...);
461 else {
462
463 /* Expand the parameter pack */
464 T* array_ptrs[] = {arrays ... };
465
466 auto set_value = [=](cs_lnum_t i_elt) {
467 for (T* array : array_ptrs)
468 for (size_t k = 0; k < stride; k++) {
469 array[i_elt*stride + k] = ref_val;
470 }
471 };
472
473 /* Loop through each index on the subset and assign values */
474# ifdef _OPENMP
475 int n_threads = cs_parall_n_threads(n_elts, CS_THR_MIN);
476# pragma omp parallel for num_threads(n_threads)
477# endif
478 for (cs_lnum_t i = 0; i < n_elts; i++)
479 set_value(elt_ids[i]);
480 }
481}
482
483/*----------------------------------------------------------------------------*/
494/*----------------------------------------------------------------------------*/
495
496template <typename T>
497void
499 const T* src,
500 T* dest)
501{
502#if defined (__CUDACC__)
503 bool is_available_on_device = cs_check_device_ptr(src)
504 && cs_check_device_ptr(dest);
505
506 if (is_available_on_device) {
507 cudaStream_t stream_ = cs_cuda_get_stream(0);
508 cs_array_copy<T>(stream_, size, src, dest);
509 return;
510 }
511#endif
512
513# ifdef _OPENMP
514 int n_threads = cs_parall_n_threads(size, CS_THR_MIN);
515# pragma omp parallel for num_threads(n_threads)
516# endif
517 for (cs_lnum_t ii = 0; ii < size; ii++)
518 dest[ii] = src[ii];
519}
520
521/*----------------------------------------------------------------------------*/
535/*----------------------------------------------------------------------------*/
536
537template <typename T>
538void
540 const T *x,
541 const T *y,
542 T *diff)
543{
544 cs_array_copy(size, x, diff);
545
546# ifdef _OPENMP
547 int n_threads = cs_parall_n_threads(size, CS_THR_MIN);
548# pragma omp parallel for num_threads(n_threads)
549# endif
550 for (cs_lnum_t ii = 0; ii < size; ii++)
551 diff[ii] -= y[ii];
552}
553
554#endif // __cplusplus
555
557
558/*----------------------------------------------------------------------------*/
559/*
560 * \brief Assign true to all elements of an array. Case of an array of booleans
561 *
562 * \param[in] size total number of elements to set
563 * \param[in, out] a array to set
564 */
565/*----------------------------------------------------------------------------*/
566
567void
569 bool a[]);
570
571/*----------------------------------------------------------------------------*/
572/*
573 * \brief Assign false to all elements of an array. Case of an array of booleans
574 *
575 * \param[in] size total number of elements to set
576 * \param[in, out] a array to set
577 */
578/*----------------------------------------------------------------------------*/
579
580void
582 bool a[]);
583
584/*----------------------------------------------------------------------------*/
585/*
586 * \brief Assign zero to all elements of an array. Case of a cs_flag_t array.
587 *
588 * \param[in] size total number of elements to set to zero
589 * \param[in, out] a array of flags to set
590 */
591/*----------------------------------------------------------------------------*/
592
593void
595 cs_flag_t a[]);
596
597/*----------------------------------------------------------------------------*/
598/*
599 * \brief Assign zero to all elements of an array. Case of a cs_lnum_t array.
600 *
601 * \param[in] size total number of elements to set to zero
602 * \param[in, out] a array to set
603 */
604/*----------------------------------------------------------------------------*/
605
606void
608 cs_lnum_t a[]);
609
610/*----------------------------------------------------------------------------*/
611/*
612 * \brief Copy values from an array of cs_lnum_t type to another of the same
613 * dimensions.
614 *
615 * \param[in] size number of elements * dimension
616 * \param[in] src source array values (size: n_elts*dim)
617 * \param[out] dest destination array values (size: n_elts*dim)
618 */
619/*----------------------------------------------------------------------------*/
620
621void
623 const cs_lnum_t src[],
624 cs_lnum_t dest[]);
625
626/*----------------------------------------------------------------------------*/
627/*
628 * \brief Assign the value "num" to all elements of an array. Case of a
629 * cs_lnum_t array.
630 *
631 * \param[in] size total number of elements to set
632 * \param[in] num value to set
633 * \param[in, out] a array to set
634 */
635/*----------------------------------------------------------------------------*/
636
637void
639 cs_lnum_t num,
640 cs_lnum_t a[]);
641
642/*----------------------------------------------------------------------------*/
643/*
644 * \brief Assign the value "num" to an array on a selected subset of elements.
645 * if elt_ids is null, then one recovers the function
646 * \ref cs_array_lnum_set_value
647 *
648 * \param[in] n_elts number of elements
649 * \param[in] elt_ids list of ids in the subset or null (size: n_elts)
650 * \param[in] num value to set
651 * \param[in, out] a array to set
652 */
653/*----------------------------------------------------------------------------*/
654
655void
657 const cs_lnum_t elt_ids[],
658 cs_lnum_t num,
659 cs_lnum_t a[]);
660/*----------------------------------------------------------------------------*/
661/*
662 * \brief Assign zero to all elements of an array. Case of a int array.
663 *
664 * \param[in] size total number of elements to set to zero
665 * \param[in, out] a array to set
666 */
667/*----------------------------------------------------------------------------*/
668
669void
671 int a[]);
672
673/*----------------------------------------------------------------------------*/
674/*
675 * \brief Assign the value "num" to all elements of an array. Case of a
676 * int array.
677 *
678 * \param[in] size total number of elements to set
679 * \param[in] num value to set
680 * \param[in, out] a array to set
681 */
682/*----------------------------------------------------------------------------*/
683
684void
686 int num,
687 int a[]);
688
689/*----------------------------------------------------------------------------*/
690/*
691 * \brief Assign the value "num" to an array on a selected subset of elements.
692 * if elt_ids is null, then one recovers the function
693 * \ref cs_array_int_set_value
694 *
695 * \param[in] n_elts number of elements
696 * \param[in] elt_ids list of ids in the subset or null (size: n_elts)
697 * \param[in] num value to set
698 * \param[in, out] a array to set
699 */
700/*----------------------------------------------------------------------------*/
701
702void
704 const cs_lnum_t elt_ids[],
705 int num,
706 int a[]);
707
708/*----------------------------------------------------------------------------*/
709/*
710 * \brief Copy an array ("ref") into another array ("dest") on possibly only a
711 * part of the array(s). Array with stride > 1 are assumed to be
712 * interlaced. The subset of elements on which working is defined by
713 * "elt_ids". The way to apply the subset is defined with the parameter
714 * "mode" as follows:
715 * - Only the "ref" array if mode = 0 (CS_ARRAY_SUBSET_IN)
716 * - Only the "dest" array if mode = 1 (CS_ARRAY_SUBSET_OUT)
717 * - Both "ref" and "dest" arrays if mode = 2 (CS_ARRAY_SUBSET_INOUT)
718 *
719 * It elt_ids is null or mode < 0 (CS_ARRAY_SUBSET_NULL), then the
720 * behavior is as \ref cs_array_real_copy
721 *
722 * One assumes that all arrays are allocated with a correct size.
723 *
724 * \param[in] n_elts number of elements in the array
725 * \param[in] stride number of values for each element
726 * \param[in] mode apply the subset ids to which array(s)
727 * \param[in] elt_ids list of ids in the subset or null (size: n_elts)
728 * \param[in] ref reference values to copy
729 * \param[in, out] dest array storing values after applying the indirection
730 */
731/*----------------------------------------------------------------------------*/
732
733void
735 int stride,
736 const cs_lnum_t elt_ids[],
737 int mode,
738 const cs_real_t ref[],
739 cs_real_t dest[]);
740
741/*----------------------------------------------------------------------------*/
742/*
743 * \brief Copy real values from an array to another of the same dimensions.
744 *
745 * \param[in] size number of elements * dimension
746 * \param[in] src source array values (size: n_elts*dim)
747 * \param[out] dest destination array values (size: n_elts*dim)
748 */
749/*----------------------------------------------------------------------------*/
750
751void
753 const cs_real_t src[],
754 cs_real_t dest[]);
755
756/*----------------------------------------------------------------------------*/
757/*
758 * \brief Multiply each value by a scaling factor s.t. dest *= scaling_factor
759 * If elt_ids is non-null, one applies an indirection.
760 * A stride can also be applied. One assumes an interlaced array.
761 *
762 * \param[in] n_elts number of elements
763 * \param[in] stride number of values for each element
764 * \param[in] elt_ids list of ids in the subset or null (size: n_elts)
765 * \param[in] scaling_factor value of the scaling factor
766 * \param[out] dest destination array values
767 */
768/*----------------------------------------------------------------------------*/
769
770void
772 int stride,
773 const cs_lnum_t *elt_ids,
774 cs_real_t scaling_factor,
775 cs_real_t dest[]);
776
777/*----------------------------------------------------------------------------*/
778/*
779 * \brief Add in place an array s.t. r += l_add
780 *
781 * \param[in] n_elts number of elements
782 * \param[in] l_add array to add
783 * \param[out] r destination array values
784 */
785/*----------------------------------------------------------------------------*/
786
787void
789 const cs_real_t l_add[],
790 cs_real_t r[]);
791
792/*----------------------------------------------------------------------------*/
793/*
794 * \brief Assign a constant value of dim "stride" to an interlaced array
795 * sharing the same stride
796 *
797 * \param[in] n_elts number of elements
798 * \param[in] stride number of values for each element
799 * \param[in] ref_val list of values to assign (size: stride)
800 * \param[in, out] a array to set (size: n_elts*stride)
801 */
802/*----------------------------------------------------------------------------*/
803
804void
806 int stride,
807 const cs_real_t ref_val[],
808 cs_real_t a[]);
809
810/*----------------------------------------------------------------------------*/
811/*
812 * \brief Assign a weighted constant value of dim "stride" to an interlaced
813 * array sharing the same stride. Apply a weight for each element. This
814 * weight is constant for each component of an element.
815 *
816 * \param[in] n_elts number of elements
817 * \param[in] stride number of values for each element
818 * \param[in] ref_val list of values to assign (size: stride)
819 * \param[in] weight values of the weight to apply (size: n_elts)
820 * \param[in, out] a array to set (size: n_elts*stride)
821 */
822/*----------------------------------------------------------------------------*/
823
824void
826 int stride,
827 const cs_real_t ref_val[],
828 const cs_real_t weight[],
829 cs_real_t a[]);
830
831/*----------------------------------------------------------------------------*/
832/*
833 * \brief Assign a constant value of dim "stride" to an interlaced array
834 * sharing the same stride. Only a subset of elements are considered.
835 * If elt_ids is null, then one recovers the function
836 * \ref cs_array_real_set_value
837 *
838 * \param[in] n_elts number of elements
839 * \param[in] stride number of values for each element
840 * \param[in] elt_ids list of ids in the subset or null (size: n_elts)
841 * \param[in] ref_val list of values to assign (size: stride)
842 * \param[in, out] a array to set (size >= n_elts * stride)
843 */
844/*----------------------------------------------------------------------------*/
845
846void
848 int stride,
849 const cs_lnum_t elt_ids[],
850 const cs_real_t ref_val[],
851 cs_real_t a[]);
852
853/*----------------------------------------------------------------------------*/
854/*
855 * \brief Assign a weighted constant value of dim "stride" to an interlaced
856 * array sharing the same stride. Only a subset of elements are
857 * considered. If elt_ids is null, then one recovers the function \ref
858 * cs_array_real_set_wvalue Apply a weight for each element. This
859 * weight is constant for each component of an element.
860 *
861 * \param[in] n_elts number of elements
862 * \param[in] stride number of values for each element
863 * \param[in] elt_ids list of ids in the subset or null (size: n_elts)
864 * \param[in] ref_val list of values to assign (size: stride)
865 * \param[in] weight values of the weight to apply (size >= n_elts)
866 * \param[in, out] a array to set (size >= n_elts*stride)
867 */
868/*----------------------------------------------------------------------------*/
869
870void
872 int stride,
873 const cs_lnum_t elt_ids[],
874 const cs_real_t ref_val[],
875 const cs_real_t weight[],
876 cs_real_t a[]);
877
878/*----------------------------------------------------------------------------*/
879/*
880 * \brief Assign a constant scalar value to an array.
881 *
882 * \param[in] n_elts number of elements
883 * \param[in] ref_val value to assign
884 * \param[in, out] a array to set
885 */
886/*----------------------------------------------------------------------------*/
887
888void
890 cs_real_t ref_val,
891 cs_real_t a[]);
892
893/*----------------------------------------------------------------------------*/
894/*
895 * \brief Assign a weighted constant scalar value to an array.
896 * The weight array has the same size as the array "a".
897 *
898 * \param[in] n_elts number of elements
899 * \param[in] ref_val value to assign
900 * \param[in] weight values of the weight to apply
901 * \param[in, out] a array to set
902 */
903/*----------------------------------------------------------------------------*/
904
905void
907 cs_real_t ref_val,
908 const cs_real_t weight[],
909 cs_real_t a[]);
910
911/*----------------------------------------------------------------------------*/
912/*
913 * \brief Assign a constant scalar value to an array on a selected subset of
914 * elements. If elt_ids is null, then one recovers the function
915 * cs_array_real_set_scalar
916 *
917 * \param[in] n_elts number of elements
918 * \param[in] elt_ids list of ids in the subset or null (size: n_elts)
919 * \param[in] ref_val value to assign
920 * \param[in, out] a array to set
921 */
922/*----------------------------------------------------------------------------*/
923
924void
926 const cs_lnum_t elt_ids[],
927 cs_real_t ref_val,
928 cs_real_t a[]);
929
930/*----------------------------------------------------------------------------*/
931/*
932 * \brief Assign a weighted constant scalar value to an array on a selected
933 * subset of elements. If elt_ids is null, then one recovers the function
934 * cs_array_real_set_wscalar
935 *
936 * \param[in] n_elts number of elements
937 * \param[in] elt_ids list of ids in the subset or null (size: n_elts)
938 * \param[in] ref_val value to assign
939 * \param[in] weight values of weights to apply
940 * \param[in, out] a array to set
941 */
942/*----------------------------------------------------------------------------*/
943
944void
946 const cs_lnum_t elt_ids[],
947 cs_real_t ref_val,
948 const cs_real_t weight[],
949 cs_real_t a[]);
950
951/*----------------------------------------------------------------------------*/
952/*
953 * \brief Assign a constant vector to an array of stride 3 which is interlaced
954 *
955 * \param[in] n_elts number of elements
956 * \param[in] ref_val vector to assign
957 * \param[in, out] a array to set
958 */
959/*----------------------------------------------------------------------------*/
960
961void
963 const cs_real_t ref_val[3],
964 cs_real_t a[]);
965
966/*----------------------------------------------------------------------------*/
967/*
968 * \brief Assign a weighted constant vector value to an interlaced array (of
969 * stride 3). The array of weights has the same size as the array "a".
970 *
971 * \param[in] n_elts number of elements
972 * \param[in] ref_val vector to assign
973 * \param[in] weight values of the weight to apply
974 * \param[in, out] a array to set
975 */
976/*----------------------------------------------------------------------------*/
977
978void
980 const cs_real_t ref_val[3],
981 const cs_real_t weight[],
982 cs_real_t a[]);
983
984/*----------------------------------------------------------------------------*/
985/*
986 * \brief Assign a constant vector to an interlaced array (of stride 3) on a
987 * selected subset of elements. If elt_ids is null, then one recovers the
988 * function cs_array_real_set_vector
989 *
990 * \param[in] n_elts number of elements
991 * \param[in] elt_ids list of ids in the subset or null (size: n_elts)
992 * \param[in] ref_val vector to assign
993 * \param[in, out] a array to set
994 */
995/*----------------------------------------------------------------------------*/
996
997void
999 const cs_lnum_t elt_ids[],
1000 const cs_real_t ref_val[3],
1001 cs_real_t a[]);
1002
1003/*----------------------------------------------------------------------------*/
1004/*
1005 * \brief Assign a weighted constant vector value to an interlaced array (of
1006 * stride 3). The subset selection is given by elt_ids. If null, then
1007 * one recovers the function \ref cs_array_real_set_wvector
1008 * The array of weights has the same size as the array "a".
1009 *
1010 * \param[in] n_elts number of elements
1011 * \param[in] elt_ids list of ids in the subset or null (size: n_elts)
1012 * \param[in] ref_val vector to assign
1013 * \param[in] weight values of the weight to apply
1014 * \param[in, out] a array to set
1015 */
1016/*----------------------------------------------------------------------------*/
1017
1018void
1020 const cs_lnum_t elt_ids[],
1021 const cs_real_t ref_val[3],
1022 const cs_real_t weight[],
1023 cs_real_t a[]);
1024
1025/*----------------------------------------------------------------------------*/
1026/*
1027 * \brief Assign a constant 3x3 tensor to an array (of stride 9) which is
1028 * interlaced
1029 *
1030 * \param[in] n_elts number of elements
1031 * \param[in] ref_tens tensor to assign
1032 * \param[in, out] a array to set
1033 */
1034/*----------------------------------------------------------------------------*/
1035
1036void
1038 const cs_real_t ref_tens[3][3],
1039 cs_real_t a[]);
1040
1041/*----------------------------------------------------------------------------*/
1042/*
1043 * \brief Assign a constant 3x3 tensor to an interlaced array (of stride 9) on
1044 * a subset of elements. If elt_ids is null, then one recovers the
1045 * function cs_array_real_set_tensor
1046 *
1047 * \param[in] n_elts number of elements
1048 * \param[in] elt_ids list of ids defining the subset or nullptr
1049 * \param[in] ref_tens tensor to assign
1050 * \param[in, out] a array to set
1051 */
1052/*----------------------------------------------------------------------------*/
1053
1054void
1056 const cs_lnum_t elt_ids[],
1057 const cs_real_t ref_tens[3][3],
1058 cs_real_t a[]);
1059
1060/*----------------------------------------------------------------------------*/
1061/*
1062 * \brief Assign zero to all elements of an array.
1063 *
1064 * \param[in] size total number of elements to set to zero
1065 * \param[in, out] a array to set
1066 */
1067/*----------------------------------------------------------------------------*/
1068
1069void
1071 cs_real_t a[]);
1072
1073/*----------------------------------------------------------------------------*/
1074/*
1075 * \brief Assign a constant value to an array (deprecated function).
1076 *
1077 * \param[in] n_elts number of associated elements
1078 * \param[in] dim associated dimension
1079 * \param[in] v value to assign
1080 * \param[out] a array values (size: n_elts*dim)
1081 */
1082/*----------------------------------------------------------------------------*/
1083
1084void
1086 cs_lnum_t dim,
1087 cs_real_t v,
1088 cs_real_t a[]);
1089
1090/*----------------------------------------------------------------------------*/
1091
1093
1094/*----------------------------------------------------------------------------*/
1095
1096#if defined(__cplusplus)
1097
1098namespace cs {
1099
1100/*----------------------------------------------------------------------------*/
1108/*----------------------------------------------------------------------------*/
1109
1110template<class T, int N = 1, layout L = layout::right>
1111class array {
1112
1113public:
1114
1115 /*--------------------------------------------------------------------------*/
1119 /*--------------------------------------------------------------------------*/
1120
1123 _extent{0},
1124 _offset{0},
1125 _size(0),
1126 _owner(true),
1127 _data(nullptr)
1128 {
1129#if !defined(__CUDA_ARCH__) && \
1130 !defined(SYCL_LANGUAGE_VERSION) && \
1131 !defined(__HIP_DEVICE_COMPILE__)
1132 _mode = cs_alloc_mode;
1133#else
1134 // use default and avoid compiler warnings
1136#endif
1137 }
1138
1139 /*--------------------------------------------------------------------------*/
1143 /*--------------------------------------------------------------------------*/
1144
1145 CS_F_HOST
1147 (
1148 cs_lnum_t size,
1149 cs_alloc_mode_t alloc_mode = cs_alloc_mode,
1150#if (defined(__GNUC__) || defined(__clang__)) && \
1151 __has_builtin(__builtin_LINE) && \
1152 __has_builtin(__builtin_FILE)
1153 const char *file_name = __builtin_FILE(),
1154 const int line_number = __builtin_LINE()
1155#else
1156 const char *file_name = __FILE__,
1157 const int line_number = __LINE__
1158#endif
1159 )
1160 :
1161 _extent{0},
1162 _offset{0},
1163 _size(size),
1164 _owner(true),
1165 _data(nullptr),
1166 _mode(alloc_mode)
1167 {
1168 /* Only usable for array */
1169 static_assert(N == 1,
1170 "Instantiation using only total size only possible for "
1171 " cs::array<T,1> or cs::array<T>");
1172 set_size_(size);
1173 allocate_(file_name, line_number);
1174 }
1175
1176 /*--------------------------------------------------------------------------*/
1180 /*--------------------------------------------------------------------------*/
1181
1182 CS_F_HOST
1184 (
1185 const cs_lnum_t(&dims)[N],
1186 cs_alloc_mode_t alloc_mode = cs_alloc_mode,
1187#if (defined(__GNUC__) || defined(__clang__)) && \
1188 __has_builtin(__builtin_LINE) && \
1189 __has_builtin(__builtin_FILE)
1190 const char *file_name = __builtin_FILE(),
1191 const int line_number = __builtin_LINE()
1192#else
1193 const char *file_name = __FILE__,
1194 const int line_number = __LINE__
1195#endif
1196 )
1197 :
1198 _extent{0},
1199 _offset{0},
1200 _owner(true),
1201 _data(nullptr),
1202 _mode(alloc_mode)
1203 {
1204 set_size_(dims);
1205 allocate_(file_name, line_number);
1206 }
1207
1208 /*--------------------------------------------------------------------------*/
1212 /*--------------------------------------------------------------------------*/
1213
1214 CS_F_HOST
1216 (
1217 cs_lnum_t size1,
1218 cs_lnum_t size2,
1219 cs_alloc_mode_t alloc_mode = cs_alloc_mode,
1220#if (defined(__GNUC__) || defined(__clang__)) && \
1221 __has_builtin(__builtin_LINE) && \
1222 __has_builtin(__builtin_FILE)
1223 const char *file_name = __builtin_FILE(),
1224 const int line_number = __builtin_LINE()
1225#else
1226 const char *file_name = __FILE__,
1227 const int line_number = __LINE__
1228#endif
1229 )
1230 :
1231 _extent{0},
1232 _offset{0},
1233 _owner(true),
1234 _data(nullptr),
1235 _mode(alloc_mode)
1236 {
1237 /* Only usable for array */
1238 static_assert(N == 2,
1239 "Instantiation using (size1, size2) only possible for "
1240 " cs::array<T,2>");
1241 cs_lnum_t tmp_size[N] = {size1, size2};
1242 set_size_(tmp_size);
1243 allocate_(file_name, line_number);
1244 }
1245
1246 /*--------------------------------------------------------------------------*/
1250 /*--------------------------------------------------------------------------*/
1251
1252 CS_F_HOST
1254 (
1255 cs_lnum_t size1,
1256 cs_lnum_t size2,
1257 cs_lnum_t size3,
1258 cs_alloc_mode_t alloc_mode = cs_alloc_mode,
1259#if (defined(__GNUC__) || defined(__clang__)) && \
1260 __has_builtin(__builtin_LINE) && \
1261 __has_builtin(__builtin_FILE)
1262 const char *file_name = __builtin_FILE(),
1263 const int line_number = __builtin_LINE()
1264#else
1265 const char *file_name = __FILE__,
1266 const int line_number = __LINE__
1267#endif
1268 )
1269 :
1270 _extent{0},
1271 _offset{0},
1272 _owner(true),
1273 _data(nullptr),
1274 _mode(alloc_mode)
1275 {
1276 /* Only usable for array */
1277 static_assert(N == 3,
1278 "Instantiation using (size1, size2, size3) only possible "
1279 " for cs::array<T,3>");
1280 cs_lnum_t tmp_size[N] = {size1, size2, size3};
1281 set_size_(tmp_size);
1282 allocate_(file_name, line_number);
1283 }
1284
1285 /*--------------------------------------------------------------------------*/
1289 /*--------------------------------------------------------------------------*/
1290
1291 CS_F_HOST
1293 (
1294 cs_lnum_t size1,
1295 cs_lnum_t size2,
1296 cs_lnum_t size3,
1297 cs_lnum_t size4,
1298 cs_alloc_mode_t alloc_mode = cs_alloc_mode,
1299#if (defined(__GNUC__) || defined(__clang__)) && \
1300 __has_builtin(__builtin_LINE) && \
1301 __has_builtin(__builtin_FILE)
1302 const char *file_name = __builtin_FILE(),
1303 const int line_number = __builtin_LINE()
1304#else
1305 const char *file_name = __FILE__,
1306 const int line_number = __LINE__
1307#endif
1308 )
1309 :
1310 _extent{0},
1311 _offset{0},
1312 _owner(true),
1313 _data(nullptr),
1314 _mode(alloc_mode)
1315 {
1316 /* Only usable for array */
1317 static_assert(N == 4,
1318 "Instantiation using (size1, size2, size3, size4) only possible "
1319 " for cs::array<T,4>");
1320 cs_lnum_t tmp_size[N] = {size1, size2, size3, size4};
1321 set_size_(tmp_size);
1322 allocate_(file_name, line_number);
1323 }
1324
1325 /*--------------------------------------------------------------------------*/
1329 /*--------------------------------------------------------------------------*/
1330
1331 template<typename... Args>
1334 (
1335 T *data,
1336 Args... indices
1337 )
1338 :
1339 _extent{0},
1340 _offset{0},
1341 _owner(false)
1342 {
1343 check_operator_args_(indices...);
1344 set_size_(indices...);
1345 _data = data;
1346 }
1347
1348 /*--------------------------------------------------------------------------*/
1352 /*--------------------------------------------------------------------------*/
1353
1356 (
1357 T *data,
1358 const cs_lnum_t(&dims)[N]
1359 )
1360 :
1361 _extent{0},
1362 _offset{0},
1363 _owner(false)
1364 {
1365#if !defined(__CUDA_ARCH__) && \
1366 !defined(SYCL_LANGUAGE_VERSION) && \
1367 !defined(__HIP_DEVICE_COMPILE__)
1368 _mode = cs_alloc_mode;
1369#else
1370 // use default and avoid compiler warnings
1372#endif
1373 set_size_(dims);
1374 _data = data;
1375 }
1376
1377 /*--------------------------------------------------------------------------*/
1386 /*--------------------------------------------------------------------------*/
1387
1390 (
1391 const array& other,
1392 bool deep_copy=false,
1393#if (defined(__GNUC__) || defined(__clang__)) && \
1394 __has_builtin(__builtin_LINE) && \
1395 __has_builtin(__builtin_FILE)
1396 const char *file_name = __builtin_FILE(),
1397 const int line_number = __builtin_LINE()
1398#else
1399 const char *file_name = __FILE__,
1400 const int line_number = __LINE__
1401#endif
1402 )
1403 {
1404 set_size_(other._extent);
1405 _mode = other._mode;
1406
1407 /* If shallow copy new instance is not owner. Otherwise same ownership
1408 * as original instance since we copy it.
1409 */
1410#if !defined(__CUDA_ARCH__) && \
1411 !defined(SYCL_LANGUAGE_VERSION) && \
1412 !defined(__HIP_DEVICE_COMPILE__)
1413 _owner = deep_copy;
1414#else
1415 _owner = false;
1416#endif
1417
1418 if (_owner) {
1419#if !defined(__CUDA_ARCH__) && \
1420 !defined(SYCL_LANGUAGE_VERSION) && \
1421 !defined(__HIP_DEVICE_COMPILE__)
1422 // Only HOST can allocate and be owner. We avoid a compiler warning
1423 // since a static test is done above.
1424 allocate_(file_name, line_number);
1425 copy_data(other._data);
1426#else
1427 CS_UNUSED(file_name);
1428 CS_UNUSED(line_number);
1429#endif
1430 }
1431 else {
1432 CS_UNUSED(file_name);
1433 CS_UNUSED(line_number);
1434 _data = other._data;
1435 }
1436 }
1437
1438 /*--------------------------------------------------------------------------*/
1442 /*--------------------------------------------------------------------------*/
1443
1444 CS_F_HOST
1446 (
1447 array&& other
1448 )
1449 : array()
1450 {
1451 swap(*this, other);
1452 }
1453
1454 /*--------------------------------------------------------------------------*/
1458 /*--------------------------------------------------------------------------*/
1459
1462 {
1463 clear();
1464 }
1465
1466 /*--------------------------------------------------------------------------*/
1470 /*--------------------------------------------------------------------------*/
1471
1472 CS_F_HOST
1473 friend void
1475 (
1476 array& first,
1477 array& second
1478 )
1479 {
1480 using std::swap;
1481
1482 /* Swap the different members between the two references. */
1483 swap(first._extent, second._extent);
1484 swap(first._offset, second._offset);
1485 swap(first._size, second._size);
1486 swap(first._owner, second._owner);
1487 swap(first._data, second._data);
1488 swap(first._mode, second._mode);
1489
1490 }
1491
1492 /*--------------------------------------------------------------------------*/
1496 /*--------------------------------------------------------------------------*/
1497
1498 CS_F_HOST
1500 {
1501 swap(*this, other);
1502
1503 return *this;
1504 }
1505
1506 /*--------------------------------------------------------------------------*/
1510 /*--------------------------------------------------------------------------*/
1511
1513 void
1515 {
1516 if (_owner) {
1517 // Avoid compiler warnings. If device, object cannot be owner...
1518#if !defined(__CUDA_ARCH__) && \
1519 !defined(SYCL_LANGUAGE_VERSION) && \
1520 !defined(__HIP_DEVICE_COMPILE__)
1521 CS_FREE(_data);
1522#endif
1523 }
1524 else {
1525 _data = nullptr;
1526 }
1527 _size = 0;
1528 for (int i = 0; i < N; i++) {
1529 _offset[i] = 0;
1530 _extent[i] = 0;
1531 }
1532 }
1533
1534 /*--------------------------------------------------------------------------*/
1540 /*--------------------------------------------------------------------------*/
1541
1545 {
1546 return mdspan<T,N,L>(_data, _extent);
1547 }
1548
1549 /*--------------------------------------------------------------------------*/
1555 /*--------------------------------------------------------------------------*/
1556
1557 template<typename... Args>
1559 mdspan<T, N-sizeof...(Args), L>
1561 (
1562 Args... indices
1563 )
1564 {
1565 check_sub_function_args_(indices...);
1566
1567 return (this->view()).sub_view(indices...);
1568 }
1569
1570 /*--------------------------------------------------------------------------*/
1578 /*--------------------------------------------------------------------------*/
1579
1580 template<int _N_, layout _L_ = L>
1584 (
1585 const cs_lnum_t(&dims)[_N_]
1586 )
1587 {
1588 /* sanity check for total size */
1589 cs_lnum_t s = 1;
1590 for (int i = 0; i < _N_; i++)
1591 s *= dims[i];
1592
1593 if (s != _size) {
1594#if !defined(__CUDA_ARCH__) && \
1595 !defined(SYCL_LANGUAGE_VERSION) && \
1596 !defined(__HIP_DEVICE_COMPILE__)
1597 bft_error(__FILE__, __LINE__, 0,
1598 _("%s: requested span has total size of %d instead of %d "
1599 "for this array.\n"),
1600 __func__, s, _size);
1601#else
1602 return mdspan<T,_N_,_L_>();
1603#endif
1604 }
1605
1606 return mdspan<T,_N_,_L_>(_data, dims);
1607 }
1608
1609 /*--------------------------------------------------------------------------*/
1617 /*--------------------------------------------------------------------------*/
1618
1619 template<layout _L_ = L, typename... Args>
1621 mdspan<T, sizeof...(Args), _L_>
1623 (
1624 Args... indices
1625 )
1626 {
1627 check_non_empty_args_(indices...);
1628
1629 constexpr int n_vals = sizeof...(Args);
1630 cs_lnum_t tmp_args[n_vals] = {indices...};
1631
1632 return get_mdspan<n_vals, _L_>(tmp_args);
1633 }
1634
1635 /*--------------------------------------------------------------------------*/
1639 /*--------------------------------------------------------------------------*/
1640
1641 CS_F_HOST
1643 (
1644 T val,
1645 const cs_lnum_t n_vals = -1
1647 )
1648 {
1649 assert(n_vals <= _size);
1650
1651 const cs_lnum_t loop_size = (n_vals == -1) ? _size : n_vals;
1652
1653 // Explicit pointer, avoid passing internal member of class to the functor
1654 T* data_ptr = _data;
1655
1657
1658 ctx.parallel_for(loop_size, CS_LAMBDA (cs_lnum_t e_id) {
1659 data_ptr[e_id] = val;
1660 });
1661
1662 ctx.wait();
1663 }
1664
1665 /*--------------------------------------------------------------------------*/
1671 /*--------------------------------------------------------------------------*/
1672
1673 CS_F_HOST
1675 (
1676 cs_dispatch_context &ctx,
1677 T val,
1678 const cs_lnum_t n_vals = -1
1680 )
1681 {
1682 assert(n_vals <= _size);
1683
1684 const cs_lnum_t loop_size = (n_vals == -1) ? _size : n_vals;
1685
1686 // Explicit pointer, avoid passing internal member of class to the functor
1687 T* data_ptr = _data;
1688
1689 /* No wait here since context is passed as argument */
1690 ctx.parallel_for(loop_size, CS_LAMBDA (cs_lnum_t e_id) {
1691 data_ptr[e_id] = val;
1692 });
1693 }
1694
1695 /*--------------------------------------------------------------------------*/
1699 /*--------------------------------------------------------------------------*/
1700
1701 CS_F_HOST
1703 (
1704 T val,
1705 const cs_lnum_t n_elts,
1706 const cs_lnum_t elt_ids[]
1707 )
1708 {
1709 assert(n_elts <= _size && n_elts >= 0);
1710
1711 if (n_elts < 1)
1712 return;
1713
1715
1716 if (elt_ids == nullptr)
1717 set_to_val(ctx, val, n_elts);
1718 else {
1719 // Explicit pointer, avoid passing internal member of class to the functor
1720 T* data_ptr = _data;
1721
1722 ctx.parallel_for(n_elts, CS_LAMBDA (cs_lnum_t e_id) {
1723 data_ptr[elt_ids[e_id]] = val;
1724 });
1725 }
1726
1727 ctx.wait();
1728 }
1729
1730 /*--------------------------------------------------------------------------*/
1736 /*--------------------------------------------------------------------------*/
1737
1738 CS_F_HOST
1740 (
1741 cs_dispatch_context &ctx,
1742 T val,
1743 const cs_lnum_t n_elts,
1744 const cs_lnum_t elt_ids[]
1746 )
1747 {
1748 assert(n_elts <= _size && n_elts >= 0);
1749
1750 if (n_elts < 1)
1751 return;
1752
1753 /* No wait here since context is passed as argument */
1754 if (elt_ids == nullptr)
1755 set_to_val(ctx, val, n_elts);
1756 else {
1757 // Explicit pointer, avoid passing internal member of class to the functor
1758 T* data_ptr = _data;
1759
1760 ctx.parallel_for(n_elts, CS_LAMBDA (cs_lnum_t e_id) {
1761 data_ptr[elt_ids[e_id]] = val;
1762 });
1763 }
1764 }
1765
1766 /*--------------------------------------------------------------------------*/
1770 /*--------------------------------------------------------------------------*/
1771
1772 CS_F_HOST
1773 void
1775 {
1777
1778 // Explicit pointer, avoid passing internal member of class to the functor
1779 T* data_ptr = _data;
1780
1781 ctx.parallel_for(_size, CS_LAMBDA (cs_lnum_t e_id) {
1782 data_ptr[e_id] = static_cast<T>(0);
1783 });
1784
1785 ctx.wait();
1786 }
1787
1788 /*--------------------------------------------------------------------------*/
1792 /*--------------------------------------------------------------------------*/
1793
1794 CS_F_HOST
1795 void
1797 (
1799 )
1800 {
1801 // Explicit pointer, avoid passing internal member of class to the functor
1802 T* data_ptr = _data;
1803
1804 ctx.parallel_for(_size, CS_LAMBDA (cs_lnum_t e_id) {
1805 data_ptr[e_id] = static_cast<T>(0);
1806 });
1807 }
1808
1809 /*--------------------------------------------------------------------------*/
1813 /*--------------------------------------------------------------------------*/
1814
1816 void
1818 {
1819 set_size_(0);
1820 _owner = false;
1821 _data = nullptr;
1822 _mode = cs_alloc_mode;
1823 }
1824
1825 /*--------------------------------------------------------------------------*/
1829 /*--------------------------------------------------------------------------*/
1830
1831 CS_F_HOST
1833 (
1834 cs_lnum_t new_size,
1835#if (defined(__GNUC__) || defined(__clang__)) && \
1836 __has_builtin(__builtin_LINE) && \
1837 __has_builtin(__builtin_FILE)
1838 const char *file_name = __builtin_FILE(),
1839 const int line_number = __builtin_LINE()
1840#else
1841 const char *file_name = __FILE__,
1842 const int line_number = __LINE__
1843#endif
1844 )
1845 {
1846 static_assert(N == 1,
1847 "Method reshape(size, ...) only possible for "
1848 " cs::array<T> or cs::array<T,1>");
1849 assert(new_size >= 0);
1850
1851 /* If same dimensions, nothing to do ... */
1852 if (new_size == _size)
1853 return;
1854
1855 if (_owner) {
1856 set_size_(new_size);
1857 reallocate_(file_name, line_number);
1858 }
1859 else
1860 bft_error(__FILE__, __LINE__, 0,
1861 _("%s: array cannot be reshaped if non owner, "
1862 "use mdspan instead.\n"), __func__);
1863
1864 }
1865
1866 /*--------------------------------------------------------------------------*/
1870 /*--------------------------------------------------------------------------*/
1871
1872 CS_F_HOST
1873 void
1875 (
1876 cs_lnum_t new_size,
1877 cs_lnum_t size_to_keep,
1878#if (defined(__GNUC__) || defined(__clang__)) && \
1879 __has_builtin(__builtin_LINE) && \
1880 __has_builtin(__builtin_FILE)
1881 const char *file_name = __builtin_FILE(),
1882 const int line_number = __builtin_LINE()
1883#else
1884 const char *file_name = __FILE__,
1885 const int line_number = __LINE__
1886#endif
1887 )
1888 {
1889 static_assert(N==1,
1890 "Method reshape_and_copy(size, ...) only possible for "
1891 " cs::array<T> or cs::array<T,1>");
1892 assert(new_size >= 0);
1893 assert(size_to_keep <= new_size && size_to_keep <= _size);
1894
1895 /* If same dimensions, nothing to do ... */
1896 if (new_size == _size)
1897 return;
1898
1899 if (!(size_to_keep <= new_size && size_to_keep <= _size))
1900 bft_error(__FILE__, __LINE__, 0,
1901 "%s: Data cannot be saved when new size is smaller "
1902 "than size to keep.\n",
1903 __func__);
1904
1905 if (_owner) {
1906 /* If new size is larger, realloc is sufficient. */
1907 set_size_(new_size);
1908 reallocate_(file_name, line_number);
1909 }
1910 else
1911 bft_error(__FILE__, __LINE__, 0,
1912 _("%s: array cannot be reshaped if non owner, "
1913 "use mdspan instead.\n"), __func__);
1914 }
1915
1916 /*--------------------------------------------------------------------------*/
1920 /*--------------------------------------------------------------------------*/
1921
1922 CS_F_HOST
1924 (
1925 const cs_lnum_t(&dims)[N],
1926#if (defined(__GNUC__) || defined(__clang__)) && \
1927 __has_builtin(__builtin_LINE) && \
1928 __has_builtin(__builtin_FILE)
1929 const char *file_name = __builtin_FILE(),
1930 const int line_number = __builtin_LINE()
1931#else
1932 const char *file_name = __FILE__,
1933 const int line_number = __LINE__
1934#endif
1935 )
1936 {
1937 /* If same dimensions, nothing to do ... */
1938 bool same_dim = true;
1939 for (int i = 0; i < N; i++)
1940 if (dims[i] != _extent[i])
1941 same_dim = false;
1942
1943 if (same_dim)
1944 return;
1945
1946 if (_owner) {
1947 set_size_(dims);
1948 reallocate_(file_name, line_number);
1949 }
1950 else
1951 bft_error(__FILE__, __LINE__, 0,
1952 _("%s: array cannot be reshaped if non owner, "
1953 "use mdspan instead.\n"), __func__);
1954 }
1955
1956 /*--------------------------------------------------------------------------*/
1960 /*--------------------------------------------------------------------------*/
1961
1962 CS_F_HOST
1964 (
1965 const cs_lnum_t(&dims)[N],
1966 cs_lnum_t size_to_keep = -1,
1967#if (defined(__GNUC__) || defined(__clang__)) && \
1968 __has_builtin(__builtin_LINE) && \
1969 __has_builtin(__builtin_FILE)
1970 const char *file_name = __builtin_FILE(),
1971 const int line_number = __builtin_LINE()
1972#else
1973 const char *file_name = __FILE__,
1974 const int line_number = __LINE__
1975#endif
1976 )
1977 {
1978 assert(size_to_keep <= dims[N-1] && size_to_keep <= _extent[N-1]);
1979
1980 /* If same dimensions, nothing to do ... */
1981 bool same_dim = true;
1982 for (int i = 0; i < N; i++)
1983 if (dims[i] != _extent[i])
1984 same_dim = false;
1985
1986 if (same_dim)
1987 return;
1988
1989 if (_owner) {
1990 /* check that size_to_keep is either >0 or == -1 (keep all) */
1991 if (size_to_keep != 0) {
1992 /* Temporary copy */
1993 array<T,N,L> tmp(*this, true);
1994
1995 /* Update instance size */
1996 set_size_(dims);
1997 reallocate_(file_name, line_number);
1998
1999 /* Local copies for the context loop */
2000 cs_lnum_t new_o[N], new_e[N], old_o[N], old_e[N], max_e[N];
2001
2002 for (int i = 0; i < N; i++) {
2003 new_o[i] = _offset[i];
2004 new_e[i] = _extent[i];
2005 old_o[i] = tmp._offset[i];
2006 old_e[i] = tmp._extent[i];
2007
2008 max_e[i] = (new_e[i] > old_e[i]) ? old_e[i] : new_e[i];
2009 }
2010 if (max_e[N-1] > size_to_keep && size_to_keep > 0)
2011 max_e[N-1] = size_to_keep;
2012
2013 T* old_data = tmp._data;
2014 T* new_data = _data;
2015
2016 cs_lnum_t loop_size = tmp._size;
2017
2018 /* Loop using dispatch */
2020
2021 ctx.parallel_for(loop_size, CS_LAMBDA (cs_lnum_t e_id) {
2022 cs_lnum_t idx[N];
2023 cs_lnum_t dummy = e_id;
2024
2025 for (int i = N-1; i >= 1; i--) {
2026 idx[i] = dummy % old_e[i];
2027 dummy = (dummy - idx[i]) / old_e[i];
2028 }
2029 idx[0] = dummy;
2030
2031 bool all_ok = true;
2032 for (int i = 0; i < N; i++) {
2033 if (idx[i] >= max_e[i]) {
2034 all_ok = false;
2035 break;
2036 }
2037 }
2038
2039 if (all_ok) {
2040 cs_lnum_t o_id = 0;
2041 cs_lnum_t n_id = 0;
2042
2043 for (int i = 0; i < N; i++) {
2044 o_id += idx[i] * old_o[i];
2045 n_id += idx[i] * new_o[i];
2046 }
2047
2048 new_data[n_id] = old_data[o_id];
2049 }
2050 });
2051
2052 ctx.wait();
2053 }
2054 else {
2055 reshape(dims, file_name, line_number);
2056 }
2057 }
2058 else
2059 bft_error(__FILE__, __LINE__, 0,
2060 _("%s: array cannot be reshaped if non owner, "
2061 "use mdspan instead.\n"), __func__);
2062 }
2063
2064 /*--------------------------------------------------------------------------*/
2068 /*--------------------------------------------------------------------------*/
2069
2070 CS_F_HOST
2071 void
2073 (
2074 const cs_lnum_t size1,
2075 const cs_lnum_t size2,
2076#if (defined(__GNUC__) || defined(__clang__)) && \
2077 __has_builtin(__builtin_LINE) && \
2078 __has_builtin(__builtin_FILE)
2079 const char *file_name = __builtin_FILE(),
2080 const int line_number = __builtin_LINE()
2081#else
2082 const char *file_name = __FILE__,
2083 const int line_number = __LINE__
2084#endif
2085 )
2086 {
2087 static_assert(N == 2,
2088 "Method reshape(size1, size2, ...) only possible for "
2089 " cs::array<T,2>");
2090 cs_lnum_t tmp_size[N] = {size1, size2};
2091
2092 reshape(tmp_size, file_name, line_number);
2093 }
2094
2095 /*--------------------------------------------------------------------------*/
2099 /*--------------------------------------------------------------------------*/
2100
2101 CS_F_HOST
2102 void
2104 (
2105 const cs_lnum_t size1,
2106 const cs_lnum_t size2,
2107 cs_lnum_t size_to_keep,
2108#if (defined(__GNUC__) || defined(__clang__)) && \
2109 __has_builtin(__builtin_LINE) && \
2110 __has_builtin(__builtin_FILE)
2111 const char *file_name = __builtin_FILE(),
2112 const int line_number = __builtin_LINE()
2113#else
2114 const char *file_name = __FILE__,
2115 const int line_number = __LINE__
2116#endif
2117 )
2118 {
2119 static_assert(N == 2,
2120 "Method reshape_and_copy(size1, size2, ...) only possible for "
2121 " cs::array<T,2>");
2122 cs_lnum_t tmp_size[N] = {size1, size2};
2123
2124 reshape_and_copy(tmp_size, size_to_keep, file_name, line_number);
2125 }
2126
2127 /*--------------------------------------------------------------------------*/
2131 /*--------------------------------------------------------------------------*/
2132
2133 CS_F_HOST
2134 void
2136 (
2137 const cs_lnum_t size1,
2138 const cs_lnum_t size2,
2139 const cs_lnum_t size3,
2140#if (defined(__GNUC__) || defined(__clang__)) && \
2141 __has_builtin(__builtin_LINE) && \
2142 __has_builtin(__builtin_FILE)
2143 const char *file_name = __builtin_FILE(),
2144 const int line_number = __builtin_LINE()
2145#else
2146 const char *file_name = __FILE__,
2147 const int line_number = __LINE__
2148#endif
2149 )
2150 {
2151 static_assert(N == 3,
2152 "Method reshape(size1, size2, size3, ...) only possible for "
2153 " cs::array<T,3>");
2154 cs_lnum_t tmp_size[N] = {size1, size2, size3};
2155
2156 reshape(tmp_size, file_name, line_number);
2157 }
2158
2159 /*--------------------------------------------------------------------------*/
2163 /*--------------------------------------------------------------------------*/
2164
2166 void
2168 (
2169 const cs_lnum_t size1,
2170 const cs_lnum_t size2,
2171 const cs_lnum_t size3,
2172 cs_lnum_t size_to_keep,
2173#if (defined(__GNUC__) || defined(__clang__)) && \
2174 __has_builtin(__builtin_LINE) && \
2175 __has_builtin(__builtin_FILE)
2176 const char *file_name = __builtin_FILE(),
2177 const int line_number = __builtin_LINE()
2178#else
2179 const char *file_name = __FILE__,
2180 const int line_number = __LINE__
2181#endif
2182 )
2183 {
2184 static_assert(N == 3,
2185 "Method reshape_and_copy(size1, size2, size3, ...) only "
2186 "possible for cs::array<T,3>");
2187 cs_lnum_t tmp_size[N] = {size1, size2, size3};
2188
2189 reshape_and_copy(tmp_size, size_to_keep, file_name, line_number);
2190 }
2191
2192 /*--------------------------------------------------------------------------*/
2196 /*--------------------------------------------------------------------------*/
2197
2200 (
2202 )
2203 {
2204 if (_mode != mode) {
2205 /* Since allocation mode is changed, deallocate data */
2206 if (_size != 0)
2207 clear();
2208
2209 _mode = mode;
2210 }
2211 }
2212
2213 /*--------------------------------------------------------------------------*/
2219 /*--------------------------------------------------------------------------*/
2220
2221 template<typename... Args>
2223 T*
2225 (
2226 Args... indices
2227 )
2228 {
2229 check_sub_function_args_(indices...);
2230
2231 return _data + contiguous_data_offset_(indices...);
2232 }
2233
2234 /*--------------------------------------------------------------------------*/
2240 /*--------------------------------------------------------------------------*/
2241
2242 template<typename... Args>
2244 T*
2246 (
2247 Args... indices
2248 ) const
2249 {
2250 check_sub_function_args_(indices...);
2251
2252 return _data + contiguous_data_offset_(indices...);
2253 }
2254
2255 /*--------------------------------------------------------------------------*/
2259 /*--------------------------------------------------------------------------*/
2260
2262 T*
2264 {
2265 return _data;
2266 }
2267
2268 /*--------------------------------------------------------------------------*/
2276 /*--------------------------------------------------------------------------*/
2277
2278 template<typename U>
2280 U*
2282 {
2283 return reinterpret_cast<U*>(_data);
2284 }
2285
2286 /*--------------------------------------------------------------------------*/
2290 /*--------------------------------------------------------------------------*/
2291
2293 const T*
2294 data() const
2295 {
2296 return _data;
2297 }
2298
2299 /*--------------------------------------------------------------------------*/
2305 /*--------------------------------------------------------------------------*/
2306
2308 inline
2309 T& operator[]
2310 (
2311 cs_lnum_t i
2312 )
2313 {
2314 return _data[i];
2315 }
2316
2317 /*--------------------------------------------------------------------------*/
2323 /*--------------------------------------------------------------------------*/
2324
2326 inline
2327 T& operator[]
2328 (
2329 cs_lnum_t i
2330 ) const
2331 {
2332 return _data[i];
2333 }
2334
2335 /*--------------------------------------------------------------------------*/
2341 /*--------------------------------------------------------------------------*/
2342
2343 template<typename Id1>
2345 inline
2346 std::enable_if_t<cs::always_true<Id1>::value && N==1, T&>
2347 operator()
2348 (
2349 Id1 i
2350 )
2351 {
2352 check_operator_args_(i);
2353 return _data[i];
2354 }
2355
2356 /*--------------------------------------------------------------------------*/
2362 /*--------------------------------------------------------------------------*/
2363
2364 template<typename Id1>
2366 inline
2367 std::enable_if_t<cs::always_true<Id1>::value && N==1, T&>
2368 operator()
2369 (
2370 Id1 i
2371 ) const
2372 {
2373 check_operator_args_(i);
2374 return _data[i];
2375 }
2376
2377 /*--------------------------------------------------------------------------*/
2383 /*--------------------------------------------------------------------------*/
2384
2385 template<typename Id1, typename Id2>
2387 inline
2388 std::enable_if_t<cs::always_true<Id1, Id2>::value && N==2, T&>
2389 operator()
2390 (
2391 Id1 i,
2392 Id2 j
2393 )
2394 {
2395 check_operator_args_(i,j);
2396
2397 if (L == layout::right)
2398 return _data[i*_offset[0] + j];
2399 else if (L == layout::left)
2400 return _data[i + j*_offset[1]];
2401 else
2402 return _data[i*_offset[0] + j*_offset[1]];
2403 }
2404
2405 /*--------------------------------------------------------------------------*/
2411 /*--------------------------------------------------------------------------*/
2412
2413 template<typename Id1, typename Id2>
2415 inline
2416 std::enable_if_t<cs::always_true<Id1, Id2>::value && N==2, T&>
2417 operator()
2418 (
2419 Id1 i,
2420 Id2 j
2421 ) const
2422 {
2423 check_operator_args_(i,j);
2424
2425 if (L == layout::right)
2426 return _data[i*_offset[0] + j];
2427 else if (L == layout::left)
2428 return _data[i + j*_offset[1]];
2429 else
2430 return _data[i*_offset[0] + j*_offset[1]];
2431 }
2432
2433 /*--------------------------------------------------------------------------*/
2439 /*--------------------------------------------------------------------------*/
2440
2441 template<typename... Args>
2443 inline
2444 std::enable_if_t<cs::always_true<Args...>::value && (N!=2) && (N!=1), T&>
2445 operator()
2446 (
2447 Args... indices
2448 )
2449 {
2450 check_operator_args_(indices...);
2451
2452 return _data[data_offset_(indices...)];
2453 }
2454
2455 /*--------------------------------------------------------------------------*/
2461 /*--------------------------------------------------------------------------*/
2462
2463 template<typename... Args>
2465 inline
2466 std::enable_if_t<cs::always_true<Args...>::value && (N!=2) && (N!=1), T&>
2467 operator()
2468 (
2469 Args... indices
2470 ) const
2471 {
2472 check_operator_args_(indices...);
2473
2474 return _data[data_offset_(indices...)];
2475 }
2476
2477 /*--------------------------------------------------------------------------*/
2483 /*--------------------------------------------------------------------------*/
2484
2486 inline
2487 cs_lnum_t
2489 (
2490 int i
2491 )
2492 {
2493 return _extent[i];
2494 }
2495
2496 /*--------------------------------------------------------------------------*/
2502 /*--------------------------------------------------------------------------*/
2503
2505 inline
2506 cs_lnum_t
2508 (
2509 int i
2510 ) const
2511 {
2512 return _extent[i];
2513 }
2514
2515 /*--------------------------------------------------------------------------*/
2521 /*--------------------------------------------------------------------------*/
2522
2524 inline
2525 cs_lnum_t
2527 (
2528 int i
2529 )
2530 {
2531 return _offset[i];
2532 }
2533
2534 /*--------------------------------------------------------------------------*/
2540 /*--------------------------------------------------------------------------*/
2541
2543 inline
2544 cs_lnum_t
2546 (
2547 int i
2548 ) const
2549 {
2550 return _offset[i];
2551 }
2552
2553 /*--------------------------------------------------------------------------*/
2559 /*--------------------------------------------------------------------------*/
2560
2562 inline
2563 cs_lnum_t
2565 {
2566 return _size;
2567 }
2568
2569 /*--------------------------------------------------------------------------*/
2575 /*--------------------------------------------------------------------------*/
2576
2578 inline
2579 cs_lnum_t
2580 size() const
2581 {
2582 return _size;
2583 }
2584
2585 /*--------------------------------------------------------------------------*/
2591 /*--------------------------------------------------------------------------*/
2592
2594 inline
2595 bool
2597 {
2598 return _owner;
2599 }
2600
2601 /*--------------------------------------------------------------------------*/
2607 /*--------------------------------------------------------------------------*/
2608
2610 inline
2611 bool
2612 owner() const
2613 {
2614 return _owner;
2615 }
2616
2617 /*--------------------------------------------------------------------------*/
2623 /*--------------------------------------------------------------------------*/
2624
2626 inline
2629 {
2630 return _mode;
2631 }
2632
2633 /*--------------------------------------------------------------------------*/
2639 /*--------------------------------------------------------------------------*/
2640
2642 inline
2644 mode() const
2645 {
2646 return _mode;
2647 }
2648
2649 /*--------------------------------------------------------------------------*/
2654 /*--------------------------------------------------------------------------*/
2655
2656 CS_F_HOST
2657 void
2659 (
2660 T *data,
2661 const cs_lnum_t n_vals = -1
2663 )
2664 {
2665 assert(n_vals <= _size);
2666
2667 const cs_lnum_t loop_size = (n_vals == -1) ? _size : n_vals;
2668
2669 // Explicit pointer, avoid passing internal member of class to the functor
2670 T* data_ptr = _data;
2671
2673
2674 ctx.parallel_for(loop_size, CS_LAMBDA (cs_lnum_t e_id) {
2675 data_ptr[e_id] = data[e_id];
2676 });
2677
2678 ctx.wait();
2679 }
2680
2681 /*--------------------------------------------------------------------------*/
2686 /*--------------------------------------------------------------------------*/
2687
2688 CS_F_HOST
2689 void
2691 (
2692 const T *data,
2693 const cs_lnum_t n_vals = -1
2695 )
2696 {
2697 assert(n_vals <= _size);
2698
2699 const cs_lnum_t loop_size = (n_vals == -1) ? _size : n_vals;
2700
2701 // Explicit pointer, avoid passing internal member of class to the functor
2702 T* data_ptr = _data;
2703
2705
2706 ctx.parallel_for(loop_size, CS_LAMBDA (cs_lnum_t e_id) {
2707 data_ptr[e_id] = data[e_id];
2708 });
2709
2710 ctx.wait();
2711 }
2712
2713 /*--------------------------------------------------------------------------*/
2718 /*--------------------------------------------------------------------------*/
2719
2720 CS_F_HOST
2721 void
2723 (
2724 array& other,
2725 const cs_lnum_t n_vals = -1
2727 )
2728 {
2729 const cs_lnum_t loop_size = (n_vals == -1) ? _size : n_vals;
2730
2731 assert(loop_size <= _size);
2732 assert(loop_size <= other._size);
2733
2734 // Explicit pointer, avoid passing internal member of class to the functor
2735 T* data_ptr = _data;
2736 T* o_data_ptr = other._data;
2737
2739
2740 ctx.parallel_for(loop_size, CS_LAMBDA (cs_lnum_t e_id) {
2741 data_ptr[e_id] = o_data_ptr[e_id];
2742 });
2743
2744 ctx.wait();
2745 }
2746
2747 /*--------------------------------------------------------------------------*/
2752 /*--------------------------------------------------------------------------*/
2753
2754 CS_F_HOST
2755 void
2757 (
2758 mdspan<T, N, L>& span,
2759 const cs_lnum_t n_vals = -1
2761 )
2762 {
2763 const cs_lnum_t loop_size = (n_vals == -1) ? _size : n_vals;
2764
2765 assert(loop_size <= _size);
2766 assert(loop_size <= span._size);
2767
2768 // Explicit pointer, avoid passing internal member of class to the functor
2769 T* data_ptr = _data;
2770 T* s_data_ptr = span._data;
2771
2773
2774 ctx.parallel_for(loop_size, CS_LAMBDA (cs_lnum_t e_id) {
2775 data_ptr[e_id] = s_data_ptr[e_id];
2776 });
2777
2778 ctx.wait();
2779 }
2780
2781 /*--------------------------------------------------------------------------*/
2787 /*--------------------------------------------------------------------------*/
2788
2789 CS_F_HOST
2790 void
2792 (
2793 cs_dispatch_context &ctx,
2794 T *data,
2795 const cs_lnum_t n_vals = -1
2797 )
2798 {
2799 assert(n_vals <= _size);
2800 const cs_lnum_t loop_size = (n_vals == -1) ? _size : n_vals;
2801
2802 // Explicit pointer, avoid passing internal member of class to the functor
2803 T* data_ptr = _data;
2804
2805 ctx.parallel_for(loop_size, CS_LAMBDA (cs_lnum_t e_id) {
2806 data_ptr[e_id] = data[e_id];
2807 });
2808 }
2809
2810 /*--------------------------------------------------------------------------*/
2816 /*--------------------------------------------------------------------------*/
2817
2818 CS_F_HOST
2819 void
2821 (
2822 cs_dispatch_context &ctx,
2823 const T *data,
2824 const cs_lnum_t n_vals = -1
2826 )
2827 {
2828 assert(n_vals <= _size);
2829 const cs_lnum_t loop_size = (n_vals == -1) ? _size : n_vals;
2830
2831 // Explicit pointer, avoid passing internal member of class to the functor
2832 T* data_ptr = _data;
2833
2834 ctx.parallel_for(loop_size, CS_LAMBDA (cs_lnum_t e_id) {
2835 data_ptr[e_id] = data[e_id];
2836 });
2837 }
2838
2839 /*--------------------------------------------------------------------------*/
2845 /*--------------------------------------------------------------------------*/
2846
2847 CS_F_HOST
2848 void
2850 (
2851 cs_dispatch_context &ctx,
2852 array &other,
2853 const cs_lnum_t n_vals = -1
2855 )
2856 {
2857 const cs_lnum_t loop_size = (n_vals == -1) ? _size : n_vals;
2858
2859 assert(loop_size <= _size);
2860 assert(loop_size <= other._size);
2861
2862 // Explicit pointer, avoid passing internal member of class to the functor
2863 T* data_ptr = _data;
2864 T* o_data_ptr = other._data;
2865
2866 ctx.parallel_for(loop_size, CS_LAMBDA (cs_lnum_t e_id) {
2867 data_ptr[e_id] = o_data_ptr[e_id];
2868 });
2869 }
2870
2871 /*--------------------------------------------------------------------------*/
2877 /*--------------------------------------------------------------------------*/
2878
2879 CS_F_HOST
2880 void
2882 (
2883 cs_dispatch_context &ctx,
2884 mdspan<T,N,L> &span,
2885 const cs_lnum_t n_vals = -1
2887 )
2888 {
2889 const cs_lnum_t loop_size = (n_vals == -1) ? _size : n_vals;
2890
2891 assert(loop_size <= _size);
2892 assert(loop_size <= span._size);
2893
2894 // Explicit pointer, avoid passing internal member of class to the functor
2895 T* data_ptr = _data;
2896 T* s_data_ptr = span._data;
2897
2898 ctx.parallel_for(loop_size, CS_LAMBDA (cs_lnum_t e_id) {
2899 data_ptr[e_id] = s_data_ptr[e_id];
2900 });
2901 }
2902
2903private:
2904
2905 /*--------------------------------------------------------------------------*/
2909 /*--------------------------------------------------------------------------*/
2910
2911 template<typename... Args>
2913 static inline
2914 void
2915 check_operator_args_
2916 (
2917 Args...
2918 )
2919 {
2920 static_assert(sizeof...(Args) == N, "Wrong number of arguments");
2921 static_assert(cs::are_integral<Args...>::value, "Non integral input arguments.");
2922 }
2923
2924 /*--------------------------------------------------------------------------*/
2928 /*--------------------------------------------------------------------------*/
2929
2930 template<typename... Args>
2932 static inline
2933 void
2934 check_sub_function_args_
2935 (
2936 Args...
2937 )
2938 {
2939 static_assert(sizeof...(Args) < N, "Too many input arguments.");
2940 static_assert(cs::are_integral<Args...>::value, "Non integral input arguments.");
2941 }
2942
2943 /*--------------------------------------------------------------------------*/
2947 /*--------------------------------------------------------------------------*/
2948
2949 template<typename... Args>
2951 static inline
2952 void
2953 check_non_empty_args_
2954 (
2955 Args...
2956 )
2957 {
2958 static_assert(sizeof...(Args) > 0, "At least one input argument is needed.");
2959 static_assert(cs::are_integral<Args...>::value, "Non integral input arguments.");
2960 }
2961
2962 /*--------------------------------------------------------------------------*/
2966 /*--------------------------------------------------------------------------*/
2967
2968 template<typename... Args>
2970 inline
2971 cs_lnum_t
2972 data_offset_
2973 (
2974 Args... indices
2975 ) const
2976 {
2977 static_assert(sizeof...(Args) <= N && sizeof...(Args) > 0,
2978 "Too many or too few input arguments.");
2979
2980 constexpr int n_idx = sizeof...(Args);
2981
2982 cs_lnum_t _indices[n_idx] = {indices...};
2983
2984 cs_lnum_t retval = 0;
2985 for (int i = 0; i < n_idx; i++)
2986 retval +=_indices[i] * _offset[i];
2987
2988 return retval;
2989 }
2990
2991 /*--------------------------------------------------------------------------*/
2995 /*--------------------------------------------------------------------------*/
2996
2997 template<typename... Args>
2999 inline
3000 cs_lnum_t
3001 contiguous_data_offset_
3002 (
3003 Args... indices
3004 ) const
3005 {
3006 static_assert(sizeof...(Args) <= N && sizeof...(Args) > 0,
3007 "Too many or too few input arguments.");
3008
3009 constexpr int n_idx = sizeof...(Args);
3010
3011 cs_lnum_t _indices[n_idx] = {indices...};
3012
3013 cs_lnum_t retval = 0;
3014 if (L == layout::right) {
3015 for (int i = 0; i < n_idx; i++)
3016 retval +=_indices[i] * _offset[i];
3017 }
3018 else if (L == layout::left) {
3019 for (int i = 0; i < n_idx; i++)
3020 retval +=_indices[i] * _offset[N-1-i];
3021 }
3022
3023 return retval;
3024 }
3025
3026 /*--------------------------------------------------------------------------*/
3030 /*--------------------------------------------------------------------------*/
3031
3033 void
3034 set_size_
3035 (
3036 const cs_lnum_t(&dims)[N]
3037 )
3038 {
3039 _size = (N > 0) ? 1 : 0;
3040 for (int i = 0; i < N; i++) {
3041 _extent[i] = dims[i];
3042 _size *= dims[i];
3043 _offset[i] = 1;
3044 }
3045
3046 /* Compute offset values for getters */
3047
3048 if (L == layout::right) {
3049 /* Version for Layout right */
3050 for (int i = 0; i < N-1; i++) {
3051 for (int j = i + 1; j < N; j++)
3052 _offset[i] *= dims[j];
3053 }
3054 }
3055 else if (L == layout::left) {
3056 for (int i = N-1; i >= 1; i--) {
3057 for (int j = i - 1; j >= 0; j--)
3058 _offset[i] *= dims[j];
3059 }
3060 }
3061 }
3062
3063 /*--------------------------------------------------------------------------*/
3067 /*--------------------------------------------------------------------------*/
3068
3069 template<typename... Args>
3071 void
3072 set_size_
3073 (
3074 Args... dims
3075 )
3076 {
3077 check_operator_args_(dims...);
3078
3079 cs_lnum_t loc_dims[N] = {dims...};
3080
3081 set_size_(loc_dims);
3082 }
3083
3084 /*--------------------------------------------------------------------------*/
3088 /*--------------------------------------------------------------------------*/
3089
3090 CS_F_HOST
3091 void
3092 allocate_
3093 (
3094 const char *file_name,
3095 const int line_number
3096 )
3097 {
3098 const char *_ptr_name = "cs::array._data";
3099 _data = static_cast<T *>(cs_mem_malloc_hd(_mode,
3100 _size,
3101 sizeof(T),
3102 _ptr_name,
3103 file_name,
3104 line_number));
3105 }
3106
3107 /*--------------------------------------------------------------------------*/
3111 /*--------------------------------------------------------------------------*/
3112
3113 CS_F_HOST
3114 void
3115 reallocate_
3116 (
3117 const char *file_name,
3118 const int line_number
3119 )
3120 {
3121 /* If not owner no-op */
3122 if (_owner) {
3123 const char *_ptr_name = "cs::array._data";
3124 _data = static_cast<T *>(cs_mem_realloc_hd(_data,
3125 _mode,
3126 _size,
3127 sizeof(T),
3128 _ptr_name,
3129 file_name,
3130 line_number));
3131 }
3132 };
3133
3134 /*===========================================================================
3135 * Private members
3136 *==========================================================================*/
3137
3138 cs_lnum_t _extent[N];
3139 cs_lnum_t _offset[N];
3140 cs_lnum_t _size;
3141 bool _owner;
3142 T* _data;
3143 cs_alloc_mode_t _mode;
3144
3145};
3146
3147} /* namespace cs */
3148
3149/*--------------------------------------------------------------------------*/
3150/* Usefull aliases without namespace */
3151/*--------------------------------------------------------------------------*/
3152
3153template<class T>
3155
3156template<class T, cs::layout L = cs::layout::right>
3158
3159template<class T, cs::layout L = cs::layout::right>
3161
3162template<class T, cs::layout L = cs::layout::right>
3164
3165template<class T, int N, cs::layout L = cs::layout::right>
3167
3168template<class T, int N>
3170
3171template<class T, int N>
3173
3174#endif /* __cplusplus */
3175
3176#endif /* __CS_ARRAY_H__ */
void bft_error(const char *const file_name, const int line_num, const int sys_error_code, const char *const format,...)
Calls the error handler (set by bft_error_handler_set() or default).
Definition: bft_error.cpp:193
Definition: cs_array.h:1111
CS_F_HOST_DEVICE mdspan< T, _N_, _L_ > get_mdspan(const cs_lnum_t(&dims)[_N_])
Get span view of array with a given array of dimensions (extent). If total size of span is different ...
Definition: cs_array.h:1584
CS_F_HOST void copy_data(cs_dispatch_context &ctx, T *data, const cs_lnum_t n_vals=-1)
Copy data from raw pointer, we suppose that data size is same as that of the array....
Definition: cs_array.h:2792
CS_F_HOST void copy_data(T *data, const cs_lnum_t n_vals=-1)
Copy data from raw pointer, we suppose that data size is same as that of the array.
Definition: cs_array.h:2659
CS_F_HOST_DEVICE void set_alloc_mode(cs_alloc_mode_t mode)
Set memory allocation mode.
Definition: cs_array.h:2200
CS_F_HOST array(cs_lnum_t size1, cs_lnum_t size2, cs_alloc_mode_t alloc_mode=cs_alloc_mode, const char *file_name=__FILE__, const int line_number=__LINE__)
Constructor method for 2D array based on sizes.
Definition: cs_array.h:1216
CS_F_HOST_DEVICE const T * data() const
Const getter to data array.
Definition: cs_array.h:2294
CS_F_HOST_DEVICE cs_alloc_mode_t mode() const
Getter function for allocation mode.
Definition: cs_array.h:2644
CS_F_HOST void set_to_val_on_subset(T val, const cs_lnum_t n_elts, const cs_lnum_t elt_ids[])
Set subset of values of the data array to a constant value.
Definition: cs_array.h:1703
CS_F_HOST void zero(cs_dispatch_context &ctx)
Set all values of the data array to 0.
Definition: cs_array.h:1797
CS_F_HOST void set_to_val(T val, const cs_lnum_t n_vals=-1)
Set all values of the data array to a constant value.
Definition: cs_array.h:1643
CS_F_HOST_DEVICE cs_lnum_t size()
Getter function for total size.
Definition: cs_array.h:2564
CS_F_HOST friend void swap(array &first, array &second)
Class swap operator used for assignment or move.
Definition: cs_array.h:1475
CS_F_HOST void copy_data(cs_dispatch_context &ctx, const T *data, const cs_lnum_t n_vals=-1)
Copy data from const raw pointer, we suppose that data size is same as that of the array....
Definition: cs_array.h:2821
CS_F_HOST_DEVICE bool owner()
Getter function for owner status.
Definition: cs_array.h:2596
CS_F_HOST void reshape(cs_lnum_t new_size, const char *file_name=__FILE__, const int line_number=__LINE__)
Resize data array.
Definition: cs_array.h:1833
CS_F_HOST array(cs_lnum_t size, cs_alloc_mode_t alloc_mode=cs_alloc_mode, const char *file_name=__FILE__, const int line_number=__LINE__)
Constructor method using only size.
Definition: cs_array.h:1147
CS_F_HOST void reshape_and_copy(cs_lnum_t new_size, cs_lnum_t size_to_keep, const char *file_name=__FILE__, const int line_number=__LINE__)
Resize data array and copy old data.
Definition: cs_array.h:1875
CS_F_HOST_DEVICE T * sub_array(Args... indices)
Get sub array based on index.
Definition: cs_array.h:2225
CS_F_HOST_DEVICE mdspan< T, N, L > view()
Get span view of array, same dimensions as array.
Definition: cs_array.h:1544
CS_F_HOST void copy_data(mdspan< T, N, L > &span, const cs_lnum_t n_vals=-1)
Copy data from a span, we suppose that data size is same as that of the array. An assert test the siz...
Definition: cs_array.h:2757
CS_F_HOST_DEVICE cs_lnum_t size() const
Getter function for total size.
Definition: cs_array.h:2580
CS_F_HOST array & operator=(array other)
Assignment operator.
Definition: cs_array.h:1499
CS_F_HOST_DEVICE cs_lnum_t offset(int i) const
Getter function for offset per dimension.
Definition: cs_array.h:2546
CS_F_HOST array(cs_lnum_t size1, cs_lnum_t size2, cs_lnum_t size3, cs_lnum_t size4, cs_alloc_mode_t alloc_mode=cs_alloc_mode, const char *file_name=__FILE__, const int line_number=__LINE__)
Constructor method for 4D array based on sizes.
Definition: cs_array.h:1293
CS_F_HOST void reshape(const cs_lnum_t(&dims)[N], const char *file_name=__FILE__, const int line_number=__LINE__)
Resize data array and do not keep old data.
Definition: cs_array.h:1924
CS_F_HOST void reshape(const cs_lnum_t size1, const cs_lnum_t size2, const char *file_name=__FILE__, const int line_number=__LINE__)
Resize data array 2D (no explicit copy, only realloc)
Definition: cs_array.h:2073
CS_F_HOST_DEVICE cs_lnum_t extent(int i)
Getter function for size per dimension (extent).
Definition: cs_array.h:2489
CS_F_HOST void zero()
Set all values of the data array to 0.
Definition: cs_array.h:1774
CS_F_HOST void set_to_val_on_subset(cs_dispatch_context &ctx, T val, const cs_lnum_t n_elts, const cs_lnum_t elt_ids[])
Set a subset of values of the data array to a constant value while providing a dispatch context....
Definition: cs_array.h:1740
CS_F_HOST_DEVICE array()
Default constructor method leading to "empty container".
Definition: cs_array.h:1122
CS_F_HOST_DEVICE cs_lnum_t offset(int i)
Getter function for offset per dimension.
Definition: cs_array.h:2527
CS_F_HOST array(cs_lnum_t size1, cs_lnum_t size2, cs_lnum_t size3, cs_alloc_mode_t alloc_mode=cs_alloc_mode, const char *file_name=__FILE__, const int line_number=__LINE__)
Constructor method for 3D array based on sizes.
Definition: cs_array.h:1254
CS_F_HOST_DEVICE U * data()
Getter for data raw pointer with recast (casts T* to U*)
Definition: cs_array.h:2281
CS_F_HOST void copy_data(cs_dispatch_context &ctx, mdspan< T, N, L > &span, const cs_lnum_t n_vals=-1)
Copy data from a span, we suppose that data size is same as that of the array. A dispatch_context is ...
Definition: cs_array.h:2882
CS_F_HOST_DEVICE cs_alloc_mode_t mode()
Getter function for allocation mode.
Definition: cs_array.h:2628
CS_F_HOST_DEVICE void clear()
Clear data (empty container).
Definition: cs_array.h:1514
CS_F_HOST_DEVICE T * data()
Getter to data array.
Definition: cs_array.h:2263
CS_F_HOST void reshape_and_copy(const cs_lnum_t size1, const cs_lnum_t size2, cs_lnum_t size_to_keep, const char *file_name=__FILE__, const int line_number=__LINE__)
Resize data array and copy old data (2D)
Definition: cs_array.h:2104
CS_F_HOST void reshape(const cs_lnum_t size1, const cs_lnum_t size2, const cs_lnum_t size3, const char *file_name=__FILE__, const int line_number=__LINE__)
Resize data array 3D (no explicit copy, only realloc)
Definition: cs_array.h:2136
CS_F_HOST array(const cs_lnum_t(&dims)[N], cs_alloc_mode_t alloc_mode=cs_alloc_mode, const char *file_name=__FILE__, const int line_number=__LINE__)
Constructor method using dimensions.
Definition: cs_array.h:1184
CS_F_HOST_DEVICE array(T *data, const cs_lnum_t(&dims)[N])
Constructor method for non owner version in multidimensional.
Definition: cs_array.h:1356
CS_F_HOST_DEVICE mdspan< T, N-sizeof...(Args), L > sub_view(Args... indices)
Get span sub-view of array, with lower dimensionality.
Definition: cs_array.h:1561
CS_F_HOST void reshape_and_copy(const cs_lnum_t(&dims)[N], cs_lnum_t size_to_keep=-1, const char *file_name=__FILE__, const int line_number=__LINE__)
Resize data array and copy old data.
Definition: cs_array.h:1964
CS_F_HOST_DEVICE ~array()
Destructor method.
Definition: cs_array.h:1461
CS_F_HOST_DEVICE cs_lnum_t extent(int i) const
Getter function for size per dimension (extent).
Definition: cs_array.h:2508
CS_F_HOST void copy_data(cs_dispatch_context &ctx, array &other, const cs_lnum_t n_vals=-1)
Copy data from raw pointer, we suppose that data size is same as that of the array....
Definition: cs_array.h:2850
CS_F_HOST void copy_data(const T *data, const cs_lnum_t n_vals=-1)
Copy data from const raw pointer, we suppose that data size is same as that of the array.
Definition: cs_array.h:2691
CS_F_HOST void set_to_val(cs_dispatch_context &ctx, T val, const cs_lnum_t n_vals=-1)
Set all values of the data array to a constant value while providing a dispatch context....
Definition: cs_array.h:1675
CS_F_HOST_DEVICE void set_empty()
Initializer method for empty containers.
Definition: cs_array.h:1817
CS_F_HOST_DEVICE array(const array &other, bool deep_copy=false, const char *file_name=__FILE__, const int line_number=__LINE__)
Constructor method using copy. May be a shallow copy.
Definition: cs_array.h:1390
CS_F_HOST void copy_data(array &other, const cs_lnum_t n_vals=-1)
Copy data from another array, we suppose that data size is same as that of the array....
Definition: cs_array.h:2723
CS_F_HOST_DEVICE bool owner() const
Getter function for owner status.
Definition: cs_array.h:2612
CS_F_HOST_DEVICE T * sub_array(Args... indices) const
Const get sub array based on index.
Definition: cs_array.h:2246
CS_F_HOST array(array &&other)
Move constructor.
Definition: cs_array.h:1446
CS_F_HOST_DEVICE array(T *data, Args... indices)
Constructor method for non owner version.
Definition: cs_array.h:1334
CS_F_HOST_DEVICE void reshape_and_copy(const cs_lnum_t size1, const cs_lnum_t size2, const cs_lnum_t size3, cs_lnum_t size_to_keep, const char *file_name=__FILE__, const int line_number=__LINE__)
Resize data array and copy old data (3D)
Definition: cs_array.h:2168
Definition: cs_mdspan.h:69
auto parallel_for(cs_lnum_t n, F &&f, Args &&... args)
Definition: cs_dispatch.h:1570
void wait(void)
Wait (synchronize) until launched computations have finished.
Definition: cs_dispatch.h:1635
Definition: cs_dispatch.h:1711
void cs_array_real_fill_zero(cs_lnum_t size, cs_real_t a[])
Assign zero to all elements of an array.
Definition: cs_array.cpp:1019
void cs_array_real_set_wvalue(cs_lnum_t n_elts, int stride, const cs_real_t ref_val[], const cs_real_t weight[], cs_real_t a[])
Assign a weighted constant value of dim "stride" to an interlaced array sharing the same stride....
Definition: cs_array.cpp:602
void cs_array_real_set_scalar_on_subset(cs_lnum_t n_elts, const cs_lnum_t elt_ids[], cs_real_t ref_val, cs_real_t a[])
Assign a constant scalar value to an array on a selected subset of elements. If elt_ids is null,...
Definition: cs_array.cpp:766
void cs_array_set_value_real(cs_lnum_t n_elts, cs_lnum_t dim, cs_real_t v, cs_real_t a[])
Assign a constant value to an array (deprecated function).
Definition: cs_array.cpp:1045
void cs_array_bool_fill_false(cs_lnum_t size, bool a[])
Assign false to all elements of an array. Case of an array of booleans.
Definition: cs_array.cpp:106
void cs_arrays_set_zero(cs_dispatch_context &ctx, const cs_lnum_t n_elts, Arrays &&... arrays)
Assign zero value to all elements of multiple arrays.
Definition: cs_array.h:343
void cs_array_real_set_tensor(cs_lnum_t n_elts, const cs_real_t ref_tens[3][3], cs_real_t a[])
Assign a constant 3x3 tensor to an array (of stride 9) which is interlaced.
Definition: cs_array.cpp:950
void cs_array_int_fill_zero(cs_lnum_t size, int a[])
Assign zero to all elements of an array. Case of a int array.
Definition: cs_array.cpp:247
void cs_array_real_set_vector(cs_lnum_t n_elts, const cs_real_t ref_val[3], cs_real_t a[])
Assign a constant vector to an array of stride 3 which is interlaced.
Definition: cs_array.cpp:830
void cs_array_flag_fill_zero(cs_lnum_t size, cs_flag_t a[])
Assign zero to all elements of an array. Case of a cs_flag_t array.
Definition: cs_array.cpp:124
void cs_array_real_scale(cs_lnum_t n_elts, int stride, const cs_lnum_t *elt_ids, cs_real_t scaling_factor, cs_real_t dest[])
Multiply each value by a scaling factor s.t. dest *= scaling_factor If elt_ids is non-null,...
Definition: cs_array.cpp:489
void cs_arrays_set_value_on_subset(const cs_lnum_t n_elts, const cs_lnum_t *elt_ids, const T *ref_val, Arrays &&... arrays)
Assign values on a selected subset of elements to multiple arrays. ref_val is input as a pointer or a...
Definition: cs_array.h:402
void cs_array_difference(cs_lnum_t size, const T *x, const T *y, T *diff)
Compute the difference diff = x - y. All arrays have the same dimension.
Definition: cs_array.h:539
void cs_array_real_set_wvector_on_subset(cs_lnum_t n_elts, const cs_lnum_t elt_ids[], const cs_real_t ref_val[3], const cs_real_t weight[], cs_real_t a[])
Assign a weighted constant vector value to an interlaced array (of stride 3). The subset selection is...
Definition: cs_array.cpp:916
void cs_array_copy(const cs_lnum_t size, const T *src, T *dest)
Copy values from an array to another of the same dimensions.
Definition: cs_array.h:498
void cs_array_real_copy_subset(cs_lnum_t n_elts, int stride, const cs_lnum_t elt_ids[], int mode, const cs_real_t ref[], cs_real_t dest[])
Copy an array ("ref") into another array ("dest") on possibly only a part of the array(s)....
Definition: cs_array.cpp:336
void cs_array_lnum_fill_zero(cs_lnum_t size, cs_lnum_t a[])
Assign zero to all elements of an array. Case of a cs_lnum_t array.
Definition: cs_array.cpp:146
void cs_array_lnum_set_value_on_subset(cs_lnum_t n_elts, const cs_lnum_t elt_ids[], cs_lnum_t num, cs_lnum_t a[])
Assign the value "num" to an array on a selected subset of elements. if elt_ids is null,...
Definition: cs_array.cpp:222
void cs_array_real_set_scalar(cs_lnum_t n_elts, cs_real_t ref_val, cs_real_t a[])
Assign a constant scalar value to an array.
Definition: cs_array.cpp:720
void cs_array_bool_fill_true(cs_lnum_t size, bool a[])
Assign true to all elements of an array. Case of an array of booleans.
Definition: cs_array.cpp:88
void cs_array_real_set_vector_on_subset(cs_lnum_t n_elts, const cs_lnum_t elt_ids[], const cs_real_t ref_val[3], cs_real_t a[])
Assign a constant vector to an interlaced array (of stride 3) on a selected subset of elements....
Definition: cs_array.cpp:881
void cs_array_real_set_value_on_subset(cs_lnum_t n_elts, int stride, const cs_lnum_t elt_ids[], const cs_real_t ref_val[], cs_real_t a[])
Assign a constant value of dim "stride" to an interlaced array sharing the same stride....
Definition: cs_array.cpp:641
void cs_array_real_padd(cs_lnum_t n_elts, const cs_real_t l_add[], cs_real_t r[])
Add in place an array s.t. r += l_add.
Definition: cs_array.cpp:541
void cs_array_lnum_copy(cs_lnum_t size, const cs_lnum_t src[], cs_lnum_t dest[])
Copy values from an array of cs_lnum_t type to another of the same dimensions.
Definition: cs_array.cpp:170
void cs_array_real_copy(cs_lnum_t size, const cs_real_t src[], cs_real_t dest[])
Copy real values from an array to another of the same dimensions.
Definition: cs_array.cpp:457
void cs_array_int_set_value(cs_lnum_t size, int num, int a[])
Assign the value "num" to all elements of an array. Case of a int array.
Definition: cs_array.cpp:271
void cs_array_real_set_value(cs_lnum_t n_elts, int stride, const cs_real_t ref_val[], cs_real_t a[])
Assign a constant value of dim "stride" to an interlaced array sharing the same stride.
Definition: cs_array.cpp:567
void cs_arrays_set_value(const cs_lnum_t n_elts, const T *ref_val, Arrays &&... arrays)
Assign values to all elements of multiple arrays. ref_val is input as a pointer or an array.
Definition: cs_array.h:103
void cs_array_real_set_wscalar_on_subset(cs_lnum_t n_elts, const cs_lnum_t elt_ids[], cs_real_t ref_val, const cs_real_t weight[], cs_real_t a[])
Assign a weighted constant scalar value to an array on a selected subset of elements....
Definition: cs_array.cpp:799
void cs_array_real_set_tensor_on_subset(cs_lnum_t n_elts, const cs_lnum_t elt_ids[], const cs_real_t ref_tens[3][3], cs_real_t a[])
Assign a constant 3x3 tensor to an interlaced array (of stride 9) on a subset of elements....
Definition: cs_array.cpp:985
void cs_array_real_set_wvector(cs_lnum_t n_elts, const cs_real_t ref_val[3], const cs_real_t weight[], cs_real_t a[])
Assign a weighted constant vector value to an interlaced array (of stride 3). The array of weights ha...
Definition: cs_array.cpp:854
void cs_array_lnum_set_value(cs_lnum_t size, cs_lnum_t num, cs_lnum_t a[])
Assign the value "num" to all elements of an array. Case of a cs_lnum_t array.
Definition: cs_array.cpp:199
void cs_array_real_set_wvalue_on_subset(cs_lnum_t n_elts, int stride, const cs_lnum_t elt_ids[], const cs_real_t ref_val[], const cs_real_t weight[], cs_real_t a[])
Assign a weighted constant value of dim "stride" to an interlaced array sharing the same stride....
Definition: cs_array.cpp:682
void cs_array_int_set_value_on_subset(cs_lnum_t n_elts, const cs_lnum_t elt_ids[], int num, int a[])
Assign the value "num" to an array on a selected subset of elements. if elt_ids is null,...
Definition: cs_array.cpp:294
void cs_array_real_set_wscalar(cs_lnum_t n_elts, cs_real_t ref_val, const cs_real_t weight[], cs_real_t a[])
Assign a weighted constant scalar value to an array. The weight array has the same size as the array ...
Definition: cs_array.cpp:742
#define __has_builtin(x)
Definition: cs_array_2dspan.h:48
#define BEGIN_C_DECLS
Definition: cs_defs.h:554
#define CS_F_HOST_DEVICE
Definition: cs_defs.h:585
double cs_real_t
Floating-point value.
Definition: cs_defs.h:357
#define _(String)
Definition: cs_defs.h:67
#define CS_LAMBDA
Definition: cs_defs.h:594
#define CS_THR_MIN
Definition: cs_defs.h:508
#define CS_UNUSED(x)
Definition: cs_defs.h:543
#define END_C_DECLS
Definition: cs_defs.h:555
int cs_lnum_t
local mesh entity id
Definition: cs_defs.h:350
unsigned short int cs_flag_t
Definition: cs_defs.h:359
#define CS_F_HOST
Definition: cs_defs.h:583
@ k
Definition: cs_field_pointer.h:72
#define CS_FREE(_ptr)
Definition: cs_mem.h:155
#define cs_alloc_mode
Definition: cs_mem.h:187
static cs_alloc_mode_t cs_check_device_ptr(const void *ptr)
Check if a pointer is associated with a device.
Definition: cs_mem.h:737
cs_alloc_mode_t
Definition: cs_mem.h:50
@ CS_ALLOC_HOST_DEVICE_SHARED
Definition: cs_mem.h:57
static int cs_parall_n_threads(cs_lnum_t n_elements, cs_lnum_t min_thread_elements)
Compute recommended number of threads for a section.
Definition: cs_parall.h:553
Definition: cs_array.h:1098
layout
Definition: cs_mdspan.h:52
Utility template which returns "true" for a given pack. This is necessary of std::enable_if_t<> when ...
Definition: cs_defs.h:832
Utility template to check if a pack of parameters is made of integral types.
Definition: cs_defs.h:809