Compadre 1.6.4
Loading...
Searching...
No Matches
Compadre_Functors.hpp
Go to the documentation of this file.
1// @HEADER
2// *****************************************************************************
3// Compadre: COMpatible PArticle Discretization and REmap Toolkit
4//
5// Copyright 2018 NTESS and the Compadre contributors.
6// SPDX-License-Identifier: BSD-2-Clause
7// *****************************************************************************
8// @HEADER
9#ifndef _COMPADRE_FUNCTORS_HPP_
10#define _COMPADRE_FUNCTORS_HPP_
11
14#include "Compadre_Basis.hpp"
16#include "Compadre_Targets.hpp"
18#include "Compadre_Functors.hpp"
20#include "KokkosBatched_Gemm_Decl.hpp"
21#include "Compadre_GMLS.hpp"
22
23namespace Compadre {
24
26
27 Kokkos::View<double**, layout_right> _source_extra_data;
28 Kokkos::View<double**, layout_right> _target_extra_data;
29 Kokkos::View<double*> _epsilons;
30 Kokkos::View<double*****, layout_right> _prestencil_weights;
31 Kokkos::View<TargetOperation*> _curvature_support_operations;
32 Kokkos::View<TargetOperation*> _operations;
33
36 int _NP;
51
65
66 // convenience variables (not from GMLS class)
68 double * RHS_data;
70 double * P_data;
74 double * Coeffs_data;
75 double * w_data;
82 double * T_data;
84 double * ref_N_data;
88
90};
91
109
110/** @name Functors
111 * Member functions that perform operations on the entire batch
112 */
113///@{
114
115//! Functor to apply target evaluation to polynomial coefficients to store in _alphas
117
119
121
122 KOKKOS_INLINE_FUNCTION
123 void operator()(const member_type& teamMember) const {
124
125 const int local_index = teamMember.league_rank();
126
127 /*
128 * Data
129 */
130
131 // Coefficients for polynomial basis have overwritten _data._RHS
138
139 applyTargetsToCoefficients(_data, teamMember, Coeffs, P_target_row);
140 }
141};
142
143//! Functor to evaluate targets operations on the basis
145
147
149
150 KOKKOS_INLINE_FUNCTION
151 void operator()(const member_type& teamMember) const {
152 /*
153 * Dimensions
154 */
155
156 const int local_index = teamMember.league_rank();
157
158 /*
159 * Data
160 */
161
165
166 scratch_vector_type delta(teamMember.thread_scratch(_data._pm.getThreadScratchLevel(0)),
168 scratch_vector_type thread_workspace(teamMember.thread_scratch(_data._pm.getThreadScratchLevel(0)),
170
171 /*
172 * Evaluate Standard Targets
173 */
174
175 computeTargetFunctionals(_data, teamMember, delta, thread_workspace, P_target_row);
176
177 }
178};
179
180//! Functor to calculate prestencil weights to apply to data to transform into a format expected by a GMLS stencil
182
184
186
187 KOKKOS_INLINE_FUNCTION
188 void operator()(const member_type& teamMember) const {
189
190 /*
191 * Dimensions
192 */
193
194 const int target_index = _data._initial_index_for_batch + teamMember.league_rank();
195 const int local_index = teamMember.league_rank();
196 const int dimensions = _data._dimensions;
197
198 /*
199 * Data
200 */
201
202 scratch_vector_type delta(teamMember.thread_scratch(_data._pm.getThreadScratchLevel(0)),
204 scratch_vector_type thread_workspace(teamMember.thread_scratch(_data._pm.getThreadScratchLevel(0)),
206
208 scratch_vector_type t1, t2;
210 tangent = scratch_matrix_right_type(teamMember.thread_scratch(_data._pm.getThreadScratchLevel(1)),
211 dimensions-1, dimensions);
212 t1 = scratch_vector_type(teamMember.team_scratch(_data._pm.getTeamScratchLevel(1)),
214 t2 = scratch_vector_type(teamMember.team_scratch(_data._pm.getTeamScratchLevel(1)),
216 for (int j = 0; j < delta.extent_int(0); ++j) {
217 delta(j) = 0;
218 }
219 for (int j = 0; j < thread_workspace.extent_int(0); ++j) {
220 thread_workspace(j) = 0;
221 }
222 }
223
224
225 // holds polynomial coefficients for curvature reconstruction
229
231 + TO_GLOBAL(target_index)*TO_GLOBAL(dimensions)*TO_GLOBAL(dimensions),
232 dimensions, dimensions);
233
234 scratch_vector_type manifold_gradient(teamMember.team_scratch(_data._pm.getTeamScratchLevel(1)),
236
237 /*
238 * Prestencil Weight Calculations
239 */
240
242 Kokkos::single(Kokkos::PerTeam(teamMember), [&] () {
243 Kokkos::single(Kokkos::PerThread(teamMember), [&] () {
244 _data._prestencil_weights(0,0,0,0,0) = -1;
245 _data._prestencil_weights(1,0,0,0,0) = 1;
246 });
247 });
249 Kokkos::single(Kokkos::PerTeam(teamMember), [&] () {
250 Kokkos::single(Kokkos::PerThread(teamMember), [&] () {
251 for (int j=0; j<dimensions; ++j) {
252 for (int k=0; k<dimensions-1; ++k) {
253 _data._prestencil_weights(0,target_index,0,k,j) = T(k,j);
254 }
255 }
256 });
257 });
260 && "StaggeredEdgeIntegralSample prestencil weight only written for manifolds.");
261 Kokkos::parallel_for(Kokkos::TeamThreadRange(teamMember,
262 _data._pc._nla.getNumberOfNeighborsDevice(target_index)), [&] (const int m) {
263 Kokkos::single(Kokkos::PerThread(teamMember), [&] () {
264 for (int quadrature = 0; quadrature<_data._qm.getNumberOfQuadraturePoints(); ++quadrature) {
265 XYZ tangent_quadrature_coord_2d;
266 for (int j=0; j<dimensions-1; ++j) {
267 tangent_quadrature_coord_2d[j] = _data._pc.getTargetCoordinate(target_index, j, &T);
268 tangent_quadrature_coord_2d[j] -= _data._pc.getNeighborCoordinate(target_index, m, j, &T);
269 }
270 double tangent_vector[3];
271 tangent_vector[0] = tangent_quadrature_coord_2d[0]*T(0,0) + tangent_quadrature_coord_2d[1]*T(1,0);
272 tangent_vector[1] = tangent_quadrature_coord_2d[0]*T(0,1) + tangent_quadrature_coord_2d[1]*T(1,1);
273 tangent_vector[2] = tangent_quadrature_coord_2d[0]*T(0,2) + tangent_quadrature_coord_2d[1]*T(1,2);
274
275 for (int j=0; j<dimensions; ++j) {
276 _data._prestencil_weights(0,target_index,m,0,j) += (1-_data._qm.getSite(quadrature,0))*tangent_vector[j]*_data._qm.getWeight(quadrature);
277 _data._prestencil_weights(1,target_index,m,0,j) += _data._qm.getSite(quadrature,0)*tangent_vector[j]*_data._qm.getWeight(quadrature);
278 }
279 }
280 });
281 });
283
284 Kokkos::parallel_for(Kokkos::TeamThreadRange(teamMember,
285 _data._pc._nla.getNumberOfNeighborsDevice(target_index)), [&] (const int m) {
286
287 Kokkos::single(Kokkos::PerThread(teamMember), [&] () {
288 calcGradientPij(_data, teamMember, delta.data(), thread_workspace.data(), target_index,
289 m, 0 /*alpha*/, 0 /*partial_direction*/, dimensions-1, _data._curvature_poly_order,
290 false /*specific order only*/, &T, ReconstructionSpace::ScalarTaylorPolynomial,
291 PointSample);
292 });
293 // reconstructs gradient at local neighbor index m
294 double grad_xi1 = 0, grad_xi2 = 0;
295 Kokkos::parallel_reduce(Kokkos::ThreadVectorRange(teamMember,
296 _data._pc._nla.getNumberOfNeighborsDevice(target_index)), [&] (const int i, double &t_grad_xi1) {
297 double alpha_ij = 0;
298 for (int l=0; l<_data.manifold_NP; ++l) {
299 alpha_ij += delta(l)*Q(l,i);
300 }
301 XYZ rel_coord = _data._pc.getRelativeCoord(target_index, i, dimensions, &T);
302 double normal_coordinate = rel_coord[dimensions-1];
303
304 // apply coefficients to sample data
305 t_grad_xi1 += alpha_ij * normal_coordinate;
306 }, grad_xi1);
307 t1(m) = grad_xi1;
308
309 Kokkos::single(Kokkos::PerThread(teamMember), [&] () {
310 calcGradientPij(_data, teamMember, delta.data(), thread_workspace.data(), target_index,
311 m, 0 /*alpha*/, 1 /*partial_direction*/, dimensions-1, _data._curvature_poly_order,
312 false /*specific order only*/, &T, ReconstructionSpace::ScalarTaylorPolynomial,
314 });
315 Kokkos::parallel_reduce(Kokkos::ThreadVectorRange(teamMember,
316 _data._pc._nla.getNumberOfNeighborsDevice(target_index)), [&] (const int i, double &t_grad_xi2) {
317 double alpha_ij = 0;
318 for (int l=0; l<_data.manifold_NP; ++l) {
319 alpha_ij += delta(l)*Q(l,i);
320 }
321 XYZ rel_coord = _data._pc.getRelativeCoord(target_index, i, dimensions, &T);
322 double normal_coordinate = rel_coord[dimensions-1];
323
324 // apply coefficients to sample data
325 if (dimensions>2) t_grad_xi2 += alpha_ij * normal_coordinate;
326 }, grad_xi2);
327 t2(m) = grad_xi2;
328 });
329
330 teamMember.team_barrier();
331
332 Kokkos::parallel_for(Kokkos::TeamThreadRange(teamMember,
333 _data._pc._nla.getNumberOfNeighborsDevice(target_index)), [&] (const int m) {
334 // constructs tangent vector at neighbor site
335 Kokkos::single(Kokkos::PerThread(teamMember), [&] () {
336 for (int j=0; j<dimensions; ++j) {
337 tangent(0,j) = t1(m)*T(dimensions-1,j) + T(0,j);
338 tangent(1,j) = t2(m)*T(dimensions-1,j) + T(1,j);
339 }
340
341 // calculate norm
342 double norm = 0;
343 for (int j=0; j<dimensions; ++j) {
344 norm += tangent(0,j)*tangent(0,j);
345 }
346
347 // normalize first vector
348 norm = std::sqrt(norm);
349 for (int j=0; j<dimensions; ++j) {
350 tangent(0,j) /= norm;
351 }
352
353 // orthonormalize next vector
354 if (dimensions-1 == 2) { // 2d manifold
355 double dot_product = tangent(0,0)*tangent(1,0)
356 + tangent(0,1)*tangent(1,1)
357 + tangent(0,2)*tangent(1,2);
358 for (int j=0; j<dimensions; ++j) {
359 tangent(1,j) -= dot_product*tangent(0,j);
360 }
361 // normalize second vector
362 norm = 0;
363 for (int j=0; j<dimensions; ++j) {
364 norm += tangent(1,j)*tangent(1,j);
365 }
366 norm = std::sqrt(norm);
367 for (int j=0; j<dimensions; ++j) {
368 tangent(1,j) /= norm;
369 }
370 }
371
372 // stores matrix of tangent and normal directions as a prestencil weight
373 for (int j=0; j<dimensions; ++j) {
374 for (int k=0; k<dimensions-1; ++k) {
375 _data._prestencil_weights(0,target_index,m,k,j) = tangent(k,j);
376 }
377 }
378 });
379 });
380 }
381 teamMember.team_barrier();
382 }
383};
384
385//! Functor to assemble the P*sqrt(weights) matrix and construct sqrt(weights)*Identity
387
389
391
392 KOKKOS_INLINE_FUNCTION
393 void operator()(const member_type& teamMember) const {
394
395 /*
396 * Dimensions
397 */
398
399 const int target_index = _data._initial_index_for_batch + teamMember.league_rank();
400 const int local_index = teamMember.league_rank();
401 const int this_num_rows = _data._sampling_multiplier*_data._pc._nla.getNumberOfNeighborsDevice(target_index);
402
403 /*
404 * Data
405 */
406
408 + TO_GLOBAL(local_index)*TO_GLOBAL(_data.P_dim_0*_data.P_dim_1),
409 _data.P_dim_0, _data.P_dim_1);
411 + TO_GLOBAL(local_index)*TO_GLOBAL(_data.RHS_dim_0*_data.RHS_dim_1),
412 _data.RHS_dim_0, _data.RHS_dim_1);
414 + TO_GLOBAL(local_index)*TO_GLOBAL(_data.max_num_rows),
415 _data.max_num_rows);
416
417 // delta, used for each thread
418 scratch_vector_type delta(teamMember.thread_scratch(_data._pm.getThreadScratchLevel(0)),
419 _data.this_num_cols);
420 scratch_vector_type thread_workspace(teamMember.thread_scratch(_data._pm.getThreadScratchLevel(0)),
422
423 /*
424 * Assemble P*sqrt(W) and sqrt(w)*Identity
425 */
426
427 // creates the matrix sqrt(W)*P
428 createWeightsAndP(_data, teamMember, delta, thread_workspace, PsqrtW, w, _data._dimensions,
429 _data._poly_order, true /*weight_p*/, NULL /*&V*/, _data._reconstruction_space,
431
432 if ((_data._constraint_type == ConstraintType::NO_CONSTRAINT) && (_data._dense_solver_type != DenseSolverType::LU)) {
433 // fill in RHS with Identity * sqrt(weights)
434 double * rhs_data = RHS.data();
435 Kokkos::parallel_for(Kokkos::TeamVectorRange(teamMember, this_num_rows), [&] (const int i) {
436 rhs_data[i] = std::sqrt(w(i));
437 });
438 } else {
439 // create global memory for matrix M = PsqrtW^T*PsqrtW
440 // don't need to cast into scratch_matrix_left_type since the matrix is symmetric
442 + TO_GLOBAL(local_index)*TO_GLOBAL(_data.RHS_dim_0*_data.RHS_dim_1),
443 _data.RHS_dim_0, _data.RHS_dim_1);
444 KokkosBatched::TeamVectorGemm<member_type,KokkosBatched::Trans::Transpose,KokkosBatched::Trans::NoTranspose,KokkosBatched::Algo::Gemm::Unblocked>
445 ::invoke(teamMember,
446 1.0,
447 PsqrtW,
448 PsqrtW,
449 0.0,
450 M);
451 teamMember.team_barrier();
452
453 // Multiply PsqrtW with sqrt(W) to get PW
454 Kokkos::parallel_for(Kokkos::TeamThreadRange(teamMember, _data.max_num_rows), [&] (const int i) {
455 for (int j=0; j < _data.this_num_cols; j++) {
456 PsqrtW(i,j) = PsqrtW(i,j)*std::sqrt(w(i));
457 }
458 });
459 teamMember.team_barrier();
460
461 // conditionally fill in rows determined by constraint type
462 if (_data._constraint_type == ConstraintType::NEUMANN_GRAD_SCALAR) {
463 // normal vector is contained in last row of T
465 + TO_GLOBAL(target_index)*TO_GLOBAL(_data._dimensions*_data._dimensions),
466 _data._dimensions, _data._dimensions);
467
468 // Get the number of neighbors for target index
469 int num_neigh_target = _data._pc._nla.getNumberOfNeighborsDevice(target_index);
470 double cutoff_p = _data._epsilons(target_index);
471
473 _data._NP, cutoff_p, _data._dimensions, num_neigh_target, &T);
474 }
475 }
476 teamMember.team_barrier();
477 }
478};
479
480//! Functor to create a coarse tangent approximation from a given neighborhood of points
482
484
485 // random number generator pool
487
489 // seed random number generator pool
490 _random_number_pool = pool_type(1);
491 }
492
493 KOKKOS_INLINE_FUNCTION
494 void operator()(const member_type& teamMember) const {
495
496 /*
497 * Dimensions
498 */
499
500 const int target_index = _data._initial_index_for_batch + teamMember.league_rank();
501 const int local_index = teamMember.league_rank();
502 const int dimensions = _data._dimensions;
503
504 /*
505 * Data
506 */
507
509 + TO_GLOBAL(local_index)*TO_GLOBAL(_data.P_dim_0*_data.P_dim_1),
510 _data.P_dim_0, _data.P_dim_1);
512 + TO_GLOBAL(local_index)*TO_GLOBAL(_data.max_num_rows),
513 _data.max_num_rows);
515 + TO_GLOBAL(target_index)*TO_GLOBAL(dimensions*dimensions),
516 dimensions, dimensions);
517
518 scratch_matrix_right_type PTP(teamMember.team_scratch(_data._pm.getTeamScratchLevel(1)),
519 dimensions, dimensions);
520
521 // delta, used for each thread
522 scratch_vector_type delta(teamMember.thread_scratch(_data._pm.getThreadScratchLevel(0)),
523 _data.this_num_cols);
524 scratch_vector_type thread_workspace(teamMember.thread_scratch(_data._pm.getThreadScratchLevel(0)),
526
527 /*
528 * Determine Coarse Approximation of Manifold Tangent Plane
529 */
530
531 // getting x y and z from which to derive a manifold
532 createWeightsAndPForCurvature(_data, teamMember, delta, thread_workspace, PsqrtW, w, dimensions, true /* only specific order */);
533
534 // create PsqrtW^T*PsqrtW
535 KokkosBatched::TeamVectorGemm<member_type,KokkosBatched::Trans::Transpose,KokkosBatched::Trans::NoTranspose,KokkosBatched::Algo::Gemm::Unblocked>
536 ::invoke(teamMember,
537 1.0,
538 PsqrtW,
539 PsqrtW,
540 0.0,
541 PTP);
542 teamMember.team_barrier();
543
544 // create coarse approximation of tangent plane in first two rows of T, with normal direction in third column
545 GMLS_LinearAlgebra::largestTwoEigenvectorsThreeByThreeSymmetric(teamMember, T, PTP, dimensions,
546 const_cast<pool_type&>(_random_number_pool));
547
548 teamMember.team_barrier();
549 }
550};
551
552//! Functor to assemble the P*sqrt(weights) matrix and construct sqrt(weights)*Identity for curvature
554
556
558
559 KOKKOS_INLINE_FUNCTION
560 void operator()(const member_type& teamMember) const {
561
562 /*
563 * Dimensions
564 */
565
566 const int target_index = _data._initial_index_for_batch + teamMember.league_rank();
567 const int local_index = teamMember.league_rank();
568 const int this_num_neighbors = _data._pc._nla.getNumberOfNeighborsDevice(target_index);
569
570 /*
571 * Data
572 */
573
574 scratch_matrix_right_type CurvaturePsqrtW(_data.P_data
575 + TO_GLOBAL(local_index)*TO_GLOBAL(_data.P_dim_0*_data.P_dim_1),
576 _data.P_dim_0, _data.P_dim_1);
578 + TO_GLOBAL(local_index)*TO_GLOBAL(_data.RHS_dim_0*_data.RHS_dim_1),
579 _data.RHS_dim_0, _data.RHS_dim_1);
581 + TO_GLOBAL(local_index)*TO_GLOBAL(_data.max_num_rows), _data.max_num_rows);
583 + TO_GLOBAL(target_index)*TO_GLOBAL(_data._dimensions*_data._dimensions),
584 _data._dimensions, _data._dimensions);
585
586 // delta, used for each thread
587 scratch_vector_type delta(teamMember.thread_scratch(_data._pm.getThreadScratchLevel(0)),
588 _data.this_num_cols);
589 scratch_vector_type thread_workspace(teamMember.thread_scratch(_data._pm.getThreadScratchLevel(0)),
591
592 //
593 // RECONSTRUCT ON THE TANGENT PLANE USING LOCAL COORDINATES
594 //
595
596 // creates the matrix sqrt(W)*P
597 createWeightsAndPForCurvature(_data, teamMember, delta, thread_workspace, CurvaturePsqrtW, w,
598 _data._dimensions-1, false /* only specific order */, &T);
599
600 // CurvaturePsqrtW is sized according to max_num_rows x this_num_cols of which in this case
601 // we are only using this_num_neighbors x manifold_NP
602 if (_data._dense_solver_type != DenseSolverType::LU) {
603 // fill in RHS with Identity * sqrt(weights)
604 double * rhs_data = RHS.data();
605 Kokkos::parallel_for(Kokkos::TeamVectorRange(teamMember, this_num_neighbors), [&] (const int i) {
606 rhs_data[i] = std::sqrt(w(i));
607 });
608 } else {
609 // create global memory for matrix M = PsqrtW^T*PsqrtW
610 // don't need to cast into scratch_matrix_left_type since the matrix is symmetric
612 + TO_GLOBAL(local_index)*TO_GLOBAL(_data.RHS_dim_0*_data.RHS_dim_1),
613 _data.RHS_dim_0, _data.RHS_dim_1);
614 // Assemble matrix M
615 KokkosBatched::TeamVectorGemm<member_type,KokkosBatched::Trans::Transpose,
616 KokkosBatched::Trans::NoTranspose,KokkosBatched::Algo::Gemm::Unblocked>
617 ::invoke(teamMember,
618 1.0,
619 CurvaturePsqrtW,
620 CurvaturePsqrtW,
621 0.0,
622 M);
623 teamMember.team_barrier();
624
625 // Multiply PsqrtW with sqrt(W) to get PW
626 Kokkos::parallel_for(Kokkos::TeamThreadRange(teamMember, this_num_neighbors), [&] (const int i) {
627 for (int j=0; j < _data.manifold_NP; j++) {
628 CurvaturePsqrtW(i, j) = CurvaturePsqrtW(i, j)*std::sqrt(w(i));
629 }
630 });
631 }
632 teamMember.team_barrier();
633 }
634};
635
636//! Functor to evaluate curvature targets and construct accurate tangent direction approximation for manifolds
638
640
642
643 KOKKOS_INLINE_FUNCTION
644 void operator()(const member_type& teamMember) const {
645
646 /*
647 * Dimensions
648 */
649
650 const int target_index = _data._initial_index_for_batch + teamMember.league_rank();
651 const int local_index = teamMember.league_rank();
652 auto dimensions = _data._dimensions;
653
654 /*
655 * Data
656 */
657
659 + TO_GLOBAL(local_index)*TO_GLOBAL(_data.Coeffs_dim_0*_data.Coeffs_dim_1),
660 _data.Coeffs_dim_0, _data.Coeffs_dim_1);
662 + TO_GLOBAL(target_index)*TO_GLOBAL(dimensions*dimensions),
663 dimensions, dimensions);
665 + TO_GLOBAL(local_index)*TO_GLOBAL(_data.P_target_row_dim_0*_data.P_target_row_dim_1),
667
668 scratch_vector_type manifold_gradient(teamMember.team_scratch(_data._pm.getTeamScratchLevel(1)),
670
671 // delta, used for each thread
672 scratch_vector_type delta(teamMember.thread_scratch(_data._pm.getThreadScratchLevel(0)),
673 _data.this_num_cols);
674 scratch_vector_type thread_workspace(teamMember.thread_scratch(_data._pm.getThreadScratchLevel(0)),
676
677 /*
678 * Manifold
679 */
680
681
682 //
683 // GET TARGET COEFFICIENTS RELATED TO GRADIENT TERMS
684 //
685 // reconstruct grad_xi1 and grad_xi2, not used for manifold_coeffs
686 computeCurvatureFunctionals(_data, teamMember, delta, thread_workspace, P_target_row, &T);
687 teamMember.team_barrier();
688
689 double grad_xi1 = 0, grad_xi2 = 0;
690 for (int i=0; i<_data._pc._nla.getNumberOfNeighborsDevice(target_index); ++i) {
691 for (int k=0; k<dimensions-1; ++k) {
692 double alpha_ij = 0;
693 int offset = _data._d_ss.getTargetOffsetIndex(0, 0, k, 0);
694 Kokkos::parallel_reduce(Kokkos::TeamThreadRange(teamMember,
695 _data.manifold_NP), [&] (const int l, double &talpha_ij) {
696 Kokkos::single(Kokkos::PerThread(teamMember), [&] () {
697 talpha_ij += P_target_row(offset,l)*Q(l,i);
698 });
699 }, alpha_ij);
700 teamMember.team_barrier();
701 Kokkos::single(Kokkos::PerTeam(teamMember), [&] () {
702 // stored staggered, grad_xi1, grad_xi2, grad_xi1, grad_xi2, ....
703 manifold_gradient(i*(dimensions-1) + k) = alpha_ij;
704 });
705 }
706 teamMember.team_barrier();
707
708 XYZ rel_coord = _data._pc.getRelativeCoord(target_index, i, dimensions, &T);
709 double normal_coordinate = rel_coord[dimensions-1];
710
711 // apply coefficients to sample data
712 grad_xi1 += manifold_gradient(i*(dimensions-1)) * normal_coordinate;
713 if (dimensions>2) grad_xi2 += manifold_gradient(i*(dimensions-1)+1) * normal_coordinate;
714 teamMember.team_barrier();
715 }
716
717 // Constructs high order orthonormal tangent space T and inverse of metric tensor
718 Kokkos::single(Kokkos::PerTeam(teamMember), [&] () {
719
720 double grad_xi[2] = {grad_xi1, grad_xi2};
721 double T_row[3];
722
723 // Construct T (high order approximation of orthonormal tangent vectors)
724 for (int i=0; i<dimensions-1; ++i) {
725 for (int j=0; j<dimensions; ++j) {
726 T_row[j] = T(i,j);
727 }
728 // build
729 for (int j=0; j<dimensions; ++j) {
730 T(i,j) = grad_xi[i]*T(dimensions-1,j);
731 T(i,j) += T_row[j];
732 }
733 }
734
735 // calculate norm
736 double norm = 0;
737 for (int j=0; j<dimensions; ++j) {
738 norm += T(0,j)*T(0,j);
739 }
740
741 // normalize first vector
742 norm = std::sqrt(norm);
743 for (int j=0; j<dimensions; ++j) {
744 T(0,j) /= norm;
745 }
746
747 // orthonormalize next vector
748 if (dimensions-1 == 2) { // 2d manifold
749 double dot_product = T(0,0)*T(1,0) + T(0,1)*T(1,1) + T(0,2)*T(1,2);
750 for (int j=0; j<dimensions; ++j) {
751 T(1,j) -= dot_product*T(0,j);
752 }
753 // normalize second vector
754 norm = 0;
755 for (int j=0; j<dimensions; ++j) {
756 norm += T(1,j)*T(1,j);
757 }
758 norm = std::sqrt(norm);
759 for (int j=0; j<dimensions; ++j) {
760 T(1,j) /= norm;
761 }
762 }
763
764 // get normal vector to first two rows of T
765 double norm_t_normal = 0;
766 if (dimensions>2) {
767 T(dimensions-1,0) = T(0,1)*T(1,2) - T(1,1)*T(0,2);
768 norm_t_normal += T(dimensions-1,0)*T(dimensions-1,0);
769 T(dimensions-1,1) = -(T(0,0)*T(1,2) - T(1,0)*T(0,2));
770 norm_t_normal += T(dimensions-1,1)*T(dimensions-1,1);
771 T(dimensions-1,2) = T(0,0)*T(1,1) - T(1,0)*T(0,1);
772 norm_t_normal += T(dimensions-1,2)*T(dimensions-1,2);
773 } else {
774 T(dimensions-1,0) = T(1,1) - T(0,1);
775 norm_t_normal += T(dimensions-1,0)*T(dimensions-1,0);
776 T(dimensions-1,1) = T(0,0) - T(1,0);
777 norm_t_normal += T(dimensions-1,1)*T(dimensions-1,1);
778 }
779 norm_t_normal = std::sqrt(norm_t_normal);
780 for (int i=0; i<dimensions-1; ++i) {
781 T(dimensions-1,i) /= norm_t_normal;
782 }
783 });
784 teamMember.team_barrier();
785 }
786};
787
788//! Functor to determine if tangent directions need reordered, and to reorder them if needed
789//! We require that the normal is consistent with a right hand rule on the tangent vectors
791
793
795
796 KOKKOS_INLINE_FUNCTION
797 void operator()(const member_type& teamMember) const {
798
799 /*
800 * Dimensions
801 */
802
803 const int target_index = _data._initial_index_for_batch + teamMember.league_rank();
804 auto dimensions = _data._dimensions;
805
806 /*
807 * Data
808 */
809
810 scratch_matrix_right_type T(_data.T_data + target_index*dimensions*dimensions, dimensions, dimensions);
811 scratch_vector_type N(_data.ref_N_data + target_index*dimensions, dimensions);
812
813 // take the dot product of the calculated normal in the tangent bundle with a given reference outward normal
814 // direction provided by the user. if the dot product is negative, flip the tangent vector ordering
815 // and flip the sign on the normal direction.
816 Kokkos::single(Kokkos::PerTeam(teamMember), [&] () {
817 compadre_kernel_assert_debug(dimensions > 1
818 && "FixTangentDirectionOrder called on manifold with a dimension of 0.");
819 double dot_product = 0;
820 for (int i=0; i<dimensions; ++i) {
821 dot_product += T(dimensions-1,i) * N(i);
822
823 }
824 if (dot_product < 0) {
825 if (dimensions==3) {
826 for (int i=0; i<dimensions; ++i) {
827 // swap tangent directions
828 double tmp = T(0,i);
829 T(0,i) = T(1,i);
830 T(1,i) = tmp;
831 }
832 }
833 for (int i=0; i<dimensions; ++i) {
834 // flip the sign of the normal vector
835 T(dimensions-1,i) *= -1;
836
837 }
838 }
839 });
840 teamMember.team_barrier();
841 }
842};
843
844//! Functor to evaluate curvature targets and apply to coefficients of curvature reconstruction
846
848
850
851 KOKKOS_INLINE_FUNCTION
852 void operator()(const member_type& teamMember) const {
853
854 /*
855 * Dimensions
856 */
857
858 const int target_index = _data._initial_index_for_batch + teamMember.league_rank();
859 const int local_index = teamMember.league_rank();
860 auto dimensions = _data._dimensions;
861
862 /*
863 * Data
864 */
865
867 + TO_GLOBAL(local_index)*TO_GLOBAL(_data.Coeffs_dim_0*_data.Coeffs_dim_1),
868 _data.Coeffs_dim_0, _data.Coeffs_dim_1);
869
871 + TO_GLOBAL(target_index)*TO_GLOBAL(dimensions*dimensions),
872 dimensions, dimensions);
873
875 + target_index*TO_GLOBAL(_data.manifold_NP), _data.manifold_NP);
876
878 + TO_GLOBAL(local_index)*TO_GLOBAL(_data.P_target_row_dim_0*_data.P_target_row_dim_1),
880
881
882 // delta, used for each thread
883 scratch_vector_type delta(teamMember.thread_scratch(_data._pm.getThreadScratchLevel(0)),
884 _data.this_num_cols);
885 scratch_vector_type thread_workspace(teamMember.thread_scratch(_data._pm.getThreadScratchLevel(0)),
887
888 /*
889 * Manifold
890 */
891
892 computeCurvatureFunctionals(_data, teamMember, delta, thread_workspace, P_target_row, &T);
893
894 Kokkos::single(Kokkos::PerTeam(teamMember), [&] () {
895 for (int j=0; j<_data.manifold_NP; ++j) { // set to zero
896 manifold_coeffs(j) = 0;
897 }
898 });
899 teamMember.team_barrier();
900 for (int i=0; i<_data._pc._nla.getNumberOfNeighborsDevice(target_index); ++i) {
901 XYZ rel_coord = _data._pc.getRelativeCoord(target_index, i, dimensions, &T);
902 double normal_coordinate = rel_coord[dimensions-1];
903
904 Kokkos::single(Kokkos::PerTeam(teamMember), [&] () {
905 Kokkos::single(Kokkos::PerTeam(teamMember), [&] () {
906 // coefficients without a target premultiplied
907 for (int j=0; j<_data.manifold_NP; ++j) {
908 manifold_coeffs(j) += Q(j,i) * normal_coordinate;
909 }
910 });
911 });
912 teamMember.team_barrier();
913 }
914 }
915};
916
917//! Functor to assemble the P*sqrt(weights) matrix and construct sqrt(weights)*Identity
919
921
923
924 KOKKOS_INLINE_FUNCTION
925 void operator()(const member_type& teamMember) const {
926
927 /*
928 * Dimensions
929 */
930
931 const int target_index = _data._initial_index_for_batch + teamMember.league_rank();
932 const int local_index = teamMember.league_rank();
933 auto dimensions = _data._dimensions;
934 const int this_num_rows = _data._sampling_multiplier*_data._pc._nla.getNumberOfNeighborsDevice(target_index);
935
936 /*
937 * Data
938 */
939
941 + TO_GLOBAL(local_index)*TO_GLOBAL(_data.P_dim_0)*TO_GLOBAL(_data.P_dim_1),
942 _data.P_dim_0, _data.P_dim_1);
944 + TO_GLOBAL(local_index)*TO_GLOBAL(_data.RHS_dim_0)*TO_GLOBAL(_data.RHS_dim_1),
945 _data.RHS_dim_0, _data.RHS_dim_1);
947 + TO_GLOBAL(local_index)*TO_GLOBAL(_data.max_num_rows), _data.max_num_rows);
949 + TO_GLOBAL(target_index)*TO_GLOBAL(dimensions)*TO_GLOBAL(dimensions), dimensions, dimensions);
950
951 // delta, used for each thread
952 scratch_vector_type delta(teamMember.thread_scratch(_data._pm.getThreadScratchLevel(0)),
953 _data.this_num_cols);
954 scratch_vector_type thread_workspace(teamMember.thread_scratch(_data._pm.getThreadScratchLevel(0)),
956
957 /*
958 * Manifold
959 */
960
961 createWeightsAndP(_data, teamMember, delta, thread_workspace, PsqrtW, w, dimensions-1,
962 _data._poly_order, true /* weight with W*/, &T, _data._reconstruction_space,
964
965 if (_data._dense_solver_type != DenseSolverType::LU) {
966 // fill in RHS with Identity * sqrt(weights)
967 double * Q_data = Q.data();
968 Kokkos::parallel_for(Kokkos::TeamVectorRange(teamMember,this_num_rows), [&] (const int i) {
969 Q_data[i] = std::sqrt(w(i));
970 });
971 } else {
972 // create global memory for matrix M = PsqrtW^T*PsqrtW
973 // don't need to cast into scratch_matrix_left_type since the matrix is symmetric
975 + TO_GLOBAL(local_index)*TO_GLOBAL(_data.RHS_dim_0*_data.RHS_dim_1),
976 _data.RHS_dim_0, _data.RHS_dim_1);
977
978 // Assemble matrix M
979 KokkosBatched::TeamVectorGemm<member_type,KokkosBatched::Trans::Transpose,KokkosBatched::Trans::NoTranspose,KokkosBatched::Algo::Gemm::Unblocked>
980 ::invoke(teamMember,
981 1.0,
982 PsqrtW,
983 PsqrtW,
984 0.0,
985 M);
986 teamMember.team_barrier();
987
988 // Multiply PsqrtW with sqrt(W) to get PW
989 Kokkos::parallel_for(Kokkos::TeamThreadRange(teamMember, _data.max_num_rows), [&] (const int i) {
990 Kokkos::parallel_for(Kokkos::ThreadVectorRange(teamMember, _data.this_num_cols), [&] (const int j) {
991 PsqrtW(i, j) = PsqrtW(i, j)*std::sqrt(w(i));
992 });
993 });
994 }
995 teamMember.team_barrier();
996 }
997};
998
999//! Functor to evaluate targets on a manifold
1001
1003
1005
1006 KOKKOS_INLINE_FUNCTION
1007 void operator()(const member_type& teamMember) const {
1008
1009 /*
1010 * Dimensions
1011 */
1012
1013 const int target_index = _data._initial_index_for_batch + teamMember.league_rank();
1014 const int local_index = teamMember.league_rank();
1015 auto dimensions = _data._dimensions;
1016
1017 /*
1018 * Data
1019 */
1020
1022 + TO_GLOBAL(target_index)*TO_GLOBAL(dimensions*dimensions), dimensions, dimensions);
1024 + target_index*TO_GLOBAL(_data.manifold_NP), _data.manifold_NP);
1026 + TO_GLOBAL(local_index)*TO_GLOBAL(_data.P_target_row_dim_0*_data.P_target_row_dim_1),
1028
1029 // delta, used for each thread
1030 scratch_vector_type delta(teamMember.thread_scratch(_data._pm.getThreadScratchLevel(0)),
1031 _data.this_num_cols);
1032 scratch_vector_type thread_workspace(teamMember.thread_scratch(_data._pm.getThreadScratchLevel(0)),
1033 _data.thread_workspace_dim);
1034
1035 /*
1036 * Apply Standard Target Evaluations to Polynomial Coefficients
1037 */
1038
1039 computeTargetFunctionalsOnManifold(_data, teamMember, delta, thread_workspace, P_target_row, T, manifold_coeffs);
1040
1041 }
1042};
1043
1044///@}
1045
1046} // Compadre
1047
1048#endif
#define compadre_kernel_assert_debug(condition)
#define TO_GLOBAL(variable)
KOKKOS_INLINE_FUNCTION int getThreadScratchLevel(const int level) const
KOKKOS_INLINE_FUNCTION int getTeamScratchLevel(const int level) const
Kokkos::Random_XorShift64_Pool pool_type
KOKKOS_INLINE_FUNCTION void computeCurvatureFunctionals(const TargetData &data, const member_type &teamMember, scratch_vector_type delta, scratch_vector_type thread_workspace, scratch_matrix_right_type P_target_row, const scratch_matrix_right_type *V, const local_index_type local_neighbor_index=-1)
Evaluates a polynomial basis for the curvature with a gradient target functional applied.
Kokkos::View< double *, Kokkos::MemoryTraits< Kokkos::Unmanaged > > scratch_vector_type
constexpr SamplingFunctional VaryingManifoldVectorPointSample
For integrating polynomial dotted with normal over an edge.
KOKKOS_INLINE_FUNCTION void computeTargetFunctionals(const TargetData &data, const member_type &teamMember, scratch_vector_type delta, scratch_vector_type thread_workspace, scratch_matrix_right_type P_target_row)
Evaluates a polynomial basis with a target functional applied to each member of the basis.
ConstraintType
Constraint type.
constexpr SamplingFunctional StaggeredEdgeIntegralSample
Samples consist of the result of integrals of a vector dotted with the tangent along edges between ne...
KOKKOS_INLINE_FUNCTION void createWeightsAndP(const BasisData &data, const member_type &teamMember, scratch_vector_type delta, scratch_vector_type thread_workspace, scratch_matrix_right_type P, scratch_vector_type w, const int dimension, int polynomial_order, bool weight_p=false, scratch_matrix_right_type *V=NULL, const ReconstructionSpace reconstruction_space=ReconstructionSpace::ScalarTaylorPolynomial, const SamplingFunctional polynomial_sampling_functional=PointSample)
Fills the _P matrix with either P or P*sqrt(w)
constexpr SamplingFunctional PointSample
Available sampling functionals.
team_policy::member_type member_type
ReconstructionSpace
Space in which to reconstruct polynomial.
@ ScalarTaylorPolynomial
Scalar polynomial basis centered at the target site and scaled by sum of basis powers e....
constexpr SamplingFunctional StaggeredEdgeAnalyticGradientIntegralSample
Analytical integral of a gradient source vector is just a difference of the scalar source at neighbor...
KOKKOS_INLINE_FUNCTION void createWeightsAndPForCurvature(const BasisData &data, const member_type &teamMember, scratch_vector_type delta, scratch_vector_type thread_workspace, scratch_matrix_right_type P, scratch_vector_type w, const int dimension, bool only_specific_order, scratch_matrix_right_type *V=NULL)
Fills the _P matrix with P*sqrt(w) for use in solving for curvature.
ProblemType
Problem type, that optionally can handle manifolds.
@ MANIFOLD
Solve GMLS problem on a manifold (will use QR or SVD to solve the resultant GMLS problem dependent on...
constexpr SamplingFunctional ManifoldVectorPointSample
Point evaluations of the entire vector source function (but on a manifold, so it includes a transform...
DenseSolverType
Dense solver type.
KOKKOS_INLINE_FUNCTION void applyTargetsToCoefficients(const SolutionData &data, const member_type &teamMember, scratch_matrix_right_type Q, scratch_matrix_right_type P_target_row)
For applying the evaluations from a target functional to the polynomial coefficients.
KOKKOS_INLINE_FUNCTION void evaluateConstraints(scratch_matrix_right_type M, scratch_matrix_right_type PsqrtW, const ConstraintType constraint_type, const ReconstructionSpace reconstruction_space, const int NP, const double cutoff_p, const int dimension, const int num_neighbors=0, scratch_matrix_right_type *T=NULL)
KOKKOS_INLINE_FUNCTION void calcGradientPij(const BasisData &data, const member_type &teamMember, double *delta, double *thread_workspace, const int target_index, int neighbor_index, const double alpha, const int partial_direction, const int dimension, const int poly_order, bool specific_order_only, const scratch_matrix_right_type *V, const ReconstructionSpace reconstruction_space, const SamplingFunctional polynomial_sampling_functional, const int evaluation_site_local_index=0)
Evaluates the gradient of a polynomial basis under the Dirac Delta (pointwise) sampling function.
KOKKOS_INLINE_FUNCTION void computeTargetFunctionalsOnManifold(const TargetData &data, const member_type &teamMember, scratch_vector_type delta, scratch_vector_type thread_workspace, scratch_matrix_right_type P_target_row, scratch_matrix_right_type V, scratch_vector_type curvature_coefficients)
Evaluates a polynomial basis with a target functional applied, using information from the manifold cu...
WeightingFunctionType
Available weighting kernel function types.
Kokkos::View< double **, layout_right, Kokkos::MemoryTraits< Kokkos::Unmanaged > > scratch_matrix_right_type
Functor to evaluate curvature targets and apply to coefficients of curvature reconstruction.
KOKKOS_INLINE_FUNCTION void operator()(const member_type &teamMember) const
Functor to apply target evaluation to polynomial coefficients to store in _alphas.
KOKKOS_INLINE_FUNCTION void operator()(const member_type &teamMember) const
ApplyTargets(GMLSSolutionData data)
Functor to assemble the P*sqrt(weights) matrix and construct sqrt(weights)*Identity for curvature.
KOKKOS_INLINE_FUNCTION void operator()(const member_type &teamMember) const
Functor to assemble the P*sqrt(weights) matrix and construct sqrt(weights)*Identity.
KOKKOS_INLINE_FUNCTION void operator()(const member_type &teamMember) const
Functor to assemble the P*sqrt(weights) matrix and construct sqrt(weights)*Identity.
KOKKOS_INLINE_FUNCTION void operator()(const member_type &teamMember) const
Functor to create a coarse tangent approximation from a given neighborhood of points.
KOKKOS_INLINE_FUNCTION void operator()(const member_type &teamMember) const
Functor to calculate prestencil weights to apply to data to transform into a format expected by a GML...
KOKKOS_INLINE_FUNCTION void operator()(const member_type &teamMember) const
Functor to evaluate targets on a manifold.
KOKKOS_INLINE_FUNCTION void operator()(const member_type &teamMember) const
Functor to evaluate targets operations on the basis.
KOKKOS_INLINE_FUNCTION void operator()(const member_type &teamMember) const
Functor to determine if tangent directions need reordered, and to reorder them if needed We require t...
KOKKOS_INLINE_FUNCTION void operator()(const member_type &teamMember) const
GMLS::point_connections_type _pc
SamplingFunctional _polynomial_sampling_functional
ReconstructionSpace _reconstruction_space
DenseSolverType _dense_solver_type
SolutionSet< device_memory_space > _d_ss
Kokkos::View< double **, layout_right > _source_extra_data
WeightingFunctionType _curvature_weighting_type
GMLS::point_connections_type _additional_pc
Kokkos::View< double *****, layout_right > _prestencil_weights
Kokkos::View< TargetOperation * > _curvature_support_operations
Kokkos::View< double * > _epsilons
Kokkos::View< TargetOperation * > _operations
WeightingFunctionType _weighting_type
Kokkos::View< double **, layout_right > _target_extra_data
SamplingFunctional _data_sampling_functional
Kokkos::View< int * > additional_number_of_neighbors_list
Kokkos::View< int * > number_of_neighbors_list
SolutionSet< device_memory_space > _d_ss
Functor to evaluate curvature targets and construct accurate tangent direction approximation for mani...
KOKKOS_INLINE_FUNCTION void operator()(const member_type &teamMember) const
KOKKOS_INLINE_FUNCTION int getNumberOfNeighborsDevice(int target_index) const
Get number of neighbors for a given target (device)
KOKKOS_INLINE_FUNCTION XYZ getRelativeCoord(const int target_index, const int neighbor_list_num, const int dimension, const scratch_matrix_right_type *V=NULL) const
Returns the relative coordinate as a vector between the target site and the neighbor site.
All vairables and functionality related to the layout and storage of GMLS solutions (alpha values)
KOKKOS_INLINE_FUNCTION int getTargetOffsetIndex(const int lro_num, const int input_component, const int output_component, const int evaluation_site_local_index=0) const
Handles offset from operation input/output + extra evaluation sites.