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 :
19 : #include "libmesh/communicator.h"
20 :
21 : class MooseMesh;
22 : class SystemBase;
23 :
24 : namespace Moose
25 : {
26 : namespace Kokkos
27 : {
28 :
29 : class NodalBCBase;
30 :
31 : /**
32 : * The Kokkos system class. Each nonlinear and auxiliary system in MOOSE has a corresponding Kokkos
33 : * system.
34 : */
35 : class System : public MeshHolder, public AssemblyHolder
36 : {
37 : public:
38 : /**
39 : * Constructor
40 : * @param system The associated MOOSE system
41 : */
42 : System(SystemBase & system);
43 :
44 : #ifdef MOOSE_KOKKOS_SCOPE
45 : /**
46 : * Synchronize the active tagged vectors and matrices between host and device
47 : * @param dir Copy direction
48 : */
49 : void sync(const MemcpyKind dir);
50 : /**
51 : * Synchronize the specified tagged vectors between host and device
52 : * @param tags The vector tags
53 : * @param dir Copy direction
54 : */
55 : ///{@
56 : void sync(const std::set<TagID> & tags, const MemcpyKind dir);
57 : void sync(const std::vector<TagID> & tags, const MemcpyKind dir);
58 : void sync(const TagID tag, const MemcpyKind dir);
59 : ///@}
60 : /**
61 : * Allocate the quadrature point vectors for active variable and tags and cache
62 : * quadrature point values
63 : */
64 : void reinit();
65 : /**
66 : * Set the active variables
67 : * @param vars The active MOOSE variables
68 : */
69 : void setActiveVariables(const std::set<MooseVariableFieldBase *> & vars);
70 : /**
71 : * Set the active solution tags
72 : * @param tags The active solution tags
73 : */
74 : void setActiveSolutionTags(const std::set<TagID> & tags);
75 : /**
76 : * Set the active residual tags
77 : * @param tags The active residual tags
78 : */
79 : void setActiveResidualTags(const std::set<TagID> & tags);
80 : /**
81 : * Set the active matrix tags
82 : * @param vars The active matrix tags
83 : */
84 : void setActiveMatrixTags(const std::set<TagID> & tags);
85 : /**
86 : * Clear the cached active variables
87 : */
88 104008 : void clearActiveVariables() { _active_variables.destroy(); }
89 : /**
90 : * Clear the cached active solution tags
91 : */
92 104008 : void clearActiveSolutionTags() { _active_solution_tags.destroy(); }
93 : /**
94 : * Clear the cached active residual tags
95 : */
96 34353 : void clearActiveResidualTags()
97 : {
98 34353 : _active_residual_tags.destroy();
99 34353 : _residual_tag_active = false;
100 34353 : }
101 : /**
102 : * Clear the cached active matrix tags
103 : */
104 5611 : void clearActiveMatrixTags()
105 : {
106 5611 : _active_matrix_tags.destroy();
107 5611 : _matrix_tag_active = false;
108 5611 : }
109 : /**
110 : * Get the MOOSE system
111 : * @returns The MOOSE system
112 : */
113 : ///{@
114 : auto & getSystem() { return _system; }
115 : const auto & getSystem() const { return _system; }
116 : ///@}
117 : /**
118 : * Get the libMesh DOF map
119 : * @returns The libMesh DOF map
120 : */
121 103099 : const auto & getDofMap() const { return _system.dofMap(); }
122 : /**
123 : * Get the libMesh communicator
124 : * @returns The libMesh communicator
125 : */
126 : const auto & getComm() const { return _comm; }
127 : /**
128 : * Get the list of local DOF indices to communicate
129 : * @returns The list of local DOF indices to communicate
130 : */
131 : const auto & getLocalCommList() const { return _local_comm_list; }
132 : /**
133 : * Get the list of ghost DOF indices to communicate
134 : * @returns The list of ghost DOF indices to communicate
135 : */
136 : const auto & getGhostCommList() const { return _ghost_comm_list; }
137 : /**
138 : * Get the sparisty pattern data
139 : * @returns The sparisty pattern data
140 : */
141 725 : const auto & getSparsity() const { return _sparsity; }
142 : /**
143 : * Get the list of off-diagonal coupled variable numbers of a variable
144 : * @param var The variable number
145 : * @returns The list of off-diagonal coupled variable numbers
146 : */
147 254340 : KOKKOS_FUNCTION const auto & getCoupling(unsigned int var) const { return _coupling[var]; }
148 : /**
149 : * Check whether a variable is active on a subdomain
150 : * @param var The variable number
151 : * @param subdomain The contiguous subdomain ID
152 : * @returns Whether the variable is active
153 : */
154 1270 : KOKKOS_FUNCTION bool isVariableActive(unsigned int var, ContiguousSubdomainID subdomain) const
155 : {
156 1270 : return _var_subdomain_active(var, subdomain);
157 : }
158 : /**
159 : * Check whether a residual tag is active
160 : * @param tag The residual tag
161 : * @returns Whether the residual tag is active
162 : */
163 15209172 : KOKKOS_FUNCTION bool isResidualTagActive(TagID tag) const { return _residual_tag_active[tag]; }
164 : /**
165 : * Check whether a matrix tag is active
166 : * @param tag The matrix tag
167 : * @returns Whether the matrix tag is active
168 : */
169 14174994 : KOKKOS_FUNCTION bool isMatrixTagActive(TagID tag) const { return _matrix_tag_active[tag]; }
170 : /**
171 : * Check whether a local DOF index is covered by a nodal BC
172 : * @param dof The local DOF index
173 : * @returns Whether the local DOF index is covered by a nodal BC
174 : */
175 22037882 : KOKKOS_FUNCTION bool hasNodalBC(dof_id_type dof) const
176 : {
177 22037882 : return _nbc_dof.isAlloc() && _nbc_dof[dof];
178 : }
179 : /**
180 : * Check whether a local DOF index is covered by a nodal BC for an extra residual tag
181 : * @param dof The local DOF index
182 : * @param tag The extra residual tag
183 : * @returns Whether the local DOF index is covered by a nodal BC
184 : */
185 687240 : KOKKOS_FUNCTION bool hasNodalBCResidualTag(dof_id_type dof, TagID tag) const
186 : {
187 687240 : return _nbc_residual_tag_dof[tag].isAlloc() && _nbc_residual_tag_dof[tag][dof];
188 : }
189 : /**
190 : * Check whether a local DOF index is associated with a nodal BC for an extra matrix tag
191 : * @param dof The local DOF index
192 : * @param tag The extra matrix tag
193 : * @returns Whether the local DOF index is covered by a nodal BC
194 : */
195 783200 : KOKKOS_FUNCTION bool hasNodalBCMatrixTag(dof_id_type dof, TagID tag) const
196 : {
197 783200 : return _nbc_matrix_tag_dof[tag].isAlloc() && _nbc_matrix_tag_dof[tag][dof];
198 : }
199 : /**
200 : * Get the FE type ID of a variable
201 : * @param var The variable number
202 : * @returns The FE type ID
203 : */
204 15114166 : KOKKOS_FUNCTION unsigned int getFETypeID(unsigned int var) const { return _var_fe_types[var]; }
205 : /**
206 : * Get the number of local DOFs
207 : * @returns The number of local DOFs
208 : */
209 6674 : KOKKOS_FUNCTION dof_id_type getNumLocalDofs() const { return _num_local_dofs; }
210 : /**
211 : * Get the number of ghost DOFs
212 : * @returns The number of ghost DOFs
213 : */
214 107106 : KOKKOS_FUNCTION dof_id_type getNumGhostDofs() const { return _num_ghost_dofs; }
215 : /**
216 : * Get the maximum number of DOFs per element of a variable
217 : * This number is local to each process
218 : * @param var The variable number
219 : * @returns The maximum number of DOFs per element
220 : */
221 : KOKKOS_FUNCTION unsigned int getMaxDofsPerElem(unsigned int var) const
222 : {
223 : return _max_dofs_per_elem[var];
224 : }
225 : /**
226 : * Get the local DOF index of a variable for an element
227 : * @param elem The contiguous element ID
228 : * @param i The element-local DOF index
229 : * @param var The variable number
230 : * @returns The local DOF index
231 : */
232 271394626 : KOKKOS_FUNCTION dof_id_type getElemLocalDofIndex(ContiguousElementID elem,
233 : unsigned int i,
234 : unsigned int var) const
235 : {
236 271394626 : return _local_elem_dof_index[var](elem, i);
237 : }
238 : /**
239 : * Get the local DOF index of a variable for a node
240 : * @param node The contiguous node ID
241 : * @param idx The node-local DOF index
242 : * @param var The variable number
243 : * @returns The local DOF index
244 : */
245 5840487 : KOKKOS_FUNCTION dof_id_type getNodeLocalDofIndex(ContiguousNodeID node,
246 : unsigned int i,
247 : unsigned int var) const
248 : {
249 5840487 : return _local_node_dof_index[var][node] + i;
250 : }
251 : /**
252 : * Get the global DOF index of a variable for an element
253 : * @param elem The contiguous element ID
254 : * @param i The element-local DOF index
255 : * @param var The variable number
256 : * @returns The global DOF index
257 : */
258 8833638 : KOKKOS_FUNCTION dof_id_type getElemGlobalDofIndex(ContiguousElementID elem,
259 : unsigned int i,
260 : unsigned int var) const
261 : {
262 8833638 : return _local_to_global_dof_index[_local_elem_dof_index[var](elem, i)];
263 : }
264 : /**
265 : * Get the global DOF index of a variable for a node
266 : * @param node The contiguous node ID
267 : * @param var The variable number
268 : * @returns The global DOF index
269 : */
270 334398 : KOKKOS_FUNCTION dof_id_type getNodeGlobalDofIndex(ContiguousNodeID node, unsigned int var) const
271 : {
272 334398 : return _local_to_global_dof_index[_local_node_dof_index[var][node]];
273 : }
274 : /**
275 : * Get the global DOF index of a local DOF index
276 : * @param dof The local DOF index
277 : * @returns The global DOF index
278 : */
279 : KOKKOS_FUNCTION dof_id_type localToGlobalDofIndex(dof_id_type dof) const
280 : {
281 : return _local_to_global_dof_index[dof];
282 : }
283 : /**
284 : * Get whether a variable is defined on a node
285 : * @param node The contiguous node ID
286 : * @param var The variable number
287 : * @returns Whether the variable is defined on the node
288 : */
289 4582962 : KOKKOS_FUNCTION bool isNodalDefined(ContiguousNodeID node, unsigned int var) const
290 : {
291 4582962 : return _local_node_dof_index[var][node] != libMesh::DofObject::invalid_id;
292 : }
293 : /**
294 : * Get a tagged Kokkos vector
295 : * @param tag The vector tag
296 : * @returns The Kokkos vector
297 : */
298 : KOKKOS_FUNCTION auto & getVector(TagID tag) const { return _vectors[tag]; }
299 : /**
300 : * Get a tagged Kokkos matrix
301 : * @param tag The matrix tag
302 : * @returns The Kokkos matrix
303 : */
304 351886 : KOKKOS_FUNCTION auto & getMatrix(TagID tag) const { return _matrices[tag]; }
305 : /**
306 : * Get the DOF value of a tagged vector
307 : * @param dof The local DOF index
308 : * @param tag The vector tag
309 : * @returns The DOF value
310 : */
311 267882485 : KOKKOS_FUNCTION Real & getVectorDofValue(dof_id_type dof, TagID tag) const
312 : {
313 267882485 : return _vectors[tag][dof];
314 : }
315 : /**
316 : * Get the quadrature value of a variable from a tagged vector
317 : * @param info The element information object
318 : * @param qp The global quadrature point index
319 : * @param var The variable number
320 : * @param tag The vector tag
321 : * @returns The quadrature value
322 : */
323 : KOKKOS_FUNCTION Real &
324 83681084 : getVectorQpValue(ElementInfo info, dof_id_type qp, unsigned int var, TagID tag) const
325 : {
326 83681084 : return _qp_solutions[tag](info.subdomain, var)[qp];
327 : }
328 : /**
329 : * Get the quadrature gradient of a variable from a tagged vector
330 : * @param info The element information object
331 : * @param qp The global quadrature point index
332 : * @param var The variable number
333 : * @param tag The vector tag
334 : * @returns The quadrature gradient
335 : */
336 : KOKKOS_FUNCTION Real3 &
337 87995180 : getVectorQpGrad(ElementInfo info, dof_id_type qp, unsigned int var, TagID tag) const
338 : {
339 87995180 : return _qp_solutions_grad[tag](info.subdomain, var)[qp];
340 : }
341 : /**
342 : * Get the face quadrature value of a variable from a tagged vector
343 : * @param info The element information object
344 : * @param side The side index
345 : * @param qp The local quadrature point index
346 : * @param var The vriable number
347 : * @param tag The vector tag
348 : * @returns The face quadrature value
349 : */
350 : KOKKOS_FUNCTION Real getVectorQpValueFace(const ElementInfo info,
351 : const unsigned int side,
352 : const unsigned int qp,
353 : const unsigned int var,
354 : const TagID tag) const;
355 : /**
356 : * Get the face quadrature gradient of a variable from a tagged vector
357 : * @param info The element information object
358 : * @param side The side index
359 : * @param jacobian The inverse Jacobian matrix
360 : * @param qp The local quadrature point index
361 : * @param var The variable number
362 : * @param tag The vector tag
363 : * @returns The face quadrature gradient
364 : */
365 : KOKKOS_FUNCTION Real3 getVectorQpGradFace(const ElementInfo info,
366 : const unsigned int side,
367 : const Real33 jacobian,
368 : const unsigned int qp,
369 : const unsigned int var,
370 : const TagID tag) const;
371 : /**
372 : * Get an entry from a tagged matrix
373 : * @param row The local row index
374 : * @param col The global column index
375 : * @param tag The matrix tag
376 : * @returns The entry from the tagged matrix
377 : */
378 8316600 : KOKKOS_FUNCTION Real & getMatrixValue(dof_id_type row, dof_id_type col, TagID tag) const
379 : {
380 8316600 : return _matrices[tag](row, col);
381 : }
382 :
383 : /**
384 : * Kokkos function for caching variable values on element quadrature points
385 : */
386 : KOKKOS_FUNCTION void operator()(const ThreadID tid) const;
387 : #endif
388 :
389 : /**
390 : * CSR format sparsity data
391 : */
392 : struct Sparsity
393 : {
394 : Array<PetscInt> col_idx;
395 : Array<PetscInt> row_idx;
396 : Array<PetscInt> row_ptr;
397 : };
398 :
399 : private:
400 : /**
401 : * Setup variable data
402 : */
403 : void setupVariables();
404 : /**
405 : * Setup DOF data
406 : */
407 : void setupDofs();
408 : /**
409 : * Setup sparsity data
410 : */
411 : void setupSparsity();
412 : /**
413 : * Check if the DOFs are covered by nodal BCs
414 : */
415 : void checkNodalBCs();
416 : /**
417 : * Get the list of DOFs covered by a nodal BC
418 : * @param nbc The Kokkos nodal BC object
419 : * @param dofs The flag whether each DOF is covered by the nodal BC
420 : */
421 : void getNodalBCDofs(const NodalBCBase * nbc, Array<bool> & dofs);
422 :
423 : /**
424 : * Kokkos thread object
425 : */
426 : Thread _thread;
427 : /**
428 : * Reference of the MOOSE system
429 : */
430 : SystemBase & _system;
431 : /**
432 : * Reference of the MOOSE mesh
433 : */
434 : const MooseMesh & _mesh;
435 : /**
436 : * Reference of the libMesh DOF map
437 : */
438 : const libMesh::DofMap & _dof_map;
439 : /**
440 : * Reference of the libMesh communicator
441 : */
442 : const Parallel::Communicator & _comm;
443 : /**
444 : * Number of variables
445 : */
446 : const unsigned int _num_vars;
447 : /**
448 : * Number of local DOFs
449 : */
450 : const dof_id_type _num_local_dofs;
451 : /**
452 : * Number of ghost DOFs
453 : */
454 : const dof_id_type _num_ghost_dofs;
455 :
456 : /**
457 : * Kokkos vectors and matrices on device
458 : */
459 : ///@{
460 : Array<Vector> _vectors;
461 : Array<Matrix> _matrices;
462 : ///@}
463 : /**
464 : * Cached elemental quadrature values and gradients
465 : */
466 : ///@{
467 : Array<Array2D<Array<Real>>> _qp_solutions;
468 : Array<Array2D<Array<Real3>>> _qp_solutions_grad;
469 : ///@}
470 : /**
471 : * Local DOF index of each variable
472 : */
473 : ///@{
474 : Array<Array2D<dof_id_type>> _local_elem_dof_index;
475 : Array<Array<dof_id_type>> _local_node_dof_index;
476 : ///@}
477 : /**
478 : * Map from local DOF index to global DOF index
479 : */
480 : Array<dof_id_type> _local_to_global_dof_index;
481 : /**
482 : * Maximum number of DOFs per element for each variable
483 : */
484 : Array<unsigned int> _max_dofs_per_elem;
485 : /**
486 : * FE type ID of each variable
487 : */
488 : Array<unsigned int> _var_fe_types;
489 : /**
490 : * Whether each variable is active on subdomains
491 : */
492 : Array2D<bool> _var_subdomain_active;
493 : /**
494 : * Off-diagonal coupled variable numbers of each variable
495 : */
496 : Array<Array<unsigned int>> _coupling;
497 :
498 : /**
499 : * List of active variable numbers
500 : */
501 : Array<unsigned int> _active_variables;
502 : /**
503 : * List of active tags
504 : */
505 : ///@{
506 : Array<TagID> _active_solution_tags;
507 : Array<TagID> _active_residual_tags;
508 : Array<TagID> _active_matrix_tags;
509 : ///@}
510 : /**
511 : * Flag whether each tag is active
512 : */
513 : ///@{
514 : Array<bool> _residual_tag_active;
515 : Array<bool> _matrix_tag_active;
516 : ///@}
517 : /**
518 : * Flag whether each DOF is covered by a nodal BC
519 : */
520 : Array<bool> _nbc_dof;
521 : /**
522 : * Flag whether each DOF is covered by a nodal BC for extra residual tags
523 : */
524 : ///@{
525 : Array<Array<bool>> _nbc_residual_tag_dof;
526 : Array<Array<bool>> _nbc_matrix_tag_dof;
527 : ///@}
528 : /**
529 : * List of DOFs to send and receive
530 : */
531 : ///@{
532 : Array<Array<dof_id_type>> _local_comm_list;
533 : Array<Array<dof_id_type>> _ghost_comm_list;
534 : ///@}
535 :
536 : /**
537 : * Matrix sparsity pattern data
538 : */
539 : Sparsity _sparsity;
540 : };
541 :
542 : #ifdef MOOSE_KOKKOS_SCOPE
543 : KOKKOS_FUNCTION inline Real
544 176196 : System::getVectorQpValueFace(
545 : ElementInfo info, unsigned int side, unsigned int qp, unsigned int var, TagID tag) const
546 : {
547 176196 : auto fe = _var_fe_types[var];
548 176196 : auto n_dofs = kokkosAssembly().getNumDofs(info.type, fe);
549 176196 : auto & phi = kokkosAssembly().getPhiFace(info.subdomain, info.type, fe)(side);
550 :
551 176196 : Real value = 0;
552 :
553 879148 : for (unsigned int i = 0; i < n_dofs; ++i)
554 702952 : value += getVectorDofValue(getElemLocalDofIndex(info.id, i, var), tag) * phi(i, qp);
555 :
556 176196 : return value;
557 : }
558 : KOKKOS_FUNCTION inline Real3
559 0 : System::getVectorQpGradFace(ElementInfo info,
560 : unsigned int side,
561 : Real33 jacobian,
562 : unsigned int qp,
563 : unsigned int var,
564 : TagID tag) const
565 : {
566 0 : auto fe = _var_fe_types[var];
567 0 : auto n_dofs = kokkosAssembly().getNumDofs(info.type, fe);
568 0 : auto & grad_phi = kokkosAssembly().getGradPhiFace(info.subdomain, info.type, fe)(side);
569 :
570 0 : Real3 grad = 0;
571 :
572 0 : for (unsigned int i = 0; i < n_dofs; ++i)
573 0 : grad += getVectorDofValue(getElemLocalDofIndex(info.id, i, var), tag) *
574 0 : (jacobian * grad_phi(i, qp));
575 :
576 0 : return grad;
577 : }
578 : #endif
579 :
580 : /**
581 : * The Kokkos interface that holds the host reference of the Kokkos systems and copies it to device
582 : * during parallel dispatch.
583 : * Maintains synchronization between host and device Kokkos systems and provides access to the
584 : * appropriate Kokkos systems depending on the architecture.
585 : */
586 : class SystemHolder
587 : {
588 : public:
589 : /**
590 : * Constructor
591 : * @param systems The Kokkos systems
592 : */
593 4571 : SystemHolder(Array<System> & systems) : _systems_host(systems), _systems_device(systems) {}
594 : /**
595 : * Copy constructor
596 : */
597 262930 : SystemHolder(const SystemHolder & holder)
598 222032 : : _systems_host(holder._systems_host), _systems_device(holder._systems_host)
599 : {
600 262930 : }
601 :
602 : #ifdef MOOSE_KOKKOS_SCOPE
603 : /**
604 : * Get the const reference of the Kokkos systems
605 : * @returns The const reference of the Kokkos systems depending on the architecture this function
606 : * is being called on
607 : */
608 8514019 : KOKKOS_FUNCTION const Array<System> & kokkosSystems() const
609 : {
610 8514019 : KOKKOS_IF_ON_HOST(return _systems_host;)
611 :
612 8514019 : return _systems_device;
613 : }
614 : /**
615 : * Get the writeable host reference of the Kokkos systems
616 : * @returns The writeable host reference of the Kokkos systems
617 : */
618 : Array<System> & kokkosSystems() { return _systems_host; }
619 : /**
620 : * Get the const reference of a Kokkos system
621 : * @param sys The system number
622 : * @returns The const reference of the Kokkos system depending on the architecture this function
623 : * is being called on
624 : */
625 29180100 : KOKKOS_FUNCTION const System & kokkosSystem(unsigned int sys) const
626 : {
627 29180100 : KOKKOS_IF_ON_HOST(return _systems_host[sys];)
628 :
629 29180100 : return _systems_device[sys];
630 : }
631 : /**
632 : * Get the writeable reference of a Kokkos system
633 : * @param sys The system number
634 : * @returns The writeable host reference of the Kokkos system
635 : */
636 3556 : System & kokkosSystem(unsigned int sys) { return _systems_host[sys]; }
637 : #endif
638 :
639 : private:
640 : /**
641 : * Host reference of the Kokkos systems
642 : */
643 : Array<System> & _systems_host;
644 : /**
645 : * Device copy of the Kokkos systems
646 : */
647 : const Array<System> _systems_device;
648 : };
649 :
650 : } // namespace Kokkos
651 : } // namespace Moose
|