Compadre 1.6.4
Loading...
Searching...
No Matches
Compadre_GMLS.cpp
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#include "Compadre_GMLS.hpp"
10#include "Compadre_Functors.hpp"
11
12namespace Compadre {
13
14void GMLS::generatePolynomialCoefficients(const int number_of_batches, const bool keep_coefficients, const bool clear_cache) {
15
16 compadre_assert_release( (keep_coefficients==false || number_of_batches==1)
17 && "keep_coefficients is set to true, but number of batches exceeds 1.");
18
19 /*
20 * Verify PointConnections are valid
21 */
22 this->verifyPointConnections(true);
24
25 // ensure that solution set has neighbor list consistent with point connections
28
29 /*
30 * Generate Quadrature
31 */
33
34 /*
35 * Generate SolutionSet on device
36 */
38
39 /*
40 * Operations to Device
41 */
42
43 // copy over operations
44 _operations = decltype(_operations)("operations", _h_ss._lro.size());
45 _host_operations = Kokkos::create_mirror_view(_operations);
46
47 // make sure at least one target operation specified
49 && "No target operations added to GMLS class before calling generatePolynomialCoefficients().");
50
51 // loop through list of linear reconstruction operations to be performed and set them on the host
52 for (size_t i=0; i<_h_ss._lro.size(); ++i) _host_operations(i) = _h_ss._lro[i];
53
54 // get copy of operations on the device
55 Kokkos::deep_copy(_operations, _host_operations);
56
57 /*
58 * Initialize Alphas and Prestencil Weights
59 */
60
61 // throw an assertion for QR solver incompatibility
62 // TODO: this is a temporary location for this check, in the future the
63 // constraint type could be an object that can check when given a dense_solver_type
66 && "Cannot solve GMLS problems with the NEUMANN_GRAD_SCALAR constraint using QR Factorization.");
67
68 // calculate the additional size for different constraint problems
71
72 // initialize all alpha values to be used for taking the dot product with data to get a reconstruction
73 try {
75 int total_added_alphas = _pc._target_coordinates.extent(0)*_d_ss._added_alpha_size;
76 _d_ss._alphas =
77 decltype(_d_ss._alphas)("alphas", (total_neighbors + TO_GLOBAL(total_added_alphas))
79 // this deep copy writes to all theoretically allocated memory,
80 // ensuring that allocation attempted was successful
81 Kokkos::deep_copy(_d_ss._alphas, 0.0);
82 } catch(std::exception &e) {
83 printf("Insufficient memory to store alphas: \n\n%s", e.what());
84 throw e;
85 }
86
87 // initialize the prestencil weights that are applied to sampling data to put it into a form
88 // that the GMLS operator will be able to operate on
90 int max_num_neighbors = _pc._nla.getMaxNumNeighbors();
91 try {
92 _prestencil_weights = decltype(_prestencil_weights)("Prestencil weights",
93 std::pow(2,sro.use_target_site_weights),
94 (sro.transform_type==DifferentEachTarget
95 || sro.transform_type==DifferentEachNeighbor) ?
97 (sro.transform_type==DifferentEachNeighbor) ?
98 max_num_neighbors : 1,
99 (sro.output_rank>0) ?
101 (sro.input_rank>0) ?
103 } catch(std::exception &e) {
104 printf("Insufficient memory to store prestencil weights: \n\n%s", e.what());
105 throw e;
106 }
107 Kokkos::fence();
108
109 /*
110 * Determine if Nonstandard Sampling Dimension or Basis Component Dimension
111 */
112
113 // calculate the dimension of the basis (a vector space on a manifold requires two components, for example)
115
116 // calculate sampling dimension
118
119 // effective number of components in the basis
121
122 // special case for using a higher order for sampling from a polynomial space that are gradients of a scalar polynomial
124 // if the reconstruction is being made with a gradient of a basis, then we want that basis to be one order higher so that
125 // the gradient is consistent with the convergence order expected.
126 _poly_order += 1;
128 }
129
130 /*
131 * Dimensions
132 */
133
134 // for tallying scratch space needed for device kernel calls
135 int team_scratch_size_a = 0;
136
137 // TEMPORARY, take to zero after conversion
138 int team_scratch_size_b = 0;
139 int thread_scratch_size_a = 0;
140 int thread_scratch_size_b = 0;
141
142 // dimensions that are relevant for each subproblem
143 int max_num_rows = _sampling_multiplier*max_num_neighbors;
144 int this_num_cols = _basis_multiplier*_NP;
145 int manifold_NP = 0;
146
147 // determines whether RHS will be stored implicitly
148 // if true, instead of RHS storing the full sqrt(W) explicitly,
149 // only the diagonal entries of sqrt(W) will be stored as a
150 // 1D array beginning at entry with matrix coordinate (0,0)
151 bool implicit_RHS = (_dense_solver_type != DenseSolverType::LU);
152
153 int basis_powers_space_multiplier = (_reconstruction_space == BernsteinPolynomial) ? 2 : 1;
155 // these dimensions already calculated differ in the case of manifolds
158 const int max_manifold_NP = (manifold_NP > _NP) ? manifold_NP : _NP;
159 this_num_cols = _basis_multiplier*max_manifold_NP;
160 const int max_poly_order = (_poly_order > _curvature_poly_order) ? _poly_order : _curvature_poly_order;
161
162 /*
163 * Calculate Scratch Space Allocations
164 */
165
166 team_scratch_size_b += scratch_matrix_right_type::shmem_size(_dimensions, _dimensions); // PTP matrix
167 team_scratch_size_b += scratch_vector_type::shmem_size( (_dimensions-1)*max_num_neighbors ); // manifold_gradient
168
169 thread_scratch_size_a += scratch_vector_type::shmem_size(this_num_cols); // delta, used for each thread
170 thread_scratch_size_a += scratch_vector_type::shmem_size(
171 (max_poly_order+1)*_global_dimensions*basis_powers_space_multiplier); // temporary space for powers in basis
173 team_scratch_size_b += scratch_vector_type::shmem_size(max_num_neighbors); // t1 work vector for prestencils
174 team_scratch_size_b += scratch_vector_type::shmem_size(max_num_neighbors); // t2 work vector for prestencils
175 thread_scratch_size_b += scratch_vector_type::shmem_size(_dimensions*_dimensions); // temporary tangent calculations, used for each thread
176 }
177
178 // allocate data on the device (initialized to zero)
180 _T = Kokkos::View<double*>("tangent approximation",_pc._target_coordinates.extent(0)*_dimensions*_dimensions);
181 Kokkos::deep_copy(_T, 0.0);
182 } else {
184 "Provided tangent_directions has number of targets different than target_coordinates");
185 }
186 _manifold_curvature_coefficients = Kokkos::View<double*>("manifold curvature coefficients",
187 _pc._target_coordinates.extent(0)*manifold_NP);
188 Kokkos::deep_copy(_manifold_curvature_coefficients, 0.0);
189
190 } else { // Standard GMLS
191
192 /*
193 * Calculate Scratch Space Allocations
194 */
195
196 thread_scratch_size_a += scratch_vector_type::shmem_size(this_num_cols); // delta, used for each thread
197 thread_scratch_size_a += scratch_vector_type::shmem_size(
198 (_poly_order+1)*_global_dimensions*basis_powers_space_multiplier); // temporary space for powers in basis
199
200 }
201 _pm.setTeamScratchSize(0, team_scratch_size_a);
202 _pm.setTeamScratchSize(1, team_scratch_size_b);
203 _pm.setThreadScratchSize(0, thread_scratch_size_a);
204 _pm.setThreadScratchSize(1, thread_scratch_size_b);
205
206 /*
207 * Calculate the size for matrix P and RHS
208 */
209
210 int RHS_dim_0, RHS_dim_1;
211 getRHSDims(_dense_solver_type, _constraint_type, _reconstruction_space, _dimensions, max_num_rows, this_num_cols, RHS_dim_0, RHS_dim_1);
212
213 int P_dim_0, P_dim_1;
214 getPDims(_dense_solver_type, _constraint_type, _reconstruction_space, _dimensions, max_num_rows, this_num_cols, P_dim_0, P_dim_1);
215
216 /*
217 * Allocate Global Device Storage of Data Needed Over Multiple Calls
218 */
219
220 global_index_type max_batch_size = (_pc._target_coordinates.extent(0) + TO_GLOBAL(number_of_batches) - 1) / TO_GLOBAL(number_of_batches);
221 try {
222 _RHS = Kokkos::View<double*>("RHS", max_batch_size*TO_GLOBAL(RHS_dim_0)*TO_GLOBAL(RHS_dim_1));
223 _P = Kokkos::View<double*>("P", max_batch_size*TO_GLOBAL(P_dim_0)*TO_GLOBAL(P_dim_1));
224 _w = Kokkos::View<double*>("w", max_batch_size*TO_GLOBAL(max_num_rows));
225 _Z = Kokkos::View<double*>("Z", max_batch_size*TO_GLOBAL(_d_ss._total_alpha_values*_d_ss._max_evaluation_sites_per_target*this_num_cols));
226
227 } catch (std::exception &e) {
228 printf("Failed to allocate space for RHS, P, and w. Consider increasing number_of_batches: \n\n%s", e.what());
229 throw e;
230 }
231 Kokkos::fence();
232
233 /*
234 * Calculate Optimal Threads Based On Levels of Parallelism
235 */
236
239 && "Normal vectors are required for solving GMLS problems with the NEUMANN_GRAD_SCALAR constraint.");
240 }
241
243 for (int batch_num=0; batch_num<number_of_batches; ++batch_num) {
244
245 auto this_batch_size = std::min(_pc._target_coordinates.extent(0)-_initial_index_for_batch, max_batch_size);
246 Kokkos::deep_copy(_RHS, 0.0);
247 Kokkos::deep_copy(_P, 0.0);
248 Kokkos::deep_copy(_w, 0.0);
249 Kokkos::deep_copy(_Z, 0.0);
250
251 auto gmls_basis_data = this->extractBasisData();
252 auto gmls_solution_data = this->extractSolutionData();
253
254
255 // even kernels that should run on other # of vector lanes do not (on GPU)
256 auto tp = _pm.TeamPolicyThreadsAndVectors(this_batch_size, _pm._default_threads, 1);
257 //auto tp = _pm.TeamPolicyThreadsAndVectors(this_batch_size, _pm._default_threads, _pm._default_vector_lanes);
258 //const auto work_item_property = Kokkos::Experimental::WorkItemProperty::HintLightWeight;
259 //const auto tp2 = Kokkos::Experimental::require(tp, work_item_property);
260
262
263 /*
264 * MANIFOLD Problems
265 */
266
267 //auto functor_name = Name(gmls_basis_data);
268 //Kokkos::parallel_for("Name", tp, functor_name);
269 if (!_orthonormal_tangent_space_provided) { // user did not specify orthonormal tangent directions, so we approximate them first
270 // coarse tangent plane approximation construction of P^T*P
271 auto functor_compute_coarse_tangent_plane = ComputeCoarseTangentPlane(gmls_basis_data);
272 Kokkos::parallel_for("ComputeCoarseTangentPlane", tp, functor_compute_coarse_tangent_plane);
273
274 // if the user provided the reference outward normal direction, then orient the computed or user provided
275 // outward normal directions in the tangent bundle
277 // use the reference outward normal direction provided by the user to orient
278 // the tangent bundle
279 auto functor_fix_tangent_direction_ordering = FixTangentDirectionOrdering(gmls_basis_data);
280 Kokkos::parallel_for("FixTangentDirectionOrdering", tp, functor_fix_tangent_direction_ordering);
281 }
282
283 // assembles the P*sqrt(weights) matrix and constructs sqrt(weights)*Identity for curvature
284 auto functor_assemble_curvature_psqrtw = AssembleCurvaturePsqrtW(gmls_basis_data);
285 Kokkos::parallel_for("AssembleCurvaturePsqrtW", tp, functor_assemble_curvature_psqrtw);
286
288 // solves P^T*P against P^T*W with LU, stored in P
289 Kokkos::Profiling::pushRegion("Curvature LU Factorization");
290 // batchLU expects layout_left matrix tiles for B
291 // by giving it layout_right matrix tiles with reverse ordered ldb and ndb
292 // it effects a transpose of _P in layout_left
293 GMLS_LinearAlgebra::batchQRPivotingSolve<layout_right,layout_left,layout_right>(_pm, _RHS.data(), RHS_dim_0, RHS_dim_1, _P.data(), P_dim_1, P_dim_0, manifold_NP, manifold_NP, max_num_neighbors, this_batch_size, implicit_RHS);
294 Kokkos::Profiling::popRegion();
295 } else {
296 // solves P*sqrt(weights) against sqrt(weights)*Identity with QR, stored in RHS
297 Kokkos::Profiling::pushRegion("Curvature QR+Pivoting Factorization");
298 GMLS_LinearAlgebra::batchQRPivotingSolve<layout_right,layout_right,layout_right>(_pm, _P.data(), P_dim_0, P_dim_1, _RHS.data(), RHS_dim_0, RHS_dim_1, max_num_neighbors, manifold_NP, max_num_neighbors, this_batch_size, implicit_RHS);
299 Kokkos::Profiling::popRegion();
300 }
301
302 // evaluates targets, applies target evaluation to polynomial coefficients for curvature
303 auto functor_get_accurate_tangent_directions = GetAccurateTangentDirections(gmls_basis_data);
304 Kokkos::parallel_for("GetAccurateTangentDirections", tp, functor_get_accurate_tangent_directions);
305
306 // Due to converting layout, entries that are assumed zeros may become non-zeros.
307 Kokkos::deep_copy(_P, 0.0);
308
309 if (batch_num==number_of_batches-1) {
310 // copy tangent bundle from device back to host
311 _host_T = Kokkos::create_mirror_view(_T);
312 Kokkos::deep_copy(_host_T, _T);
313 }
314 }
315
316 // this time assembling curvature PsqrtW matrix is using a highly accurate approximation of the tangent, previously calculated
317 // assembles the P*sqrt(weights) matrix and constructs sqrt(weights)*Identity for curvature
318 auto functor_assemble_curvature_psqrtw = AssembleCurvaturePsqrtW(gmls_basis_data);
319 Kokkos::parallel_for("AssembleCurvaturePsqrtW", tp, functor_assemble_curvature_psqrtw);
320
322 // solves P^T*P against P^T*W with LU, stored in P
323 Kokkos::Profiling::pushRegion("Curvature LU Factorization");
324 GMLS_LinearAlgebra::batchQRPivotingSolve<layout_right,layout_left,layout_right>(_pm, _RHS.data(), RHS_dim_0, RHS_dim_1, _P.data(), P_dim_1, P_dim_0, manifold_NP, manifold_NP, max_num_neighbors, this_batch_size, implicit_RHS);
325 Kokkos::Profiling::popRegion();
326 } else {
327 // solves P*sqrt(weights) against sqrt(weights)*Identity, stored in RHS
328 Kokkos::Profiling::pushRegion("Curvature QR+Pivoting Factorization");
329 GMLS_LinearAlgebra::batchQRPivotingSolve<layout_right,layout_right,layout_right>(_pm, _P.data(), P_dim_0, P_dim_1, _RHS.data(), RHS_dim_0, RHS_dim_1, max_num_neighbors, manifold_NP, max_num_neighbors, this_batch_size, implicit_RHS);
330 Kokkos::Profiling::popRegion();
331 }
332
333 // evaluates targets, applies target evaluation to polynomial coefficients for curvature
334 auto functor_apply_curvature_targets = ApplyCurvatureTargets(gmls_basis_data);
335 Kokkos::parallel_for("ApplyCurvatureTargets", tp, functor_apply_curvature_targets);
336 Kokkos::fence();
337
338 // prestencil weights calculated here. appropriate because:
339 // precedes polynomial reconstruction from data (replaces contents of _RHS)
340 // follows reconstruction of geometry
341 // calculate prestencil weights
342 auto functor_compute_prestencil_weights = ComputePrestencilWeights(gmls_basis_data);
343 Kokkos::parallel_for("ComputePrestencilWeights", tp, functor_compute_prestencil_weights);
344
345 // Due to converting layout, entried that are assumed zeros may become non-zeros.
346 Kokkos::deep_copy(_P, 0.0);
347
348 // assembles the P*sqrt(weights) matrix and constructs sqrt(weights)*Identity
349 auto functor_assemble_manifold_psqrtw = AssembleManifoldPsqrtW(gmls_basis_data);
350 Kokkos::parallel_for("AssembleManifoldPsqrtW", tp, functor_assemble_manifold_psqrtw);
351
352 // solves P*sqrt(weights) against sqrt(weights)*Identity, stored in RHS
354 Kokkos::Profiling::pushRegion("Manifold LU Factorization");
355 GMLS_LinearAlgebra::batchQRPivotingSolve<layout_right,layout_left,layout_right>(_pm, _RHS.data(), RHS_dim_0, RHS_dim_1, _P.data(), P_dim_1, P_dim_0, this_num_cols, this_num_cols, max_num_rows, this_batch_size, implicit_RHS);
356 Kokkos::Profiling::popRegion();
357 } else {
358 Kokkos::Profiling::pushRegion("Manifold QR+Pivoting Factorization");
359 GMLS_LinearAlgebra::batchQRPivotingSolve<layout_right,layout_right,layout_right>(_pm, _P.data(), P_dim_0, P_dim_1, _RHS.data(), RHS_dim_0, RHS_dim_1, max_num_rows, this_num_cols, max_num_rows, this_batch_size, implicit_RHS);
360 Kokkos::Profiling::popRegion();
361 }
362 Kokkos::fence();
363
364 } else {
365
366 /*
367 * STANDARD GMLS Problems
368 */
369
370 // assembles the P*sqrt(weights) matrix and constructs sqrt(weights)*Identity
371 auto functor_assemble_standard_psqrtw = AssembleStandardPsqrtW(gmls_basis_data);
372 //printf("size of assemble: %lu\n", sizeof(functor_assemble_standard_psqrtw));
373 Kokkos::parallel_for("AssembleStandardPsqrtW", tp, functor_assemble_standard_psqrtw);
374 Kokkos::fence();
375
376 // solves P*sqrt(weights) against sqrt(weights)*Identity, stored in RHS
378 Kokkos::Profiling::pushRegion("LU Factorization");
379 GMLS_LinearAlgebra::batchQRPivotingSolve<layout_right,layout_left,layout_right>(_pm, _RHS.data(), RHS_dim_0, RHS_dim_1, _P.data(), P_dim_1, P_dim_0, this_num_cols + added_coeff_size, this_num_cols + added_coeff_size, max_num_rows + _d_ss._added_alpha_size, this_batch_size, implicit_RHS);
380 Kokkos::Profiling::popRegion();
381 } else {
382 Kokkos::Profiling::pushRegion("QR+Pivoting Factorization");
384 GMLS_LinearAlgebra::batchQRPivotingSolve<layout_right,layout_right,layout_right>(_pm, _RHS.data(), RHS_dim_0, RHS_dim_1, _P.data(), P_dim_1, P_dim_0, this_num_cols + added_coeff_size, this_num_cols + added_coeff_size, max_num_rows + _d_ss._added_alpha_size, this_batch_size, implicit_RHS);
385 } else {
386 GMLS_LinearAlgebra::batchQRPivotingSolve<layout_right,layout_right,layout_right>(_pm, _P.data(), P_dim_0, P_dim_1, _RHS.data(), RHS_dim_0, RHS_dim_1, max_num_rows, this_num_cols, max_num_rows, this_batch_size, implicit_RHS);
387 }
388 Kokkos::Profiling::popRegion();
389 }
390
391 auto functor_compute_prestencil_weights = ComputePrestencilWeights(gmls_basis_data);
392 Kokkos::parallel_for("ComputePrestencilWeights", tp, functor_compute_prestencil_weights);
393 Kokkos::fence();
394 }
395
396 /*
397 * Calculate Optimal Threads Based On Levels of Parallelism
398 */
399
400
402
403 /*
404 * MANIFOLD Problems
405 */
406
407 // evaluates targets, applies target evaluation to polynomial coefficients to store in _alphas
408 auto functor_evaluate_manifold_targets = EvaluateManifoldTargets(gmls_basis_data);
409 Kokkos::parallel_for("EvaluateManifoldTargets", tp, functor_evaluate_manifold_targets);
410
411 } else {
412
413 /*
414 * STANDARD GMLS Problems
415 */
416
417 // evaluates targets, applies target evaluation to polynomial coefficients to store in _alphas
418 auto functor_evaluate_standard_targets = EvaluateStandardTargets(gmls_basis_data);
419 Kokkos::parallel_for("EvaluateStandardTargets", tp, functor_evaluate_standard_targets);
420 }
421
422
423 // fine grain control over applying target (most expensive part after QR solve)
426 const auto work_item_property = Kokkos::Experimental::WorkItemProperty::HintLightWeight;
427 const auto tp2 = Kokkos::Experimental::require(tp, work_item_property);
428 auto functor_apply_targets = ApplyTargets(gmls_solution_data);
429 //printf("size of apply: %lu\n", sizeof(functor_apply_targets));
430 Kokkos::parallel_for("ApplyTargets", tp2, functor_apply_targets);
431
432
433 _initial_index_for_batch += this_batch_size;
434 if ((size_t)_initial_index_for_batch == _pc._target_coordinates.extent(0)) break;
435 } // end of batch loops
436
437 if (clear_cache) {
438 // deallocate _P and _w
439 _w = Kokkos::View<double*>("w",0);
440 _Z = Kokkos::View<double*>("Z",0);
441 if (number_of_batches > 1) { // no reason to keep coefficients if they aren't all in memory
442 _RHS = Kokkos::View<double*>("RHS",0);
443 _P = Kokkos::View<double*>("P",0);
445 } else {
447 _RHS = Kokkos::View<double*>("RHS", 0);
448 if (!keep_coefficients) _P = Kokkos::View<double*>("P", 0);
449 } else {
451 _P = Kokkos::View<double*>("P", 0);
452 if (!keep_coefficients) _RHS = Kokkos::View<double*>("RHS", 0);
453 } else {
454 _RHS = Kokkos::View<double*>("RHS", 0);
455 if (!keep_coefficients) _P = Kokkos::View<double*>("P", 0);
456 }
457 }
458 if (keep_coefficients) _store_PTWP_inv_PTW = true;
459 }
460 }
461
462 /*
463 * Device to Host Copy Of Solution
464 */
465 // copy computed alphas back to the host
466 this->_d_ss._contains_valid_alphas = true;
469 _host_prestencil_weights = Kokkos::create_mirror_view(_prestencil_weights);
471 }
472
473}
474
475void GMLS::generateAlphas(const int number_of_batches, const bool keep_coefficients, const bool clear_cache) {
476
477 this->generatePolynomialCoefficients(number_of_batches, keep_coefficients, clear_cache);
478
479}
480
482 auto gmls = *this;
483 auto data = GMLSSolutionData();
484 data._sampling_multiplier = gmls._sampling_multiplier;
485 data._initial_index_for_batch = gmls._initial_index_for_batch;
486 data._d_ss = gmls._d_ss;
487
488 // store results of calculation in struct
489 const int max_num_neighbors = gmls._pc._nla.getMaxNumNeighbors();
490 const int max_num_rows = gmls._sampling_multiplier*max_num_neighbors;
491
492 // applyTargetsToCoefficients currently uses data.this_num_cols for the
493 // dot product range. Even for manifold problems, it is still appropriate to
494 // use gmls._basis_multiplier*gmls._NP. If GMLSSolutionData was ever to be
495 // used for applying solution coefficients for the curvature reconstruction,
496 // the manifolf_NP would have to be used for that application (requiring an
497 // extra argument to applyTargetsToCoefficients)
498 data.this_num_cols = gmls._basis_multiplier*gmls._NP;
499
500 int RHS_dim_0, RHS_dim_1;
501 getRHSDims(gmls._dense_solver_type, gmls._constraint_type, gmls._reconstruction_space,
502 gmls._dimensions, max_num_rows, data.this_num_cols, RHS_dim_0, RHS_dim_1);
503 int P_dim_0, P_dim_1;
504 getPDims(gmls._dense_solver_type, gmls._constraint_type, gmls._reconstruction_space,
505 gmls._dimensions, max_num_rows, data.this_num_cols, P_dim_0, P_dim_1);
506
507 if ((gmls._constraint_type == ConstraintType::NO_CONSTRAINT) && (gmls._dense_solver_type != DenseSolverType::LU)) {
508 data.Coeffs_data = gmls._RHS.data();
509 data.Coeffs_dim_0 = RHS_dim_0;
510 data.Coeffs_dim_1 = RHS_dim_1;
511 } else {
512 data.Coeffs_data = gmls._P.data();
513 data.Coeffs_dim_0 = P_dim_1;
514 data.Coeffs_dim_1 = P_dim_0;
515 }
516
517 data.P_target_row_dim_0 = gmls._d_ss._total_alpha_values*gmls._d_ss._max_evaluation_sites_per_target;
518 data.P_target_row_dim_1 = data.this_num_cols;
519 data.P_target_row_data = gmls._Z.data();
520
521 data.operations_size = gmls._operations.size();
522
523 data.number_of_neighbors_list = gmls._pc._nla._number_of_neighbors_list;
524 data.additional_number_of_neighbors_list = gmls._additional_pc._nla._number_of_neighbors_list;
525
526 return data;
527}
528
530 auto gmls = *this;
531 auto data = GMLSBasisData();
532 data._source_extra_data = gmls._source_extra_data;
533 data._target_extra_data = gmls._target_extra_data;
534 data._pc = gmls._pc;
535 data._epsilons = gmls._epsilons ;
536 data._prestencil_weights = gmls._prestencil_weights ;
537 data._additional_pc = gmls._additional_pc;
538 data._poly_order = gmls._poly_order ;
539 data._curvature_poly_order = gmls._curvature_poly_order;
540 data._NP = gmls._NP;
541 data._global_dimensions = gmls._global_dimensions;
542 data._local_dimensions = gmls._local_dimensions;
543 data._dimensions = gmls._dimensions;
544 data._reconstruction_space = gmls._reconstruction_space;
545 data._reconstruction_space_rank = gmls._reconstruction_space_rank;
546 data._dense_solver_type = gmls._dense_solver_type;
547 data._problem_type = gmls._problem_type;
548 data._constraint_type = gmls._constraint_type;
549 data._polynomial_sampling_functional = gmls._polynomial_sampling_functional;
550 data._data_sampling_functional = gmls._data_sampling_functional;
551 data._curvature_support_operations = gmls._curvature_support_operations;
552 data._operations = gmls._operations;
553 data._weighting_type = gmls._weighting_type;
554 data._curvature_weighting_type = gmls._curvature_weighting_type;
555 data._weighting_p = gmls._weighting_p;
556 data._weighting_n = gmls._weighting_n;
557 data._curvature_weighting_p = gmls._curvature_weighting_p;
558 data._curvature_weighting_n = gmls._curvature_weighting_n;
559 data._basis_multiplier = gmls._basis_multiplier;
560 data._sampling_multiplier = gmls._sampling_multiplier;
561 data._data_sampling_multiplier = gmls._data_sampling_multiplier;
562 data._initial_index_for_batch = gmls._initial_index_for_batch;
563 data._pm = gmls._pm;
564 data._order_of_quadrature_points = gmls._order_of_quadrature_points;
565 data._dimension_of_quadrature_points = gmls._dimension_of_quadrature_points;
566 data._qm = gmls._qm;
567 data._d_ss = gmls._d_ss;
568
569 data.max_num_neighbors = gmls._pc._nla.getMaxNumNeighbors();
570 data.max_num_rows = gmls._sampling_multiplier*data.max_num_neighbors;
571 int basis_powers_space_multiplier = (gmls._reconstruction_space == BernsteinPolynomial) ? 2 : 1;
572 if (gmls._problem_type == ProblemType::MANIFOLD) {
573 data.manifold_NP = GMLS::getNP(gmls._curvature_poly_order, gmls._dimensions-1,
575 data.max_manifold_NP = (data.manifold_NP > gmls._NP) ? data.manifold_NP : gmls._NP;
576 data.this_num_cols = gmls._basis_multiplier*data.max_manifold_NP;
577 data.max_poly_order = (gmls._poly_order > gmls._curvature_poly_order) ?
578 gmls._poly_order : gmls._curvature_poly_order;
579
580 data.ref_N_data = gmls._ref_N.data();
581 data.ref_N_dim = gmls._dimensions;
582
583 data.thread_workspace_dim = (data.max_poly_order+1)*gmls._global_dimensions*basis_powers_space_multiplier;
584 data.manifold_gradient_dim = (gmls._dimensions-1)*data.max_num_neighbors;
585
586 data.manifold_curvature_coefficients_data = gmls._manifold_curvature_coefficients.data();
587
588 } else {
589 data.manifold_NP = 0;
590 data.this_num_cols = gmls._basis_multiplier*gmls._NP;
591 data.thread_workspace_dim = (gmls._poly_order+1)*gmls._global_dimensions*basis_powers_space_multiplier;
592 data.manifold_gradient_dim = 0;
593 }
594
595
596 data.T_data = gmls._T.data();
597
598 data.P_target_row_dim_0 = gmls._d_ss._total_alpha_values*gmls._d_ss._max_evaluation_sites_per_target;
599 data.P_target_row_dim_1 = data.this_num_cols;
600 data.P_target_row_data = gmls._Z.data();
601
602 data.RHS_data = gmls._RHS.data();
603 getRHSDims(gmls._dense_solver_type, gmls._constraint_type, gmls._reconstruction_space,
604 gmls._dimensions, data.max_num_rows, data.this_num_cols, data.RHS_dim_0, data.RHS_dim_1);
605
606 data.P_data = gmls._P.data();
607 getPDims(gmls._dense_solver_type, gmls._constraint_type, gmls._reconstruction_space,
608 gmls._dimensions, data.max_num_rows, data.this_num_cols, data.P_dim_0, data.P_dim_1);
609
610 data.w_data = gmls._w.data();
611
612 if ((gmls._constraint_type == ConstraintType::NO_CONSTRAINT) && (gmls._dense_solver_type != DenseSolverType::LU)) {
613 data.Coeffs_data = gmls._RHS.data();
614 data.Coeffs_dim_0 = data.RHS_dim_0;
615 data.Coeffs_dim_1 = data.RHS_dim_1;
616 } else {
617 data.Coeffs_data = gmls._P.data();
618 data.Coeffs_dim_0 = data.P_dim_1;
619 data.Coeffs_dim_1 = data.P_dim_0;
620 }
621
622 return data;
623}
624
625} // Compadre
#define TO_GLOBAL(variable)
#define compadre_assert_release(condition)
compadre_assert_release is used for assertions that should always be checked, but generally are not e...
SolutionSet< device_memory_space > _d_ss
int _sampling_multiplier
actual dimension of the sampling functional e.g.
neighbor_lists_type * getNeighborLists() const
Get neighbor list accessor.
std::string _quadrature_type
quadrature rule type
Kokkos::View< double * > _T
Rank 3 tensor for high order approximation of tangent vectors for all problems.
point_connections_type _additional_pc
(OPTIONAL) connections between additional points and neighbors
bool _reference_outward_normal_direction_provided
whether or not the reference outward normal directions were provided by the user.
int _basis_multiplier
dimension of the reconstructed function e.g.
ReconstructionSpace _reconstruction_space
reconstruction space for GMLS problems, set at GMLS class instantiation
const GMLSSolutionData extractSolutionData() const
Get GMLS solution data.
int _order_of_quadrature_points
order of exact polynomial integration for quadrature rule
int _global_dimensions
spatial dimension of the points, set at class instantiation only
void generatePolynomialCoefficients(const int number_of_batches=1, const bool keep_coefficients=false, const bool clear_cache=true)
Generates polynomial coefficients by setting up and solving least squares problems Sets up the batch ...
ProblemType _problem_type
problem type for GMLS problem, can also be set to STANDARD for normal or MANIFOLD for manifold proble...
int _data_sampling_multiplier
effective dimension of the data sampling functional e.g.
DenseSolverType _dense_solver_type
solver type for GMLS problem - can be QR, SVD or LU
Kokkos::View< TargetOperation * >::host_mirror_type _host_operations
vector containing target functionals to be applied for reconstruction problem (host)
Kokkos::View< double * > _Z
stores evaluations of targets applied to basis
int _initial_index_for_batch
initial index for current batch
int _NP
dimension of basis for polynomial reconstruction
void generateAlphas(const int number_of_batches=1, const bool keep_coefficients=false, const bool clear_cache=true)
Meant to calculate target operations and apply the evaluations to the previously constructed polynomi...
int _poly_order
order of basis for polynomial reconstruction
bool verifyPointConnections(bool assert_valid=false)
Verify whether _pc is valid.
Kokkos::View< double * > _manifold_curvature_coefficients
curvature polynomial coefficients for all problems
SamplingFunctional _data_sampling_functional
generally the same as _polynomial_sampling_functional, but can differ if specified at
SolutionSet< host_memory_space > _h_ss
Solution Set (contains all alpha values from solution and alpha layout methods)
ConstraintType _constraint_type
constraint type for GMLS problem
int _dimensions
dimension of the problem, set at class instantiation only
int _curvature_poly_order
order of basis for curvature reconstruction
bool _use_reference_outward_normal_direction_provided_to_orient_surface
whether or not to use reference outward normal directions to orient the surface in a manifold problem...
Kokkos::View< double * >::host_mirror_type _host_T
tangent vectors information (host)
ParallelManager _pm
determines scratch level spaces and is used to call kernels
Kokkos::View< TargetOperation * > _operations
vector containing target functionals to be applied for reconstruction problem (device)
bool _orthonormal_tangent_space_provided
whether or not the orthonormal tangent directions were provided by the user.
Kokkos::View< double * > _RHS
sqrt(w)*Identity matrix for all problems, later holds polynomial coefficients for all problems
SamplingFunctional _polynomial_sampling_functional
polynomial sampling functional used to construct P matrix, set at GMLS class instantiation NOTE: ca...
Kokkos::View< double * > _P
P*sqrt(w) matrix for all problems.
static KOKKOS_INLINE_FUNCTION int getNP(const int m, const int dimension=3, const ReconstructionSpace r_space=ReconstructionSpace::ScalarTaylorPolynomial)
Returns size of the basis for a given polynomial order and dimension General to dimension 1....
const GMLSBasisData extractBasisData() const
Get GMLS basis data.
Kokkos::View< constdouble *****, layout_right >::host_mirror_type _host_prestencil_weights
generated weights for nontraditional samples required to transform data into expected sampling functi...
point_connections_type _pc
connections between points and neighbors
Kokkos::View< double *****, layout_right > _prestencil_weights
generated weights for nontraditional samples required to transform data into expected sampling functi...
Kokkos::View< double * > _w
contains weights for all problems
bool _entire_batch_computed_at_once
whether entire calculation was computed at once the alternative is that it was broken up over many sm...
int _local_dimensions
dimension of the problem, set at class instantiation only. For manifolds, generally _global_dimension...
bool verifyAdditionalPointConnections(bool assert_valid=false)
Verify whether _additional_pc is valid.
bool _store_PTWP_inv_PTW
whether polynomial coefficients were requested to be stored (in a state not yet applied to data)
int _dimension_of_quadrature_points
dimension of quadrature rule
Quadrature _qm
manages and calculates quadrature
void setThreadScratchSize(const int level, const int value)
void setTeamScratchSize(const int level, const int value)
Kokkos::TeamPolicy< device_execution_space > TeamPolicyThreadsAndVectors(const global_index_type batch_size, const int threads_per_team=-1, const int vector_lanes_per_thread=-1) const
Creates a team policy for a parallel_for parallel_for will break out over loops over teams with each ...
template void batchQRPivotingSolve< layout_right, layout_left, layout_right >(ParallelManager, double *, int, int, double *, int, int, int, int, int, const int, const bool)
template void batchQRPivotingSolve< layout_right, layout_right, layout_right >(ParallelManager, double *, int, int, double *, int, int, int, int, int, const int, const bool)
constexpr SamplingFunctional VaryingManifoldVectorPointSample
For integrating polynomial dotted with normal over an edge.
KOKKOS_INLINE_FUNCTION int getAdditionalCoeffSizeFromConstraintAndSpace(DenseSolverType dense_solver_type, ConstraintType constraint_type, ReconstructionSpace reconstruction_space, const int dimension)
KOKKOS_INLINE_FUNCTION void getRHSDims(DenseSolverType dense_solver_type, ConstraintType constraint_type, ReconstructionSpace reconstruction_space, const int dimension, const int M, const int N, int &RHS_row, int &RHS_col)
@ NEUMANN_GRAD_SCALAR
Neumann Gradient Scalar Type.
@ NO_CONSTRAINT
No constraint.
constexpr SamplingFunctional PointSample
Available sampling functionals.
@ BernsteinPolynomial
Bernstein polynomial basis.
@ ScalarTaylorPolynomial
Scalar polynomial basis centered at the target site and scaled by sum of basis powers e....
KOKKOS_INLINE_FUNCTION int calculateBasisMultiplier(const ReconstructionSpace rs, const int local_dimensions)
Calculate basis_multiplier.
KOKKOS_INLINE_FUNCTION int calculateSamplingMultiplier(const ReconstructionSpace rs, const SamplingFunctional sro, const int local_dimensions)
Calculate sampling_multiplier.
constexpr SamplingFunctional StaggeredEdgeAnalyticGradientIntegralSample
Analytical integral of a gradient source vector is just a difference of the scalar source at neighbor...
KOKKOS_INLINE_FUNCTION int getAdditionalAlphaSizeFromConstraint(DenseSolverType dense_solver_type, ConstraintType constraint_type)
@ MANIFOLD
Solve GMLS problem on a manifold (will use QR or SVD to solve the resultant GMLS problem dependent on...
KOKKOS_INLINE_FUNCTION void getPDims(DenseSolverType dense_solver_type, ConstraintType constraint_type, ReconstructionSpace reconstruction_space, const int dimension, const int M, const int N, int &out_row, int &out_col)
@ LU
LU factorization performed on P^T*W*P matrix.
@ QR
QR+Pivoting factorization performed on P*sqrt(w) matrix.
@ DifferentEachNeighbor
Each target applies a different transform for each neighbor.
@ DifferentEachTarget
Each target applies a different data transform, but the same to each neighbor.
std::size_t global_index_type
KOKKOS_INLINE_FUNCTION int getOutputDimensionOfSampling(SamplingFunctional sro, const int local_dimensions)
Dimensions ^ output rank for sampling operation (always in local chart if on a manifold,...
Functor to evaluate curvature targets and apply to coefficients of curvature reconstruction.
Functor to apply target evaluation to polynomial coefficients to store in _alphas.
Functor to assemble the P*sqrt(weights) matrix and construct sqrt(weights)*Identity for curvature.
Functor to assemble the P*sqrt(weights) matrix and construct sqrt(weights)*Identity.
Functor to assemble the P*sqrt(weights) matrix and construct sqrt(weights)*Identity.
Functor to create a coarse tangent approximation from a given neighborhood of points.
Functor to calculate prestencil weights to apply to data to transform into a format expected by a GML...
Functor to evaluate targets on a manifold.
Functor to evaluate targets operations on the basis.
Functor to determine if tangent directions need reordered, and to reorder them if needed We require t...
Functor to evaluate curvature targets and construct accurate tangent direction approximation for mani...
global_index_type getTotalNeighborsOverAllListsHost() const
Get the sum of the number of neighbors of all targets' neighborhoods (host)
KOKKOS_INLINE_FUNCTION int getMaxNumNeighbors() const
Get the maximum number of neighbors of all targets' neighborhoods (host/device)
KOKKOS_INLINE_FUNCTION int getNumberOfTargets() const
Get number of total targets having neighborhoods (host/device).
device_mirror_target_view_type _target_coordinates
All vairables and functionality related to the layout and storage of GMLS solutions (alpha values)
Kokkos::View< TargetOperation *, memory_space > _lro
vector of user requested target operations
NeighborLists< Kokkos::View< int * > > _neighbor_lists
Accessor to get neighbor list data, offset data, and number of neighbors per target.
int _max_evaluation_sites_per_target
maximum number of evaluation sites for each target (includes target site)
Kokkos::View< double *, layout_right, memory_space > _alphas
generated alpha coefficients (device)
int _added_alpha_size
additional alpha coefficients due to constraints
bool _contains_valid_alphas
whether internal alpha values are valid (set externally on a solve)
int _total_alpha_values
used for sizing P_target_row and the _alphas view