Line data Source code
1 : //* This file is part of the MOOSE framework
2 : //* https://www.mooseframework.org
3 : //*
4 : //* All rights reserved, see COPYRIGHT for full restrictions
5 : //* https://github.com/idaholab/moose/blob/master/COPYRIGHT
6 : //*
7 : //* Licensed under LGPL 2.1, please see LICENSE for details
8 : //* https://www.gnu.org/licenses/lgpl-2.1.html
9 :
10 : #pragma once
11 :
12 : #include "KokkosTypes.h"
13 : #include "KokkosVector.h"
14 : #include "KokkosMatrix.h"
15 : #include "KokkosAssembly.h"
16 :
17 : #include "SystemBase.h"
18 : #include "PerfGraphInterface.h"
19 :
20 : #include "libmesh/communicator.h"
21 :
22 : class MooseMesh;
23 : class SystemBase;
24 :
25 : namespace Moose::Kokkos
26 : {
27 :
28 : class NodalBCBase;
29 :
30 : /**
31 : * The Kokkos system class. Each nonlinear and auxiliary system in MOOSE has a corresponding Kokkos
32 : * system.
33 : */
34 : class System : public PerfGraphInterface, public MeshHolder, public AssemblyHolder
35 : {
36 : public:
37 : /**
38 : * Constructor
39 : * @param system The associated MOOSE system
40 : */
41 : System(SystemBase & system);
42 :
43 : #ifdef MOOSE_KOKKOS_SCOPE
44 : /**
45 : * Synchronize the active tagged vectors and matrices between host and device
46 : * @param dir Copy direction
47 : */
48 : void sync(const MemcpyType dir);
49 : /**
50 : * Synchronize the specified tagged vectors between host and device
51 : * @param tags The vector tags
52 : * @param dir Copy direction
53 : */
54 : ///@{
55 : void sync(const std::set<TagID> & tags, const MemcpyType dir);
56 : void sync(const std::vector<TagID> & tags, const MemcpyType dir);
57 : void sync(const TagID tag, const MemcpyType dir);
58 : ///@}
59 : /**
60 : * Allocate the quadrature point vectors for active variable and tags and cache
61 : * quadrature point values
62 : */
63 : void reinit();
64 : /**
65 : * Set the active variables
66 : * @param vars The active MOOSE variables
67 : */
68 : void setActiveVariables(const std::set<MooseVariableFieldBase *> & vars);
69 : /**
70 : * Set the active solution tags
71 : * @param tags The active solution tags
72 : */
73 : void setActiveSolutionTags(const std::set<TagID> & tags);
74 : /**
75 : * Set the active residual tags
76 : * @param tags The active residual tags
77 : */
78 : void setActiveResidualTags(const std::set<TagID> & tags);
79 : /**
80 : * Set the active matrix tags
81 : * @param vars The active matrix tags
82 : */
83 : void setActiveMatrixTags(const std::set<TagID> & tags);
84 : /**
85 : * Clear the cached active variables
86 : */
87 190102 : void clearActiveVariables() { _active_variables.destroy(); }
88 : /**
89 : * Clear the cached active solution tags
90 : */
91 319784 : void clearActiveSolutionTags() { _active_solution_tags.destroy(); }
92 : /**
93 : * Clear the cached active residual tags
94 : */
95 131257 : void clearActiveResidualTags()
96 : {
97 131257 : _active_residual_tags.destroy();
98 131257 : _residual_tag_active = false;
99 131257 : }
100 : /**
101 : * Clear the cached active matrix tags
102 : */
103 15791 : void clearActiveMatrixTags()
104 : {
105 15791 : _active_matrix_tags.destroy();
106 15791 : _matrix_tag_active = false;
107 15791 : }
108 : /**
109 : * Get the MOOSE system
110 : * @returns The MOOSE system
111 : */
112 : ///@{
113 : auto & getSystem() { return _system; }
114 : const auto & getSystem() const { return _system; }
115 : ///@}
116 : /**
117 : * Get the libMesh DOF map
118 : * @returns The libMesh DOF map
119 : */
120 376947 : const auto & getDofMap() const { return _system.dofMap(); }
121 : /**
122 : * Get the libMesh communicator
123 : * @returns The libMesh communicator
124 : */
125 : const auto & getComm() const { return _comm; }
126 : /**
127 : * Get the list of local DOF indices to communicate
128 : * @returns The list of local DOF indices to communicate
129 : */
130 : const auto & getLocalCommList() const { return _local_comm_list; }
131 : /**
132 : * Get the list of ghost DOF indices to communicate
133 : * @returns The list of ghost DOF indices to communicate
134 : */
135 : const auto & getGhostCommList() const { return _ghost_comm_list; }
136 : /**
137 : * Get the sparisty pattern data
138 : * @returns The sparisty pattern data
139 : */
140 1736 : const auto & getSparsity() const { return _sparsity; }
141 : /**
142 : * Get the list of off-diagonal coupled variable numbers of a variable
143 : * @param var The variable number
144 : * @returns The list of off-diagonal coupled variable numbers
145 : */
146 259992 : KOKKOS_FUNCTION const auto & getCoupling(unsigned int var) const { return _coupling[var]; }
147 : /**
148 : * Check whether a variable is active on a subdomain
149 : * @param var The variable number
150 : * @param subdomain The contiguous subdomain ID
151 : * @returns Whether the variable is active
152 : */
153 496600 : KOKKOS_FUNCTION bool isVariableActive(unsigned int var, ContiguousSubdomainID subdomain) const
154 : {
155 496600 : return _var_subdomain_active(var, subdomain);
156 : }
157 : /**
158 : * Check whether a residual tag is active
159 : * @param tag The residual tag
160 : * @returns Whether the residual tag is active
161 : */
162 18658064 : KOKKOS_FUNCTION bool isResidualTagActive(TagID tag) const { return _residual_tag_active[tag]; }
163 : /**
164 : * Check whether a matrix tag is active
165 : * @param tag The matrix tag
166 : * @returns Whether the matrix tag is active
167 : */
168 17535073 : KOKKOS_FUNCTION bool isMatrixTagActive(TagID tag) const { return _matrix_tag_active[tag]; }
169 : /**
170 : * Check whether a local DOF index is associated with a nodal BC for an extra matrix tag
171 : * @param dof The local DOF index
172 : * @param tag The extra matrix tag
173 : * @returns Whether the local DOF index is covered by a nodal BC
174 : */
175 11463133 : KOKKOS_FUNCTION bool hasNodalBCMatrixTag(dof_id_type dof, TagID tag) const
176 : {
177 11463133 : return _nbc_matrix_tag_dof[tag].isAlloc() && _nbc_matrix_tag_dof[tag][dof];
178 : }
179 : /**
180 : * Get the FE type ID of a variable
181 : * @param var The variable number
182 : * @returns The FE type ID
183 : */
184 48591736 : KOKKOS_FUNCTION unsigned int getFETypeID(unsigned int var) const { return _var_fe_types[var]; }
185 : /**
186 : * Get the number of local DOFs
187 : * @returns The number of local DOFs
188 : */
189 16793 : KOKKOS_FUNCTION dof_id_type getNumLocalDofs() const { return _num_local_dofs; }
190 : /**
191 : * Get the number of ghost DOFs
192 : * @returns The number of ghost DOFs
193 : */
194 387639 : KOKKOS_FUNCTION dof_id_type getNumGhostDofs() const { return _num_ghost_dofs; }
195 : /**
196 : * Get the maximum number of DOFs per element of a variable
197 : * This number is local to each process
198 : * @param var The variable number
199 : * @returns The maximum number of DOFs per element
200 : */
201 : KOKKOS_FUNCTION unsigned int getMaxDofsPerElem(unsigned int var) const
202 : {
203 : return _max_dofs_per_elem[var];
204 : }
205 : /**
206 : * Get the local DOF index of a variable for an element
207 : * @param elem The contiguous element ID
208 : * @param i The element-local DOF index
209 : * @param var The variable number
210 : * @returns The local DOF index
211 : */
212 350587812 : KOKKOS_FUNCTION dof_id_type getElemLocalDofIndex(ContiguousElementID elem,
213 : unsigned int i,
214 : unsigned int var) const
215 : {
216 350587812 : return _local_elem_dof_index[var](elem, i);
217 : }
218 : /**
219 : * Get the local DOF index of a variable for a node
220 : * @param node The contiguous node ID
221 : * @param idx The node-local DOF index
222 : * @param var The variable number
223 : * @returns The local DOF index
224 : */
225 6825647 : KOKKOS_FUNCTION dof_id_type getNodeLocalDofIndex(ContiguousNodeID node,
226 : unsigned int i,
227 : unsigned int var) const
228 : {
229 6825647 : return _local_node_dof_index[var][node] + i;
230 : }
231 : /**
232 : * Get the global DOF index of a variable for an element
233 : * @param elem The contiguous element ID
234 : * @param i The element-local DOF index
235 : * @param var The variable number
236 : * @returns The global DOF index
237 : */
238 9967885 : KOKKOS_FUNCTION dof_id_type getElemGlobalDofIndex(ContiguousElementID elem,
239 : unsigned int i,
240 : unsigned int var) const
241 : {
242 9967885 : return _local_to_global_dof_index[_local_elem_dof_index[var](elem, i)];
243 : }
244 : /**
245 : * Get the global DOF index of a variable for a node
246 : * @param node The contiguous node ID
247 : * @param var The variable number
248 : * @returns The global DOF index
249 : */
250 342407 : KOKKOS_FUNCTION dof_id_type getNodeGlobalDofIndex(ContiguousNodeID node, unsigned int var) const
251 : {
252 342407 : return _local_to_global_dof_index[_local_node_dof_index[var][node]];
253 : }
254 : /**
255 : * Get the global DOF index of a local DOF index
256 : * @param dof The local DOF index
257 : * @returns The global DOF index
258 : */
259 : KOKKOS_FUNCTION dof_id_type localToGlobalDofIndex(dof_id_type dof) const
260 : {
261 : return _local_to_global_dof_index[dof];
262 : }
263 : /**
264 : * Get whether a variable is defined on a node
265 : * @param node The contiguous node ID
266 : * @param var The variable number
267 : * @returns Whether the variable is defined on the node
268 : */
269 4918893 : KOKKOS_FUNCTION bool isNodalDefined(ContiguousNodeID node, unsigned int var) const
270 : {
271 4918893 : return _local_node_dof_index[var][node] != libMesh::DofObject::invalid_id;
272 : }
273 : /**
274 : * Get a tagged Kokkos vector
275 : * @param tag The vector tag
276 : * @returns The Kokkos vector
277 : */
278 : KOKKOS_FUNCTION auto & getVector(TagID tag) const { return _vectors[tag]; }
279 : /**
280 : * Get a tagged Kokkos matrix
281 : * @param tag The matrix tag
282 : * @returns The Kokkos matrix
283 : */
284 426675 : KOKKOS_FUNCTION auto & getMatrix(TagID tag) const { return _matrices[tag]; }
285 : /**
286 : * Get the DOF value of a tagged vector
287 : * @param dof The local DOF index
288 : * @param tag The vector tag
289 : * @returns The DOF value
290 : */
291 329597939 : KOKKOS_FUNCTION Real & getVectorDofValue(const dof_id_type dof, const TagID tag) const
292 : {
293 329597939 : return _vectors[tag][dof];
294 : }
295 : /**
296 : * Get the DOF value of a tagged vector for automatic differentiation (AD)
297 : * @param dof The local DOF index
298 : * @param tag The vector tag
299 : * @param seed The derivative seed
300 : * @returns The DOF AD value with optional seed derivative
301 : */
302 : KOKKOS_FUNCTION ADReal getVectorDofADValue(const dof_id_type dof,
303 : const TagID tag,
304 : const Real seed) const;
305 : /**
306 : * Get the quadrature value of a variable from a tagged vector
307 : * @param info The element information object
308 : * @param qp The global quadrature point index
309 : * @param var The variable number
310 : * @param tag The vector tag
311 : * @returns The quadrature value
312 : */
313 101597468 : KOKKOS_FUNCTION Real & getVectorQpValue(const ElementInfo info,
314 : const dof_id_type qp,
315 : const unsigned int var,
316 : const TagID tag) const
317 : {
318 101597468 : return _qp_solutions[tag](info.subdomain, var)[qp];
319 : }
320 : /**
321 : * Get the quadrature value of a variable from a tagged vector for automatic differentiation (AD)
322 : * @param info The element information object
323 : * @param offset The offset into the global quadrature point index
324 : * @param qp The local quadrature point index
325 : * @param var The variable number
326 : * @param tag The vector tag
327 : * @param seed The derivative seed
328 : * @returns The quadrature AD value
329 : */
330 : KOKKOS_FUNCTION ADReal getVectorQpADValue(const ElementInfo info,
331 : const dof_id_type offset,
332 : const dof_id_type qp,
333 : const unsigned int var,
334 : const TagID tag,
335 : const Real seed) const;
336 : /**
337 : * Get the quadrature gradient of a variable from a tagged vector
338 : * @param info The element information object
339 : * @param qp The global quadrature point index
340 : * @param var The variable number
341 : * @param tag The vector tag
342 : * @returns The quadrature gradient
343 : */
344 108611240 : KOKKOS_FUNCTION Real3 & getVectorQpGrad(const ElementInfo info,
345 : const dof_id_type qp,
346 : const unsigned int var,
347 : const TagID tag) const
348 : {
349 108611240 : return _qp_solutions_grad[tag](info.subdomain, var)[qp];
350 : }
351 : /**
352 : * Get the quadrature gradient of a variable from a tagged vector for automatic differentiation
353 : * (AD)
354 : * @param info The element information object
355 : * @param jacobian The inverse Jacobian matrix
356 : * @param offset The offset into the global quadrature point index
357 : * @param qp The local quadrature point index
358 : * @param var The variable number
359 : * @param tag The vector tag
360 : * @param seed The derivative seed
361 : * @returns The quadrature AD gradient
362 : */
363 : KOKKOS_FUNCTION ADReal3 getVectorQpADGrad(const ElementInfo info,
364 : const Real33 jacobian,
365 : const dof_id_type offset,
366 : const dof_id_type qp,
367 : const unsigned int var,
368 : const TagID tag,
369 : const Real seed) const;
370 : /**
371 : * Get the face quadrature value of a variable from a tagged vector
372 : * @param info The element information object
373 : * @param side The side index
374 : * @param qp The local quadrature point index
375 : * @param var The vriable number
376 : * @param tag The vector tag
377 : * @returns The face quadrature value
378 : */
379 : KOKKOS_FUNCTION Real getVectorQpValueFace(const ElementInfo info,
380 : const unsigned int side,
381 : const unsigned int qp,
382 : const unsigned int var,
383 : const TagID tag) const;
384 : /**
385 : * Get the face quadrature value of a variable from a tagged vector for automatic differentiation
386 : * (AD)
387 : * @param info The element information object
388 : * @param side The side index
389 : * @param qp The local quadrature point index
390 : * @param var The vriable number
391 : * @param tag The vector tag
392 : * @param seed The derivative seed
393 : * @returns The face quadrature AD value
394 : */
395 : KOKKOS_FUNCTION ADReal getVectorQpADValueFace(const ElementInfo info,
396 : const unsigned int side,
397 : const unsigned int qp,
398 : const unsigned int var,
399 : const TagID tag,
400 : const Real seed) const;
401 : /**
402 : * Get the face quadrature gradient of a variable from a tagged vector
403 : * @param info The element information object
404 : * @param side The side index
405 : * @param jacobian The inverse Jacobian matrix
406 : * @param qp The local quadrature point index
407 : * @param var The variable number
408 : * @param tag The vector tag
409 : * @returns The face quadrature gradient
410 : */
411 : KOKKOS_FUNCTION Real3 getVectorQpGradFace(const ElementInfo info,
412 : const unsigned int side,
413 : const Real33 jacobian,
414 : const unsigned int qp,
415 : const unsigned int var,
416 : const TagID tag) const;
417 : /**
418 : * Get the face quadrature gradient of a variable from a tagged vector for automatic
419 : * differentiation (AD)
420 : * @param info The element information object
421 : * @param side The side index
422 : * @param jacobian The inverse Jacobian matrix
423 : * @param qp The local quadrature point index
424 : * @param var The variable number
425 : * @param tag The vector tag
426 : * @param seed The derivative seed
427 : * @returns The face quadrature AD gradient
428 : */
429 : KOKKOS_FUNCTION ADReal3 getVectorQpADGradFace(const ElementInfo info,
430 : const unsigned int side,
431 : const Real33 jacobian,
432 : const unsigned int qp,
433 : const unsigned int var,
434 : const TagID tag,
435 : const Real seed) const;
436 : /**
437 : * Get an entry from a tagged matrix
438 : * @param row The local row index
439 : * @param col The global column index
440 : * @param tag The matrix tag
441 : * @returns The entry from the tagged matrix
442 : */
443 13276073 : KOKKOS_FUNCTION Real & getMatrixValue(dof_id_type row, dof_id_type col, TagID tag) const
444 : {
445 13276073 : return _matrices[tag](row, col);
446 : }
447 :
448 : /**
449 : * Kokkos function for caching variable values on element quadrature points
450 : */
451 : KOKKOS_FUNCTION void operator()(const ThreadID tid) const;
452 : #endif
453 :
454 : /**
455 : * CSR format sparsity data
456 : */
457 : struct Sparsity
458 : {
459 : Array<PetscInt> col_idx;
460 : Array<PetscInt> row_idx;
461 : Array<PetscInt> row_ptr;
462 : };
463 :
464 : private:
465 : /**
466 : * Setup variable data
467 : */
468 : void setupVariables();
469 : /**
470 : * Setup DOF data
471 : */
472 : void setupDofs();
473 : /**
474 : * Setup sparsity data
475 : */
476 : void setupSparsity();
477 : /**
478 : * Mark the DOFs covered by nodal BCs
479 : */
480 : void setupNodalBCDofs();
481 : /**
482 : * Get the list of DOFs covered by a nodal BC
483 : * @param nbc The Kokkos nodal BC object
484 : * @param dofs The flag whether each DOF is covered by the nodal BC
485 : */
486 : void getNodalBCDofs(const NodalBCBase * nbc, Array<bool> & dofs);
487 :
488 : /**
489 : * Kokkos thread object
490 : */
491 : Thread<> _thread;
492 : /**
493 : * Reference of the MOOSE system
494 : */
495 : SystemBase & _system;
496 : /**
497 : * Reference of the MOOSE mesh
498 : */
499 : const MooseMesh & _mesh;
500 : /**
501 : * Reference of the libMesh DOF map
502 : */
503 : const libMesh::DofMap & _dof_map;
504 : /**
505 : * Reference of the libMesh communicator
506 : */
507 : const Parallel::Communicator & _comm;
508 : /**
509 : * Number of variables
510 : */
511 : const unsigned int _num_vars;
512 : /**
513 : * Number of local DOFs
514 : */
515 : const dof_id_type _num_local_dofs;
516 : /**
517 : * Number of ghost DOFs
518 : */
519 : const dof_id_type _num_ghost_dofs;
520 :
521 : /**
522 : * Kokkos vectors and matrices on device
523 : */
524 : ///@{
525 : Array<Vector> _vectors;
526 : Array<Matrix> _matrices;
527 : ///@}
528 : /**
529 : * Cached elemental quadrature values and gradients
530 : */
531 : ///@{
532 : Array<Array2D<Array<Real>>> _qp_solutions;
533 : Array<Array2D<Array<Real3>>> _qp_solutions_grad;
534 : ///@}
535 : /**
536 : * Local DOF index of each variable
537 : */
538 : ///@{
539 : Array<Array2D<dof_id_type>> _local_elem_dof_index;
540 : Array<Array<dof_id_type>> _local_node_dof_index;
541 : ///@}
542 : /**
543 : * Map from local DOF index to global DOF index
544 : */
545 : Array<dof_id_type> _local_to_global_dof_index;
546 : /**
547 : * Maximum number of DOFs per element for each variable
548 : */
549 : Array<unsigned int> _max_dofs_per_elem;
550 : /**
551 : * FE type ID of each variable
552 : */
553 : Array<unsigned int> _var_fe_types;
554 : /**
555 : * Whether each variable is active on subdomains
556 : */
557 : Array2D<bool> _var_subdomain_active;
558 : /**
559 : * Off-diagonal coupled variable numbers of each variable
560 : */
561 : Array<Array<unsigned int>> _coupling;
562 :
563 : /**
564 : * List of active variable numbers
565 : */
566 : Array<unsigned int> _active_variables;
567 : /**
568 : * List of active tags
569 : */
570 : ///@{
571 : Array<TagID> _active_solution_tags;
572 : Array<TagID> _active_residual_tags;
573 : Array<TagID> _active_matrix_tags;
574 : ///@}
575 : /**
576 : * Flag whether each tag is active
577 : */
578 : ///@{
579 : Array<bool> _residual_tag_active;
580 : Array<bool> _matrix_tag_active;
581 : ///@}
582 : /**
583 : * Flag whether each DOF is covered by a nodal BC for each matrix tag
584 : */
585 : Array<Array<bool>> _nbc_matrix_tag_dof;
586 : /**
587 : * List of DOFs to send and receive
588 : */
589 : ///@{
590 : Array<Array<dof_id_type>> _local_comm_list;
591 : Array<Array<dof_id_type>> _ghost_comm_list;
592 : ///@}
593 :
594 : /**
595 : * Matrix sparsity pattern data
596 : */
597 : Sparsity _sparsity;
598 : };
599 :
600 : #ifdef MOOSE_KOKKOS_SCOPE
601 : KOKKOS_FUNCTION inline ADReal
602 17390940 : System::getVectorDofADValue(const dof_id_type dof, TagID tag, const Real seed) const
603 : {
604 17390940 : ADReal value = _vectors[tag][dof];
605 :
606 17390940 : if (seed != 0)
607 17298214 : value.derivatives().insert(_local_to_global_dof_index[dof]) = seed;
608 :
609 17390940 : return value;
610 0 : }
611 :
612 : KOKKOS_FUNCTION inline ADReal
613 2270400 : System::getVectorQpADValue(const ElementInfo info,
614 : const dof_id_type offset,
615 : const dof_id_type qp,
616 : const unsigned int var,
617 : const TagID tag,
618 : const Real seed) const
619 : {
620 2270400 : ADReal value = 0;
621 :
622 2270400 : if (seed == 0)
623 16800 : value = getVectorQpValue(info, offset + qp, var, tag);
624 : else
625 : {
626 2253600 : auto fe = _var_fe_types[var];
627 2253600 : auto n_dofs = kokkosAssembly().getNumDofs(info.type, fe);
628 2253600 : auto & phi = kokkosAssembly().getPhi(info.subdomain, info.type, fe);
629 :
630 11268000 : for (unsigned int i = 0; i < n_dofs; ++i)
631 9014400 : value += getVectorDofADValue(getElemLocalDofIndex(info.id, i, var), tag, seed) * phi(i, qp);
632 : }
633 :
634 2270400 : return value;
635 0 : }
636 :
637 : KOKKOS_FUNCTION inline ADReal3
638 2369344 : System::getVectorQpADGrad(const ElementInfo info,
639 : const Real33 jacobian,
640 : const dof_id_type offset,
641 : const dof_id_type qp,
642 : const unsigned int var,
643 : const TagID tag,
644 : const Real seed) const
645 : {
646 2369344 : ADReal3 grad;
647 :
648 2369344 : if (seed == 0)
649 276560 : grad = getVectorQpGrad(info, offset + qp, var, tag);
650 : else
651 : {
652 2092784 : auto fe = _var_fe_types[var];
653 2092784 : auto n_dofs = kokkosAssembly().getNumDofs(info.type, fe);
654 2092784 : auto & grad_phi = kokkosAssembly().getGradPhi(info.subdomain, info.type, fe);
655 :
656 10341680 : for (unsigned int i = 0; i < n_dofs; ++i)
657 8248896 : grad += getVectorDofADValue(getElemLocalDofIndex(info.id, i, var), tag, seed) *
658 16497792 : (jacobian * grad_phi(i, qp));
659 : }
660 :
661 2369344 : return grad;
662 0 : }
663 :
664 : KOKKOS_FUNCTION inline Real
665 265208 : System::getVectorQpValueFace(const ElementInfo info,
666 : const unsigned int side,
667 : const unsigned int qp,
668 : const unsigned int var,
669 : const TagID tag) const
670 : {
671 265208 : auto fe = _var_fe_types[var];
672 265208 : auto n_dofs = kokkosAssembly().getNumDofs(info.type, fe);
673 265208 : auto & phi = kokkosAssembly().getPhiFace(info.subdomain, info.type, fe)(side);
674 :
675 265208 : Real value = 0;
676 :
677 1334616 : for (unsigned int i = 0; i < n_dofs; ++i)
678 1069408 : value += getVectorDofValue(getElemLocalDofIndex(info.id, i, var), tag) * phi(i, qp);
679 :
680 265208 : return value;
681 : }
682 :
683 : KOKKOS_FUNCTION inline ADReal
684 15156 : System::getVectorQpADValueFace(const ElementInfo info,
685 : const unsigned int side,
686 : const unsigned int qp,
687 : const unsigned int var,
688 : const TagID tag,
689 : const Real seed) const
690 : {
691 15156 : auto fe = _var_fe_types[var];
692 15156 : auto n_dofs = kokkosAssembly().getNumDofs(info.type, fe);
693 15156 : auto & phi = kokkosAssembly().getPhiFace(info.subdomain, info.type, fe)(side);
694 :
695 15156 : ADReal value = 0;
696 :
697 73948 : for (unsigned int i = 0; i < n_dofs; ++i)
698 58792 : value += getVectorDofADValue(getElemLocalDofIndex(info.id, i, var), tag, seed) * phi(i, qp);
699 :
700 15156 : return value;
701 0 : }
702 :
703 : KOKKOS_FUNCTION inline Real3
704 0 : System::getVectorQpGradFace(const ElementInfo info,
705 : const unsigned int side,
706 : const Real33 jacobian,
707 : const unsigned int qp,
708 : const unsigned int var,
709 : const TagID tag) const
710 : {
711 0 : auto fe = _var_fe_types[var];
712 0 : auto n_dofs = kokkosAssembly().getNumDofs(info.type, fe);
713 0 : auto & grad_phi = kokkosAssembly().getGradPhiFace(info.subdomain, info.type, fe)(side);
714 :
715 0 : Real3 grad = 0;
716 :
717 0 : for (unsigned int i = 0; i < n_dofs; ++i)
718 0 : grad += getVectorDofValue(getElemLocalDofIndex(info.id, i, var), tag) *
719 0 : (jacobian * grad_phi(i, qp));
720 :
721 0 : return grad;
722 : }
723 :
724 : KOKKOS_FUNCTION inline ADReal3
725 0 : System::getVectorQpADGradFace(const ElementInfo info,
726 : const unsigned int side,
727 : const Real33 jacobian,
728 : const unsigned int qp,
729 : const unsigned int var,
730 : const TagID tag,
731 : const Real seed) const
732 : {
733 0 : auto fe = _var_fe_types[var];
734 0 : auto n_dofs = kokkosAssembly().getNumDofs(info.type, fe);
735 0 : auto & grad_phi = kokkosAssembly().getGradPhiFace(info.subdomain, info.type, fe)(side);
736 :
737 0 : ADReal3 grad = ADReal(0);
738 :
739 0 : for (unsigned int i = 0; i < n_dofs; ++i)
740 0 : grad += getVectorDofADValue(getElemLocalDofIndex(info.id, i, var), tag, seed) *
741 0 : (jacobian * grad_phi(i, qp));
742 :
743 0 : return grad;
744 0 : }
745 : #endif
746 :
747 : /**
748 : * The Kokkos interface that holds the host reference of the Kokkos systems and copies it to device
749 : * during parallel dispatch.
750 : * Maintains synchronization between host and device Kokkos systems and provides access to the
751 : * appropriate Kokkos systems depending on the architecture.
752 : */
753 : class SystemHolder
754 : {
755 : public:
756 : /**
757 : * Constructor
758 : * @param systems The Kokkos systems
759 : */
760 12122 : SystemHolder(Array<System> & systems) : _systems_host(systems), _systems_device(systems) {}
761 : /**
762 : * Copy constructor
763 : */
764 508913 : SystemHolder(const SystemHolder & holder)
765 289687 : : _systems_host(holder._systems_host), _systems_device(holder._systems_host)
766 : {
767 508913 : }
768 :
769 : #ifdef MOOSE_KOKKOS_SCOPE
770 : /**
771 : * Get the const reference of the Kokkos systems
772 : * @returns The const reference of the Kokkos systems depending on the architecture this function
773 : * is being called on
774 : */
775 26196050 : KOKKOS_FUNCTION const Array<System> & kokkosSystems() const
776 : {
777 26196050 : KOKKOS_IF_ON_HOST(return _systems_host;)
778 :
779 26196050 : return _systems_device;
780 : }
781 : /**
782 : * Get the writeable host reference of the Kokkos systems
783 : * @returns The writeable host reference of the Kokkos systems
784 : */
785 : Array<System> & kokkosSystems() { return _systems_host; }
786 : /**
787 : * Get the const reference of a Kokkos system
788 : * @param sys The system number
789 : * @returns The const reference of the Kokkos system depending on the architecture this function
790 : * is being called on
791 : */
792 35506248 : KOKKOS_FUNCTION const System & kokkosSystem(unsigned int sys) const
793 : {
794 35506248 : KOKKOS_IF_ON_HOST(return _systems_host[sys];)
795 :
796 35506248 : return _systems_device[sys];
797 : }
798 : /**
799 : * Get the writeable reference of a Kokkos system
800 : * @param sys The system number
801 : * @returns The writeable host reference of the Kokkos system
802 : */
803 5398 : System & kokkosSystem(unsigned int sys) { return _systems_host[sys]; }
804 : #endif
805 :
806 : private:
807 : /**
808 : * Host reference of the Kokkos systems
809 : */
810 : Array<System> & _systems_host;
811 : /**
812 : * Device copy of the Kokkos systems
813 : */
814 : const Array<System> _systems_device;
815 : };
816 :
817 : } // namespace Moose::Kokkos
|