milupHPC documentation
  • include
  • subdomain_key_tree
subdomain.cuh
Go to the documentation of this file.
1
10#ifndef MILUPHPC_DOMAIN_CUH
11#define MILUPHPC_DOMAIN_CUH
12
13#include "../parameter.h"
14#include "../cuda_utils/cuda_runtime.h"
15#include "../cuda_utils/cuda_utilities.cuh"
16#include "../helper_handler.h"
17#include "tree.cuh"
18
19//TODO: necessary/reasonable to have Key class?
20//class Key {
21 //keyType *keys;
22 //integer *maxLevel;
23 //CUDA_CALLABLE_MEMBER Key();
24 //CUDA_CALLABLE_MEMBER Key(keyType *keys, integer *maxLevel);
25 //CUDA_CALLABLE_MEMBER ~Key();
26 //CUDA_CALLABLE_MEMBER void set(keyType *keys, integer *maxLevel);*/
27//};
28
29// Forward declaration of DomainList class
30class DomainList;
31// Forward declaration of SubDomainKeyTree class
32class SubDomainKeyTree;
33
35namespace KeyNS {
36
44 CUDA_CALLABLE_MEMBER void key2Char(keyType key, integer maxLevel, char *keyAsChar);
45
56 CUDA_CALLABLE_MEMBER integer key2proc(keyType key, SubDomainKeyTree *subDomainKeyTree/*, Curve::Type curveType=Curve::lebesgue*/);
57}
58
62class SubDomainKeyTree {
63
64public:
66 integer rank;
68 integer numProcesses;
70 keyType *range;
71 //keyType *hilbertRange;
73 integer *procParticleCounter;
74
78 CUDA_CALLABLE_MEMBER SubDomainKeyTree();
87 CUDA_CALLABLE_MEMBER SubDomainKeyTree(integer rank, integer numProcesses, keyType *range,
88 integer *procParticleCounter);
92 CUDA_CALLABLE_MEMBER ~SubDomainKeyTree();
93
102 CUDA_CALLABLE_MEMBER void set(integer rank, integer numProcesses, keyType *range, integer *procParticleCounter);
103
110 CUDA_CALLABLE_MEMBER integer key2proc(keyType key/*, Curve::Type curveType=Curve::lebesgue*/);
111
121 CUDA_CALLABLE_MEMBER bool isDomainListNode(keyType key, integer maxLevel, integer level,
122 Curve::Type curveType=Curve::lebesgue);
123
124};
125
127namespace SubDomainKeyTreeNS {
128
130 namespace Kernel {
131
143 __global__ void set(SubDomainKeyTree *subDomainKeyTree, integer rank, integer numProcesses, keyType *range,
144 integer *procParticleCounter);
145
153 __global__ void test(SubDomainKeyTree *subDomainKeyTree);
154
180 __global__ void buildDomainTree(Tree *tree, Particles *particles, DomainList *domainList, integer n, integer m);
181
182 __global__ void buildDomainTree(SubDomainKeyTree *subDomainKeyTree, Tree *tree, Particles *particles, DomainList *domainList, integer n, integer m, integer level);
183
197 __global__ void getParticleKeys(SubDomainKeyTree *subDomainKeyTree, Tree *tree,
198 Particles *particles, keyType *keys, integer maxLevel,
199 integer n, Curve::Type curveType = Curve::lebesgue);
200
213 __global__ void particlesPerProcess(SubDomainKeyTree *subDomainKeyTree, Tree *tree, Particles *particles,
214 integer n, integer m, Curve::Type curveType=Curve::lebesgue);
215
233 __global__ void markParticlesProcess(SubDomainKeyTree *subDomainKeyTree, Tree *tree, Particles *particles,
234 integer n, integer m, integer *sortArray,
235 Curve::Type curveType=Curve::lebesgue);
236
246 __global__ void zeroDomainListNodes(Particles *particles, DomainList *domainList,
247 DomainList *lowestDomainList);
248
259 template <typename T>
260 __global__ void prepareLowestDomainExchange(Particles *particles, DomainList *lowestDomainList,
261 T *buffer, Entry::Name entry);
262
274 template <typename T>
275 __global__ void updateLowestDomainListNodes(Particles *particles, DomainList *lowestDomainList,
276 T *buffer, Entry::Name entry);
277
287 __global__ void compLowestDomainListNodes(Tree *tree, Particles *particles, DomainList *lowestDomainList);
288
299 __global__ void compLocalPseudoParticles(Tree *tree, Particles *particles, DomainList *domainList, int n);
300
313 __global__ void compDomainListPseudoParticlesPerLevel(Tree *tree, Particles *particles, DomainList *domainList,
314 DomainList *lowestDomainList, int n, int level);
315
316 __global__ void compDomainListPseudoParticles(Tree *tree, Particles *particles, DomainList *domainList,
317 DomainList *lowestDomainList, int n);
318
340 __global__ void repairTree(SubDomainKeyTree *subDomainKeyTree, Tree *tree, Particles *particles,
341 DomainList *domainList, DomainList *lowestDomainList,
342 int n, int m, Curve::Type curveType);
343
344 // notes:
345 // - using Helper::keyTypeBuffer as keyHistRanges
346 __global__ void createKeyHistRanges(Helper *helper, integer bins);
347
348 // notes:
349 // - using Helper::keyTypeBuffer as keyHistRanges
350 // - using Helper::integerBuffer as keyHistCounts
351 __global__ void keyHistCounter(Tree *tree, Particles *particles, SubDomainKeyTree *subDomainKeyTree,
352 Helper *helper,
353 /*keyType *keyHistRanges, integer *keyHistCounts,*/ int bins, int n,
354 Curve::Type curveType=Curve::lebesgue);
355
356 // notes:
357 // - using Helper::keyTypeBuffer as keyHistRanges
358 // - using Helper::integerBuffer as keyHistCounts
359 __global__ void calculateNewRange(SubDomainKeyTree *subDomainKeyTree, Helper *helper,
360 /*keyType *keyHistRanges, integer *keyHistCounts,*/ int bins, int n,
361 Curve::Type curveType=Curve::lebesgue);
362
363
365 namespace Launch {
366
370 void set(SubDomainKeyTree *subDomainKeyTree, integer rank, integer numProcesses, keyType *range,
371 integer *procParticleCounter);
372
376 void test(SubDomainKeyTree *subDomainKeyTree);
377
383 real buildDomainTree(Tree *tree, Particles *particles, DomainList *domainList, integer n, integer m);
384
385 real buildDomainTree(SubDomainKeyTree *subDomainKeyTree, Tree *tree, Particles *particles, DomainList *domainList, integer n, integer m, integer level);
386
392 real getParticleKeys(SubDomainKeyTree *subDomainKeyTree, Tree *tree,
393 Particles *particles, keyType *keys, integer maxLevel,
394 integer n, Curve::Type curveType = Curve::lebesgue);
395
401 real particlesPerProcess(SubDomainKeyTree *subDomainKeyTree, Tree *tree, Particles *particles,
402 integer n, integer m, Curve::Type curveType=Curve::lebesgue);
403
409 real markParticlesProcess(SubDomainKeyTree *subDomainKeyTree, Tree *tree, Particles *particles,
410 integer n, integer m, integer *sortArray,
411 Curve::Type curveType=Curve::lebesgue);
412
418 real zeroDomainListNodes(Particles *particles, DomainList *domainList, DomainList *lowestDomainList);
419
425 template <typename T>
426 real prepareLowestDomainExchange(Particles *particles, DomainList *lowestDomainList,
427 T *buffer, Entry::Name entry);
428
429
435 template <typename T>
436 real updateLowestDomainListNodes(Particles *particles, DomainList *lowestDomainList,
437 T *buffer, Entry::Name entry);
438
444 real compLowestDomainListNodes(Tree *tree, Particles *particles, DomainList *lowestDomainList);
445
451 real compLocalPseudoParticles(Tree *tree, Particles *particles, DomainList *domainList, int n);
452
458 real compDomainListPseudoParticlesPerLevel(Tree *tree, Particles *particles, DomainList *domainList,
459 DomainList *lowestDomainList, int n, int level);
460
466 real compDomainListPseudoParticles(Tree *tree, Particles *particles, DomainList *domainList,
467 DomainList *lowestDomainList, int n);
468
474 real repairTree(SubDomainKeyTree *subDomainKeyTree, Tree *tree, Particles *particles,
475 DomainList *domainList, DomainList *lowestDomainList,
476 int n, int m, Curve::Type curveType);
477
478 real createKeyHistRanges(Helper *helper, integer bins);
479
480 real keyHistCounter(Tree *tree, Particles *particles, SubDomainKeyTree *subDomainKeyTree,
481 Helper *helper, int bins, int n, Curve::Type curveType=Curve::lebesgue);
482
483 real calculateNewRange(SubDomainKeyTree *subDomainKeyTree, Helper *helper, int bins, int n,
484 Curve::Type curveType=Curve::lebesgue);
485 }
486
487 }
488
489}
490
494class DomainList {
495
496public:
497
499 integer *domainListIndices;
501 integer *domainListLevels;
503 integer *domainListIndex;
505 integer *domainListCounter;
507 keyType *domainListKeys;
509 keyType *sortedDomainListKeys;
511 integer *relevantDomainListIndices;
513 integer *relevantDomainListLevels;
515 integer *relevantDomainListProcess;
516
517 integer *relevantDomainListOriginalIndex;
518
519 real *borders;
520
524 CUDA_CALLABLE_MEMBER DomainList();
525
537 CUDA_CALLABLE_MEMBER DomainList(integer *domainListIndices, integer *domainListLevels, integer *domainListIndex,
538 integer *domainListCounter, keyType *domainListKeys, keyType *sortedDomainListKeys,
539 integer *relevantDomainListIndices, integer *relevantDomainListLevels,
540 integer *relevantDomainListProcess);
552 CUDA_CALLABLE_MEMBER void set(integer *domainListIndices, integer *domainListLevels, integer *domainListIndex,
553 integer *domainListCounter, keyType *domainListKeys, keyType *sortedDomainListKeys,
554 integer *relevantDomainListIndices, integer *relevantDomainListLevels,
555 integer *relevantDomainListProcess);
556
557 CUDA_CALLABLE_MEMBER void setBorders(real *borders, integer *relevantDomainListOriginalIndex);
558
562 CUDA_CALLABLE_MEMBER ~DomainList();
563};
564
565namespace DomainListNS {
566
567 namespace Kernel {
568
583 __global__ void set(DomainList *domainList, integer *domainListIndices, integer *domainListLevels,
584 integer *domainListIndex, integer *domainListCounter, keyType *domainListKeys,
585 keyType *sortedDomainListKeys, integer *relevantDomainListIndices,
586 integer *relevantDomainListLevels, integer *relevantDomainListProcess);
587
588 __global__ void setBorders(DomainList *domainList, real *borders, integer *relevantDomainListOriginalIndex);
589
598 __global__ void info(Particles *particles, DomainList *domainList);
599
609 __global__ void info(Particles *particles, DomainList *domainList, DomainList *lowestDomainList);
610
626 __global__ void createDomainList(SubDomainKeyTree *subDomainKeyTree, DomainList *domainList,
627 integer maxLevel, Curve::Type curveType = Curve::lebesgue);
628
645 __global__ void lowestDomainList(SubDomainKeyTree *subDomainKeyTree, Tree *tree, Particles *particles,
646 DomainList *domainList, DomainList *lowestDomainList, integer n, integer m);
647
648
649 namespace Launch {
650
654 void set(DomainList *domainList, integer *domainListIndices, integer *domainListLevels,
655 integer *domainListIndex, integer *domainListCounter, keyType *domainListKeys,
656 keyType *sortedDomainListKeys, integer *relevantDomainListIndices,
657 integer *relevantDomainListLevels, integer *relevantDomainListProcess);
658
659 void setBorders(DomainList *domainList, real *borders, integer *relevantDomainListOriginalIndex);
660
666 real info(Particles *particles, DomainList *domainList);
667
668 real info(Particles *particles, DomainList *domainList, DomainList *lowestDomainList);
669
675 real createDomainList(SubDomainKeyTree *subDomainKeyTree, DomainList *domainList,
676 integer maxLevel, Curve::Type curveType = Curve::lebesgue);
677
683 real lowestDomainList(SubDomainKeyTree *subDomainKeyTree, Tree *tree, Particles *particles,
684 DomainList *domainList, DomainList *lowestDomainList, integer n, integer m);
685
686 }
687 }
688
689}
690
692namespace ParticlesNS {
693
704 __device__ bool applySphericalCriterion(SubDomainKeyTree *subDomainKeyTree, Tree *tree, Particles *particles,
705 real d, int index);
706
717 __device__ bool applyCubicCriterion(SubDomainKeyTree *subDomainKeyTree, Tree *tree, Particles *particles,
718 real d, int index);
719
720 namespace Kernel {
721 __global__ void mark2remove(SubDomainKeyTree *subDomainKeyTree, Tree *tree, Particles *particles,
722 int *particles2remove, int *counter, int criterion, real d,
723 int numParticles);
724
725 namespace Launch {
726 real mark2remove(SubDomainKeyTree *subDomainKeyTree, Tree *tree, Particles *particles,
727 int *particles2remove, int *counter, int criterion, real d,
728 int numParticles);
729 }
730 }
731}
732
733namespace CudaUtils {
734 namespace Kernel {
735 template<typename T, typename U>
736 __global__ void markDuplicatesTemp(Tree *tree, DomainList *domainList, T *array, U *entry1, U *entry2, U *entry3, integer *duplicateCounter, integer *child, int length);
737
738 template <typename T, unsigned int blockSize>
739 __global__ void reduceBlockwise(T *array, T *outputData, int n);
740
741 template <typename T, unsigned int blockSize>
742 __global__ void blockReduction(const T *indata, T *outdata);
743
744 namespace Launch {
745 template<typename T, typename U>
746 real markDuplicatesTemp(Tree *tree, DomainList *domainList, T *array, U *entry1, U *entry2, U *entry3, integer *duplicateCounter, integer *child, int length);
747
748 template <typename T, unsigned int blockSize>
749 real reduceBlockwise(T *array, T *outputData, int n);
750
751 template <typename T, unsigned int blockSize>
752 real blockReduction(const T *indata, T *outdata);
753 }
754 }
755}
756
758namespace Physics {
759 namespace Kernel {
760
771 template <unsigned int blockSize>
772 __global__ void calculateAngularMomentumBlockwise(Particles *particles, real *outputData, int n);
773
783 template <unsigned int blockSize>
784 __global__ void sumAngularMomentum(const real *indata, real *outdata);
785
794 __global__ void kineticEnergy(Particles *particles, int n);
795
796 namespace Launch {
802 template <unsigned int blockSize>
803 real calculateAngularMomentumBlockwise(Particles *particles, real *outputData, int n);
804
810 template <unsigned int blockSize>
811 real sumAngularMomentum(const real *indata, real *outdata);
812
818 real kineticEnergy(Particles *particles, int n);
819 }
820 }
821}
822
823#endif //MILUPHPC_DOMAIN_CUH
DomainList
Definition: subdomain.cuh:494
DomainList::domainListKeys
keyType * domainListKeys
domain list node keys
Definition: subdomain.cuh:507
DomainList::domainListIndex
integer * domainListIndex
domain list node index, thus amount of domain list nodes
Definition: subdomain.cuh:503
DomainList::setBorders
CUDA_CALLABLE_MEMBER void setBorders(real *borders, integer *relevantDomainListOriginalIndex)
Definition: subdomain.cu:1714
DomainList::sortedDomainListKeys
keyType * sortedDomainListKeys
sorted domain list node keys, usable as output for sorting the keys
Definition: subdomain.cuh:509
DomainList::DomainList
CUDA_CALLABLE_MEMBER DomainList()
Constructor.
Definition: subdomain.cu:1674
DomainList::relevantDomainListOriginalIndex
integer * relevantDomainListOriginalIndex
Definition: subdomain.cuh:517
DomainList::domainListIndices
integer * domainListIndices
domain list node indices
Definition: subdomain.cuh:499
DomainList::domainListLevels
integer * domainListLevels
domain list node levels
Definition: subdomain.cuh:501
DomainList::relevantDomainListIndices
integer * relevantDomainListIndices
concentrate domain list nodes, usable to reduce domain list indices in respect to some criterion
Definition: subdomain.cuh:511
DomainList::domainListCounter
integer * domainListCounter
domain list node counter, usable as buffer
Definition: subdomain.cuh:505
DomainList::relevantDomainListProcess
integer * relevantDomainListProcess
Definition: subdomain.cuh:515
DomainList::set
CUDA_CALLABLE_MEMBER void set(integer *domainListIndices, integer *domainListLevels, integer *domainListIndex, integer *domainListCounter, keyType *domainListKeys, keyType *sortedDomainListKeys, integer *relevantDomainListIndices, integer *relevantDomainListLevels, integer *relevantDomainListProcess)
Setter, passing pointer to member variables.
Definition: subdomain.cu:1696
DomainList::borders
real * borders
Definition: subdomain.cuh:519
DomainList::~DomainList
CUDA_CALLABLE_MEMBER ~DomainList()
Destructor.
Definition: subdomain.cu:1692
DomainList::relevantDomainListLevels
integer * relevantDomainListLevels
Definition: subdomain.cuh:513
Helper
Definition: helper.cuh:24
Particles
Particle(s) class based on SoA (Structur of Arrays).
Definition: particles.cuh:50
SubDomainKeyTree
SubDomainKeyTree class handling rank, number of processes and ranges.
Definition: subdomain.cuh:62
SubDomainKeyTree::range
keyType * range
Space-filling curve ranges, mapping key ranges/borders to MPI processes.
Definition: subdomain.cuh:70
SubDomainKeyTree::~SubDomainKeyTree
CUDA_CALLABLE_MEMBER ~SubDomainKeyTree()
Destructor.
Definition: subdomain.cu:32
SubDomainKeyTree::key2proc
CUDA_CALLABLE_MEMBER integer key2proc(keyType key)
Compute particle's MPI process belonging by it's key.
Definition: subdomain.cu:68
SubDomainKeyTree::isDomainListNode
CUDA_CALLABLE_MEMBER bool isDomainListNode(keyType key, integer maxLevel, integer level, Curve::Type curveType=Curve::lebesgue)
Check whether key, thus particle, represents a domain list node.
Definition: subdomain.cu:79
SubDomainKeyTree::set
CUDA_CALLABLE_MEMBER void set(integer rank, integer numProcesses, keyType *range, integer *procParticleCounter)
Setter.
Definition: subdomain.cu:36
SubDomainKeyTree::procParticleCounter
integer * procParticleCounter
particle counter in dependence of MPI process(es)
Definition: subdomain.cuh:73
SubDomainKeyTree::rank
integer rank
MPI rank.
Definition: subdomain.cuh:66
SubDomainKeyTree::SubDomainKeyTree
CUDA_CALLABLE_MEMBER SubDomainKeyTree()
Default Constructor.
Definition: subdomain.cu:21
SubDomainKeyTree::numProcesses
integer numProcesses
MPI number of processes.
Definition: subdomain.cuh:68
Tree
Tree class.
Definition: tree.cuh:140
CUDA_CALLABLE_MEMBER
#define CUDA_CALLABLE_MEMBER
Definition: cuda_utilities.cuh:30
CudaUtils::Kernel::Launch::markDuplicatesTemp
real markDuplicatesTemp(Tree *tree, DomainList *domainList, T *array, U *entry1, U *entry2, U *entry3, integer *duplicateCounter, integer *child, int length)
Definition: subdomain.cu:2397
CudaUtils::Kernel::Launch::reduceBlockwise
real reduceBlockwise(T *array, T *outputData, int n)
Definition: subdomain.cu:2405
CudaUtils::Kernel::Launch::blockReduction
real blockReduction(const T *indata, T *outdata)
Definition: subdomain.cu:2412
CudaUtils::Kernel::markDuplicatesTemp
__global__ void markDuplicatesTemp(Tree *tree, DomainList *domainList, T *array, U *entry1, U *entry2, U *entry3, integer *duplicateCounter, integer *child, int length)
Definition: subdomain.cu:2249
CudaUtils::Kernel::reduceBlockwise
__global__ void reduceBlockwise(T *array, T *outputData, int n)
Definition: subdomain.cu:2325
CudaUtils::Kernel::blockReduction
__global__ void blockReduction(const T *indata, T *outdata)
Definition: subdomain.cu:2380
CudaUtils
Definition: cuda_utilities.cuh:100
DomainListNS::Kernel::Launch::info
real info(Particles *particles, DomainList *domainList)
Wrapper for DomainListNS::Kernel::info().
Definition: subdomain.cu:2131
DomainListNS::Kernel::Launch::set
void set(DomainList *domainList, integer *domainListIndices, integer *domainListLevels, integer *domainListIndex, integer *domainListCounter, keyType *domainListKeys, keyType *sortedDomainListKeys, integer *relevantDomainListIndices, integer *relevantDomainListLevels, integer *relevantDomainListProcess)
Wrapper for DomainListNS::Kernel::set().
Definition: subdomain.cu:2115
DomainListNS::Kernel::Launch::createDomainList
real createDomainList(SubDomainKeyTree *subDomainKeyTree, DomainList *domainList, integer maxLevel, Curve::Type curveType=Curve::lebesgue)
Wrapper for DomainListNS::Kernel::createDomainList().
Definition: subdomain.cu:2141
DomainListNS::Kernel::Launch::lowestDomainList
real lowestDomainList(SubDomainKeyTree *subDomainKeyTree, Tree *tree, Particles *particles, DomainList *domainList, DomainList *lowestDomainList, integer n, integer m)
Wrapper for ::DomainListNS::Kernel::lowestDoainList().
Definition: subdomain.cu:2149
DomainListNS::Kernel::Launch::setBorders
void setBorders(DomainList *domainList, real *borders, integer *relevantDomainListOriginalIndex)
Definition: subdomain.cu:2125
DomainListNS::Kernel::lowestDomainList
__global__ void lowestDomainList(SubDomainKeyTree *subDomainKeyTree, Tree *tree, Particles *particles, DomainList *domainList, DomainList *lowestDomainList, integer n, integer m)
Kernel to create the lowest domain list.
Definition: subdomain.cu:2039
DomainListNS::Kernel::set
__global__ void set(DomainList *domainList, integer *domainListIndices, integer *domainListLevels, integer *domainListIndex, integer *domainListCounter, keyType *domainListKeys, keyType *sortedDomainListKeys, integer *relevantDomainListIndices, integer *relevantDomainListLevels, integer *relevantDomainListProcess)
Kernel call to setter.
Definition: subdomain.cu:1723
DomainListNS::Kernel::setBorders
__global__ void setBorders(DomainList *domainList, real *borders, integer *relevantDomainListOriginalIndex)
Definition: subdomain.cu:1733
DomainListNS::Kernel::createDomainList
__global__ void createDomainList(SubDomainKeyTree *subDomainKeyTree, DomainList *domainList, integer maxLevel, Curve::Type curveType=Curve::lebesgue)
Kernel to create the domain list.
Definition: subdomain.cu:1962
DomainListNS::Kernel::info
__global__ void info(Particles *particles, DomainList *domainList)
Info kernel (for debugging purposes).
Definition: subdomain.cu:1737
DomainListNS
Definition: subdomain.cuh:565
HelperNS::sortArray
real sortArray(A *arrayToSort, A *sortedArray, B *keyIn, B *keyOut, integer n)
Definition: helper.cu:151
Kernel
Definition: device_rhs.cuh:7
KeyNS
Key (keyType) related functions and kernels.
Definition: subdomain.cuh:35
KeyNS::key2proc
CUDA_CALLABLE_MEMBER integer key2proc(keyType key, SubDomainKeyTree *subDomainKeyTree)
Convert the key to the corresponding process.
Definition: subdomain.cu:17
KeyNS::key2Char
CUDA_CALLABLE_MEMBER void key2Char(keyType key, integer maxLevel, char *keyAsChar)
Convert a key to a char for printing.
Definition: subdomain.cu:5
ParticlesNS::Kernel::Launch::mark2remove
real mark2remove(SubDomainKeyTree *subDomainKeyTree, Tree *tree, Particles *particles, int *particles2remove, int *counter, int criterion, real d, int numParticles)
Definition: subdomain.cu:2233
ParticlesNS::Kernel::mark2remove
__global__ void mark2remove(SubDomainKeyTree *subDomainKeyTree, Tree *tree, Particles *particles, int *particles2remove, int *counter, int criterion, real d, int numParticles)
Definition: subdomain.cu:2200
ParticlesNS
Particle class related functions and kernels.
Definition: particles.cuh:552
ParticlesNS::applySphericalCriterion
__device__ bool applySphericalCriterion(SubDomainKeyTree *subDomainKeyTree, Tree *tree, Particles *particles, real d, int index)
Check whether particle(s) are within sphere from simulation center.
Definition: subdomain.cu:2161
ParticlesNS::applyCubicCriterion
__device__ bool applyCubicCriterion(SubDomainKeyTree *subDomainKeyTree, Tree *tree, Particles *particles, real d, int index)
Check whether particle(s) are within cube from simulation center.
Definition: subdomain.cu:2179
Physics::Kernel::Launch::kineticEnergy
real kineticEnergy(Particles *particles, int n)
Wrapper for: Physics::Kernel::kineticEnergy().
Definition: subdomain.cu:2638
Physics::Kernel::Launch::sumAngularMomentum
real sumAngularMomentum(const real *indata, real *outdata)
Wrapper for: Physics::Kernel::sumAngularMomentum().
Definition: subdomain.cu:2631
Physics::Kernel::Launch::calculateAngularMomentumBlockwise
real calculateAngularMomentumBlockwise(Particles *particles, real *outputData, int n)
Wrapper for: Physics::Kernel::calculateAngularMomentumBlockwise().
Definition: subdomain.cu:2623
Physics::Kernel::calculateAngularMomentumBlockwise
__global__ void calculateAngularMomentumBlockwise(Particles *particles, real *outputData, int n)
Calculate angular momentum for all particles (per block).
Definition: subdomain.cu:2428
Physics::Kernel::kineticEnergy
__global__ void kineticEnergy(Particles *particles, int n)
Calculate kinetic energy.
Definition: subdomain.cu:2598
Physics::Kernel::sumAngularMomentum
__global__ void sumAngularMomentum(const real *indata, real *outdata)
Calculate angular momentum: sum over blocks.
Definition: subdomain.cu:2576
Physics
Physics related functions and kernels.
Definition: subdomain.cuh:758
ProfilerIds::Time::tree
const char *const tree
Definition: h5profiler.h:57
ProfilerIds::numParticles
const char *const numParticles
Definition: h5profiler.h:29
SubDomainKeyTreeNS::Kernel::Launch::compLocalPseudoParticles
real compLocalPseudoParticles(Tree *tree, Particles *particles, DomainList *domainList, int n)
Wrapper for SubDomainKeyTreeNS::Kernel::compLocalPseudoParticles().
Definition: subdomain.cu:1622
SubDomainKeyTreeNS::Kernel::Launch::particlesPerProcess
real particlesPerProcess(SubDomainKeyTree *subDomainKeyTree, Tree *tree, Particles *particles, integer n, integer m, Curve::Type curveType=Curve::lebesgue)
Wrapper for SubDomainKeyTreeNS::Kernel::particlesPerProcess().
Definition: subdomain.cu:1572
SubDomainKeyTreeNS::Kernel::Launch::updateLowestDomainListNodes
real updateLowestDomainListNodes(Particles *particles, DomainList *lowestDomainList, T *buffer, Entry::Name entry)
Wrapper for SubDomainKeyTreeNS::Kernel::updateLowestDomainListNodes().
Definition: subdomain.cu:1605
SubDomainKeyTreeNS::Kernel::Launch::compDomainListPseudoParticlesPerLevel
real compDomainListPseudoParticlesPerLevel(Tree *tree, Particles *particles, DomainList *domainList, DomainList *lowestDomainList, int n, int level)
Wrapper for SubDomainKeyTreeNS::Kernel::compDomainListPseudoParticlesPerLevel().
Definition: subdomain.cu:1628
SubDomainKeyTreeNS::Kernel::Launch::compDomainListPseudoParticles
real compDomainListPseudoParticles(Tree *tree, Particles *particles, DomainList *domainList, DomainList *lowestDomainList, int n)
Wrapper for SubDomainKeyTreeNS::Kernel::compDomainListPseudoParticles().
Definition: subdomain.cu:1635
SubDomainKeyTreeNS::Kernel::Launch::prepareLowestDomainExchange
real prepareLowestDomainExchange(Particles *particles, DomainList *lowestDomainList, T *buffer, Entry::Name entry)
Wrapper for SubDomainKeyTreeNS::Kernel::prepareLowestDomainExchange().
Definition: subdomain.cu:1594
SubDomainKeyTreeNS::Kernel::Launch::zeroDomainListNodes
real zeroDomainListNodes(Particles *particles, DomainList *domainList, DomainList *lowestDomainList)
Wrapper for SubDomainKeyTreeNS::Kernel::zeroDomainListNodes().
Definition: subdomain.cu:1587
SubDomainKeyTreeNS::Kernel::Launch::buildDomainTree
real buildDomainTree(Tree *tree, Particles *particles, DomainList *domainList, integer n, integer m)
Wrapper for SubDomainKeyTreeNS::Kernel::buildDomainTree().
Definition: subdomain.cu:1550
SubDomainKeyTreeNS::Kernel::Launch::createKeyHistRanges
real createKeyHistRanges(Helper *helper, integer bins)
Definition: subdomain.cu:1650
SubDomainKeyTreeNS::Kernel::Launch::repairTree
real repairTree(SubDomainKeyTree *subDomainKeyTree, Tree *tree, Particles *particles, DomainList *domainList, DomainList *lowestDomainList, int n, int m, Curve::Type curveType)
Wrapper for SubDomainKeyTreeNS::Kernel::repairTree().
Definition: subdomain.cu:1642
SubDomainKeyTreeNS::Kernel::Launch::getParticleKeys
real getParticleKeys(SubDomainKeyTree *subDomainKeyTree, Tree *tree, Particles *particles, keyType *keys, integer maxLevel, integer n, Curve::Type curveType=Curve::lebesgue)
Wrapper for SubDomainKeyTreeNS::Kernel::getParticleKeys().
Definition: subdomain.cu:1564
SubDomainKeyTreeNS::Kernel::Launch::keyHistCounter
real keyHistCounter(Tree *tree, Particles *particles, SubDomainKeyTree *subDomainKeyTree, Helper *helper, int bins, int n, Curve::Type curveType=Curve::lebesgue)
Definition: subdomain.cu:1655
SubDomainKeyTreeNS::Kernel::Launch::markParticlesProcess
real markParticlesProcess(SubDomainKeyTree *subDomainKeyTree, Tree *tree, Particles *particles, integer n, integer m, integer *sortArray, Curve::Type curveType=Curve::lebesgue)
Wrapper for ::SubDomainKeyTreeNS::Kernel::markParticlesPerProcess().
Definition: subdomain.cu:1579
SubDomainKeyTreeNS::Kernel::Launch::compLowestDomainListNodes
real compLowestDomainListNodes(Tree *tree, Particles *particles, DomainList *lowestDomainList)
Wrapper for SubDomainKeyTreeNS::Kernel::compLowestDomainListNodes().
Definition: subdomain.cu:1616
SubDomainKeyTreeNS::Kernel::Launch::set
void set(SubDomainKeyTree *subDomainKeyTree, integer rank, integer numProcesses, keyType *range, integer *procParticleCounter)
Wrapper for SubDomainKeyTreeNS::Kernel::set().
Definition: subdomain.cu:1538
SubDomainKeyTreeNS::Kernel::Launch::test
void test(SubDomainKeyTree *subDomainKeyTree)
Wrapper for SubDomainKeyTreeNS::Kernel::test().
Definition: subdomain.cu:1545
SubDomainKeyTreeNS::Kernel::Launch::calculateNewRange
real calculateNewRange(SubDomainKeyTree *subDomainKeyTree, Helper *helper, int bins, int n, Curve::Type curveType=Curve::lebesgue)
Definition: subdomain.cu:1662
SubDomainKeyTreeNS::Kernel::markParticlesProcess
__global__ void markParticlesProcess(SubDomainKeyTree *subDomainKeyTree, Tree *tree, Particles *particles, integer n, integer m, integer *sortArray, Curve::Type curveType=Curve::lebesgue)
Kernel to mark particle's belonging.
Definition: subdomain.cu:921
SubDomainKeyTreeNS::Kernel::test
__global__ void test(SubDomainKeyTree *subDomainKeyTree)
Test kernel call (for debugging/testing purposes).
Definition: subdomain.cu:142
SubDomainKeyTreeNS::Kernel::compDomainListPseudoParticles
__global__ void compDomainListPseudoParticles(Tree *tree, Particles *particles, DomainList *domainList, DomainList *lowestDomainList, int n)
Definition: subdomain.cu:1291
SubDomainKeyTreeNS::Kernel::getParticleKeys
__global__ void getParticleKeys(SubDomainKeyTree *subDomainKeyTree, Tree *tree, Particles *particles, keyType *keys, integer maxLevel, integer n, Curve::Type curveType=Curve::lebesgue)
Kernel to get all particle keys (and additional information for debugging purposes).
Definition: subdomain.cu:849
SubDomainKeyTreeNS::Kernel::repairTree
__global__ void repairTree(SubDomainKeyTree *subDomainKeyTree, Tree *tree, Particles *particles, DomainList *domainList, DomainList *lowestDomainList, int n, int m, Curve::Type curveType)
Repair tree by removing received and inserted (pseudo-)particles.
Definition: subdomain.cu:1355
SubDomainKeyTreeNS::Kernel::set
__global__ void set(SubDomainKeyTree *subDomainKeyTree, integer rank, integer numProcesses, keyType *range, integer *procParticleCounter)
Kernel call to setter.
Definition: subdomain.cu:137
SubDomainKeyTreeNS::Kernel::compLowestDomainListNodes
__global__ void compLowestDomainListNodes(Tree *tree, Particles *particles, DomainList *lowestDomainList)
Compute/Find lowest domain list nodes.
Definition: subdomain.cu:1087
SubDomainKeyTreeNS::Kernel::prepareLowestDomainExchange
__global__ void prepareLowestDomainExchange(Particles *particles, DomainList *lowestDomainList, T *buffer, Entry::Name entry)
Prepare lowest domain exchange via MPI by copying to contiguous memory.
Definition: subdomain.cu:995
SubDomainKeyTreeNS::Kernel::updateLowestDomainListNodes
__global__ void updateLowestDomainListNodes(Particles *particles, DomainList *lowestDomainList, T *buffer, Entry::Name entry)
Update lowest domain list nodes.
Definition: subdomain.cu:1034
SubDomainKeyTreeNS::Kernel::particlesPerProcess
__global__ void particlesPerProcess(SubDomainKeyTree *subDomainKeyTree, Tree *tree, Particles *particles, integer n, integer m, Curve::Type curveType=Curve::lebesgue)
Kernel to check particle's belonging and count in dependence of belonging/process.
Definition: subdomain.cu:897
SubDomainKeyTreeNS::Kernel::keyHistCounter
__global__ void keyHistCounter(Tree *tree, Particles *particles, SubDomainKeyTree *subDomainKeyTree, Helper *helper, int bins, int n, Curve::Type curveType=Curve::lebesgue)
Definition: subdomain.cu:1479
SubDomainKeyTreeNS::Kernel::compLocalPseudoParticles
__global__ void compLocalPseudoParticles(Tree *tree, Particles *particles, DomainList *domainList, int n)
Compute local (tree) pseudo particles.
Definition: subdomain.cu:1139
SubDomainKeyTreeNS::Kernel::calculateNewRange
__global__ void calculateNewRange(SubDomainKeyTree *subDomainKeyTree, Helper *helper, int bins, int n, Curve::Type curveType=Curve::lebesgue)
Definition: subdomain.cu:1507
SubDomainKeyTreeNS::Kernel::compDomainListPseudoParticlesPerLevel
__global__ void compDomainListPseudoParticlesPerLevel(Tree *tree, Particles *particles, DomainList *domainList, DomainList *lowestDomainList, int n, int level)
Compute domain list pseudo particles (per level).
Definition: subdomain.cu:1175
SubDomainKeyTreeNS::Kernel::buildDomainTree
__global__ void buildDomainTree(Tree *tree, Particles *particles, DomainList *domainList, integer n, integer m)
Kernel to build the domain tree.
Definition: subdomain.cu:148
SubDomainKeyTreeNS::Kernel::zeroDomainListNodes
__global__ void zeroDomainListNodes(Particles *particles, DomainList *domainList, DomainList *lowestDomainList)
Zero domain list nodes.
Definition: subdomain.cu:946
SubDomainKeyTreeNS::Kernel::createKeyHistRanges
__global__ void createKeyHistRanges(Helper *helper, integer bins)
Definition: subdomain.cu:1459
SubDomainKeyTreeNS
SubDomainKeyTree related functions and kernels.
Definition: subdomain.cuh:127
real
double real
Definition: parameter.h:15
keyType
unsigned long keyType
Definition: parameter.h:18
integer
int integer
Definition: parameter.h:17
Curve::Type
Type
Definition: parameter.h:206
Curve::lebesgue
@ lebesgue
Definition: parameter.h:207
Entry::Name
Name
Definition: parameter.h:255
tree.cuh
Tree related classes, kernels and functions.

milupHPC - include/subdomain_key_tree/subdomain.cuh Source File
Generated on Wed Aug 31 2022 12:16:52 by Doxygen 1.9.3