32 OutputViewType _output;
33 const inputLeftViewType _inputLeft;
34 const inputRightViewType _inputRight;
36 KOKKOS_INLINE_FUNCTION
38 inputLeftViewType inputLeft_,
39 inputRightViewType inputRight_)
41 _inputLeft(inputLeft_),
42 _inputRight(inputRight_) {}
45 KOKKOS_INLINE_FUNCTION
46 void operator()(
const ordinal_type cl,
const ordinal_type pt)
const
48 const auto val = _inputLeft(cl, pt % _inputLeft.extent(1));
50 if constexpr (equalRank) {
51 if constexpr (reciprocal)
52 _output(cl, pt) = _inputRight(cl, pt) / val;
54 _output(cl, pt) = val * _inputRight(cl, pt);
57 if constexpr (reciprocal)
58 _output(cl, pt) = _inputRight(pt) / val;
60 _output(cl, pt) = val * _inputRight(pt);
65 KOKKOS_INLINE_FUNCTION
66 void operator()(
const ordinal_type cl,
const ordinal_type bf,
const ordinal_type pt)
const
68 const auto val = _inputLeft(cl, pt % _inputLeft.extent(1));
70 if constexpr (equalRank) {
71 if constexpr (reciprocal)
72 _output(cl, bf, pt) = _inputRight(cl, bf, pt) / val;
74 _output(cl, bf, pt) = val * _inputRight(cl, bf, pt);
77 if constexpr (reciprocal)
78 _output(cl, bf, pt) = _inputRight(bf, pt) / val;
80 _output(cl, bf, pt) = val * _inputRight(bf, pt);
92 OutputViewType _output;
93 const inputLeftViewType _inputLeft;
94 const inputRightViewType _inputRight;
96 KOKKOS_INLINE_FUNCTION
98 inputLeftViewType inputLeft_,
99 inputRightViewType inputRight_)
101 _inputLeft(inputLeft_),
102 _inputRight(inputRight_) {}
105 KOKKOS_INLINE_FUNCTION
106 void operator()(
const ordinal_type cl,
const ordinal_type pt)
const
108 const auto val = _inputLeft(cl, pt % _inputLeft.extent(1));
110 if constexpr (equalRank) {
111 if constexpr (reciprocal) {
112 for (ordinal_type i = 0; i < _output.extent_int(2); ++i)
113 _output(cl, pt, i) = _inputRight(cl, pt, i) / val;
116 for (ordinal_type i = 0; i < _output.extent_int(2); ++i)
117 _output(cl, pt, i) = val * _inputRight(cl, pt, i);
121 if constexpr (reciprocal) {
122 for (ordinal_type i = 0; i < _output.extent_int(2); ++i)
123 _output(cl, pt, i) = _inputRight(pt, i) / val;
126 for (ordinal_type i = 0; i < _output.extent_int(2); ++i)
127 _output(cl, pt, i) = val * _inputRight(pt, i);
133 KOKKOS_INLINE_FUNCTION
134 void operator()(
const ordinal_type cl,
const ordinal_type bf,
const ordinal_type pt)
const
136 const auto val = _inputLeft(cl, pt % _inputLeft.extent(1));
138 if constexpr (equalRank) {
139 if constexpr (reciprocal) {
140 for (ordinal_type i = 0; i < _output.extent_int(3); ++i)
141 _output(cl, bf, pt, i) = _inputRight(cl, bf, pt, i) / val;
144 for (ordinal_type i = 0; i < _output.extent_int(3); ++i)
145 _output(cl, bf, pt, i) = val * _inputRight(cl, bf, pt, i);
149 if constexpr (reciprocal) {
150 for (ordinal_type i = 0; i < _output.extent_int(3); ++i)
151 _output(cl, bf, pt, i) = _inputRight(bf, pt, i) / val;
154 for (ordinal_type i = 0; i < _output.extent_int(3); ++i)
155 _output(cl, bf, pt, i) = val * _inputRight(bf, pt, i);
168 OutputViewType _output;
169 const inputLeftViewType _inputLeft;
170 const inputRightViewType _inputRight;
172 KOKKOS_INLINE_FUNCTION
174 inputLeftViewType inputLeft_,
175 inputRightViewType inputRight_)
177 _inputLeft(inputLeft_),
178 _inputRight(inputRight_) {}
181 KOKKOS_INLINE_FUNCTION
182 void operator()(
const ordinal_type cl,
const ordinal_type pt)
const
184 const auto val = _inputLeft(cl, pt % _inputLeft.extent(1));
186 if constexpr (equalRank) {
187 if constexpr (reciprocal) {
188 for (ordinal_type i = 0; i < _output.extent_int(2); ++i)
189 for (ordinal_type j = 0; j < _output.extent_int(3); ++j)
190 _output(cl, pt, i, j) = _inputRight(cl, pt, i, j) / val;
193 for (ordinal_type i = 0; i < _output.extent_int(2); ++i)
194 for (ordinal_type j = 0; j < _output.extent_int(3); ++j)
195 _output(cl, pt, i, j) = val * _inputRight(cl, pt, i, j);
199 if constexpr (reciprocal) {
200 for (ordinal_type i = 0; i < _output.extent_int(2); ++i)
201 for (ordinal_type j = 0; j < _output.extent_int(3); ++j)
202 _output(cl, pt, i, j) = _inputRight(pt, i, j) / val;
205 for (ordinal_type i = 0; i < _output.extent_int(2); ++i)
206 for (ordinal_type j = 0; j < _output.extent_int(3); ++j)
207 _output(cl, pt, i, j) = val * _inputRight(pt, i, j);
213 KOKKOS_INLINE_FUNCTION
214 void operator()(
const ordinal_type cl,
const ordinal_type bf,
const ordinal_type pt)
const
216 const auto val = _inputLeft(cl, pt % _inputLeft.extent(1));
218 if constexpr (equalRank) {
219 if constexpr (reciprocal) {
220 for (ordinal_type i = 0; i < _output.extent_int(3); ++i)
221 for (ordinal_type j = 0; j < _output.extent_int(4); ++j)
222 _output(cl, bf, pt, i, j) = _inputRight(cl, bf, pt, i, j) / val;
225 for (ordinal_type i = 0; i < _output.extent_int(3); ++i)
226 for (ordinal_type j = 0; j < _output.extent_int(4); ++j)
227 _output(cl, bf, pt, i, j) = val * _inputRight(cl, bf, pt, i, j);
231 if constexpr (reciprocal) {
232 for (ordinal_type i = 0; i < _output.extent_int(3); ++i)
233 for (ordinal_type j = 0; j < _output.extent_int(4); ++j)
234 _output(cl, bf, pt, i, j) = _inputRight(bf, pt, i, j) / val;
237 for (ordinal_type i = 0; i < _output.extent_int(3); ++i)
238 for (ordinal_type j = 0; j < _output.extent_int(4); ++j)
239 _output(cl, bf, pt, i, j) = val * _inputRight(bf, pt, i, j);
253 const Kokkos::DynRankView<inputDataValueType, inputDataProperties...> inputData,
254 const Kokkos::DynRankView<inputFieldValueType, inputFieldProperties...> inputFields,
255 const bool reciprocal ) {
257#ifdef HAVE_INTREPID2_DEBUG
259 INTREPID2_TEST_FOR_EXCEPTION( inputData.rank() != 2, std::invalid_argument,
260 ">>> ERROR (ArrayTools::scalarMultiplyDataField): Input data container must have rank 2.");
262 if (outputFields.rank() <= inputFields.rank()) {
263 INTREPID2_TEST_FOR_EXCEPTION( inputFields.rank() < 3 || inputFields.rank() > 5, std::invalid_argument,
264 ">>> ERROR (ArrayTools::scalarMultiplyDataField): Input fields container must have rank 3, 4, or 5.");
265 INTREPID2_TEST_FOR_EXCEPTION( outputFields.rank() != inputFields.rank(), std::invalid_argument,
266 ">>> ERROR (ArrayTools::scalarMultiplyDataField): Input and output fields containers must have the same rank.");
267 INTREPID2_TEST_FOR_EXCEPTION( inputFields.extent(0) != inputData.extent(0), std::invalid_argument,
268 ">>> ERROR (ArrayTools::scalarMultiplyDataField): Zeroth dimensions (number of integration domains) of the fields and data input containers must agree!");
269 INTREPID2_TEST_FOR_EXCEPTION( inputData.extent(1) != inputFields.extent(2) &&
270 inputData.extent(1) != 1, std::invalid_argument,
271 ">>> ERROR (ArrayTools::scalarMultiplyDataField): Second dimension of the fields input container and first dimension of data input container (number of integration points) must agree, or first data dimension must be 1!");
272 for (size_type i=0;i<inputFields.rank();++i) {
273 INTREPID2_TEST_FOR_EXCEPTION( inputFields.extent(i) != outputFields.extent(i), std::invalid_argument,
274 ">>> ERROR (ArrayTools::scalarMultiplyDataField): inputFields dimension (i) does not match to the dimension (i) of outputFields");
278 INTREPID2_TEST_FOR_EXCEPTION( inputFields.rank() < 2 || inputFields.rank() > 4, std::invalid_argument,
279 ">>> ERROR (ArrayTools::scalarMultiplyDataField): Input fields container must have rank 2, 3, or 4.");
280 INTREPID2_TEST_FOR_EXCEPTION( outputFields.rank() != (inputFields.rank()+1), std::invalid_argument,
281 ">>> ERROR (ArrayTools::scalarMultiplyDataField): The rank of the input fields container must be one less than the rank of the output fields container.");
282 INTREPID2_TEST_FOR_EXCEPTION( inputData.extent(1) != inputFields.extent(1) &&
283 inputData.extent(1) != 1, std::invalid_argument,
284 ">>> ERROR (ArrayTools::scalarMultiplyDataField): First dimensions of fields input container and data input container (number of integration points) must agree or first data dimension must be 1!");
285 INTREPID2_TEST_FOR_EXCEPTION( inputData.extent(0) != outputFields.extent(0), std::invalid_argument,
286 ">>> ERROR (ArrayTools::scalarMultiplyDataField): Zeroth dimensions of fields output container and data input containers (number of integration domains) must agree!");
287 for (size_type i=0;i<inputFields.rank();++i) {
288 INTREPID2_TEST_FOR_EXCEPTION( inputFields.extent(i) != outputFields.extent(i+1), std::invalid_argument,
289 ">>> ERROR (ArrayTools::scalarMultiplyDataField): inputFields dimension (i) does not match to the dimension (i+1) of outputFields");
295 typedef Kokkos::DynRankView<outputFieldValueType,outputFieldProperties...> outputFieldViewType;
296 typedef Kokkos::DynRankView<inputDataValueType,inputDataProperties...> inputDataViewType;
297 typedef Kokkos::DynRankView<inputFieldValueType,inputFieldProperties...> inputFieldViewType;
299 using range_policy_type = Kokkos::MDRangePolicy
300 < ExecSpaceType, Kokkos::Rank<3>, Kokkos::IndexType<ordinal_type> >;
302 const range_policy_type policy( { 0, 0, 0 },
303 { outputFields.extent(0), outputFields.extent(1), outputFields.extent(2) } );
305 const bool equalRank = ( outputFields.rank() == inputFields.rank() );
307 if (outputFields.rank() == 3) {
311 Kokkos::parallel_for( policy, FunctorType(outputFields, inputData, inputFields) );
314 Kokkos::parallel_for( policy, FunctorType(outputFields, inputData, inputFields) );
320 Kokkos::parallel_for( policy, FunctorType(outputFields, inputData, inputFields) );
323 Kokkos::parallel_for( policy, FunctorType(outputFields, inputData, inputFields) );
327 else if (outputFields.rank() == 4) {
331 Kokkos::parallel_for( policy, FunctorType(outputFields, inputData, inputFields) );
334 Kokkos::parallel_for( policy, FunctorType(outputFields, inputData, inputFields) );
340 Kokkos::parallel_for( policy, FunctorType(outputFields, inputData, inputFields) );
343 Kokkos::parallel_for( policy, FunctorType(outputFields, inputData, inputFields) );
351 Kokkos::parallel_for( policy, FunctorType(outputFields, inputData, inputFields) );
354 Kokkos::parallel_for( policy, FunctorType(outputFields, inputData, inputFields) );
360 Kokkos::parallel_for( policy, FunctorType(outputFields, inputData, inputFields) );
363 Kokkos::parallel_for( policy, FunctorType(outputFields, inputData, inputFields) );
377 const Kokkos::DynRankView<inputDataLeftValueType, inputDataLeftProperties...> inputDataLeft,
378 const Kokkos::DynRankView<inputDataRightValueType,inputDataRightProperties...> inputDataRight,
379 const bool reciprocal ) {
381#ifdef HAVE_INTREPID2_DEBUG
383 INTREPID2_TEST_FOR_EXCEPTION( inputDataLeft.rank() != 2, std::invalid_argument,
384 ">>> ERROR (ArrayTools::scalarMultiplyDataData): Left input data container must have rank 2.");
386 if (outputData.rank() <= inputDataRight.rank()) {
387 INTREPID2_TEST_FOR_EXCEPTION( inputDataRight.rank() < 2 || inputDataRight.rank() > 4, std::invalid_argument,
388 ">>> ERROR (ArrayTools::scalarMultiplyDataData): Right input data container must have rank 2, 3, or 4.");
389 INTREPID2_TEST_FOR_EXCEPTION( outputData.rank() != inputDataRight.rank(), std::invalid_argument,
390 ">>> ERROR (ArrayTools::scalarMultiplyDataData): Right input and output data containers must have the same rank.");
391 INTREPID2_TEST_FOR_EXCEPTION( inputDataRight.extent(0) != inputDataLeft.extent(0), std::invalid_argument,
392 ">>> ERROR (ArrayTools::scalarMultiplyDataData): Zeroth dimensions (number of integration domains) of the left and right data input containers must agree!");
393 INTREPID2_TEST_FOR_EXCEPTION( inputDataLeft.extent(1) != inputDataRight.extent(1) &&
394 inputDataLeft.extent(1) != 1, std::invalid_argument,
395 ">>> ERROR (ArrayTools::scalarMultiplyDataData): First dimensions of the left and right data input containers (number of integration points) must agree, or first dimension of the left data input container must be 1!");
396 for (size_type i=0;i<inputDataRight.rank();++i) {
397 INTREPID2_TEST_FOR_EXCEPTION( inputDataRight.extent(i) != outputData.extent(i), std::invalid_argument,
398 ">>> ERROR (ArrayTools::scalarMultiplyDataData): inputDataRight dimension (i) does not match to the dimension (i) of outputData");
401 INTREPID2_TEST_FOR_EXCEPTION( inputDataRight.rank() < 1 || inputDataRight.rank() > 3, std::invalid_argument,
402 ">>> ERROR (ArrayTools::scalarMultiplyDataData): Right input data container must have rank 1, 2, or 3.");
403 INTREPID2_TEST_FOR_EXCEPTION( outputData.rank() != (inputDataRight.rank()+1), std::invalid_argument,
404 ">>> ERROR (ArrayTools::scalarMultiplyDataData): The rank of the right input data container must be one less than the rank of the output data container.");
405 INTREPID2_TEST_FOR_EXCEPTION( inputDataLeft.extent(1) != inputDataRight.extent(0) &&
406 inputDataLeft.extent(1) != 1, std::invalid_argument,
407 ">>> ERROR (ArrayTools::scalarMultiplyDataData): Zeroth dimension of the right input data container and first dimension of the left data input container (number of integration points) must agree or first dimension of the left data input container must be 1!");
408 INTREPID2_TEST_FOR_EXCEPTION( inputDataLeft.extent(0) != outputData.extent(0), std::invalid_argument,
409 ">>> ERROR (ArrayTools::scalarMultiplyDataData): Zeroth dimensions of data output and left data input containers (number of integration domains) must agree!");
410 for (size_type i=0;i<inputDataRight.rank();++i) {
411 INTREPID2_TEST_FOR_EXCEPTION( inputDataRight.extent(i) != outputData.extent(i+1), std::invalid_argument,
412 ">>> ERROR (ArrayTools::scalarMultiplyDataData): inputDataRight dimension (i) does not match to the dimension (i+1) of outputData");
418 typedef Kokkos::DynRankView<outputDataValueType,outputDataProperties...> outputDataViewType;
419 typedef Kokkos::DynRankView<inputDataLeftValueType,inputDataLeftProperties...> inputDataLeftViewType;
420 typedef Kokkos::DynRankView<inputDataRightValueType,inputDataRightProperties...> inputDataRightViewType;
422 using range_policy_type = Kokkos::MDRangePolicy
423 < ExecSpaceType, Kokkos::Rank<2>, Kokkos::IndexType<ordinal_type> >;
425 const range_policy_type policy( { 0, 0 },
426 { outputData.extent(0), outputData.extent(1) } );
428 const bool equalRank = ( outputData.rank() == inputDataRight.rank() );
429 if (outputData.rank() == 2) {
433 Kokkos::parallel_for( policy, FunctorType(outputData, inputDataLeft, inputDataRight) );
436 Kokkos::parallel_for( policy, FunctorType(outputData, inputDataLeft, inputDataRight) );
442 Kokkos::parallel_for( policy, FunctorType(outputData, inputDataLeft, inputDataRight) );
445 Kokkos::parallel_for( policy, FunctorType(outputData, inputDataLeft, inputDataRight) );
449 else if (outputData.rank() == 3) {
453 Kokkos::parallel_for( policy, FunctorType(outputData, inputDataLeft, inputDataRight) );
456 Kokkos::parallel_for( policy, FunctorType(outputData, inputDataLeft, inputDataRight) );
462 Kokkos::parallel_for( policy, FunctorType(outputData, inputDataLeft, inputDataRight) );
465 Kokkos::parallel_for( policy, FunctorType(outputData, inputDataLeft, inputDataRight) );
473 Kokkos::parallel_for( policy, FunctorType(outputData, inputDataLeft, inputDataRight) );
476 Kokkos::parallel_for( policy, FunctorType(outputData, inputDataLeft, inputDataRight) );
482 Kokkos::parallel_for( policy, FunctorType(outputData, inputDataLeft, inputDataRight) );
485 Kokkos::parallel_for( policy, FunctorType(outputData, inputDataLeft, inputDataRight) );