Skip to content
GitLab
Explore
Sign in
Primary navigation
Search or go to…
Project
openfpm_pdata
Manage
Activity
Members
Labels
Plan
Issues
Issue boards
Milestones
Wiki
Code
Merge requests
Repository
Branches
Commits
Tags
Repository graph
Compare revisions
Build
Pipelines
Jobs
Pipeline schedules
Artifacts
Deploy
Releases
Model registry
Operate
Environments
Monitor
Incidents
Analyze
Value stream analytics
Contributor analytics
CI/CD analytics
Repository analytics
Model experiments
Help
Help
Support
GitLab documentation
Compare GitLab plans
Community forum
Contribute to GitLab
Provide feedback
Keyboard shortcuts
?
Snippets
Groups
Projects
Show more breadcrumbs
Sbalzarini Lab
Software
Parallel Computing
OpenFPM
openfpm_pdata
Commits
a2ce93f1
Commit
a2ce93f1
authored
3 years ago
by
Pietro Incardona
Browse files
Options
Downloads
Patches
Plain Diff
Memory BW added
parent
682a215e
No related branches found
Branches containing commit
No related tags found
No related merge requests found
Pipeline
#4013
passed
3 years ago
Stage: build
Stage: test
Changes
2
Pipelines
1
Hide whitespace changes
Inline
Side-by-side
Showing
2 changed files
example/Performance/memBW/Makefile
+62
-0
62 additions, 0 deletions
example/Performance/memBW/Makefile
example/Performance/memBW/main.cu
+62
-0
62 additions, 0 deletions
example/Performance/memBW/main.cu
with
124 additions
and
0 deletions
example/Performance/memBW/Makefile
0 → 100644
+
62
−
0
View file @
a2ce93f1
include
../../example.mk
### This is a trick to avoid "Command not found if you no not have NVCC compiler". In practice the normal C++ compiler is used
### internally the example disable with the preprocessor its code if not compiled with nvcc
CUDA_CC
=
CUDA_CC_LINK
=
CC
=
mpic++
ifdef
HIP
CUDA_CC
=
hipcc
CUDA_OPTIONS
=
-D__NVCC__
-D__HIP__
-DCUDART_VERSION
=
11000
-D__CUDACC__
-D__CUDACC_VER_MAJOR__
=
11
-D__CUDACC_VER_MINOR__
=
0
-D__CUDACC_VER_BUILD__
=
0
LIBS_SELECT
=
$(
LIBS
)
CC
=
hipcc
CUDA_CC_LINK
=
hipcc
else
ifdef
CUDA_ON_CPU
CUDA_CC
=
mpic++
-x
c++
$(
INCLUDE_PATH
)
INCLUDE_PATH_NVCC
=
CUDA_CC_LINK
=
mpic++
CUDA_OPTIONS
=
-D__NVCC__
-DCUDART_VERSION
=
11000
LIBS_SELECT
=
$(
LIBS
)
else
ifeq
(, $(shell which nvcc))
CUDA_CC
=
mpic++
-x
c++
$(
INCLUDE_PATH
)
INCLUDE_PATH_NVCC
=
CUDA_CC_LINK
=
mpic++
LIBS_SELECT
=
$(
LIBS
)
else
CUDA_CC
=
nvcc
-ccbin
=
mpic++
CUDA_CC_LINK
=
nvcc
-ccbin
=
mpic++
LIBS_SELECT
=
$(
LIBS_NVCC
)
endif
endif
endif
CC
=
mpic++
LDIR
=
OBJ
=
main.o
miniBUDE
:
%.o
:
%.cu
$(
CUDA_CC
)
-g
-O3
$(
CUDA_OPTIONS
)
$(
OPT
)
-c
--std
=
c++14
-o
$@
$<
$(
INCLUDE_PATH_NVCC
)
%.o
:
%.cpp
$(
CC
)
-g
-O3
$(
OPT
)
-g
-c
--std
=
c++14
-o
$@
$<
$(
INCLUDE_PATH
)
miniBUDE
:
$(OBJ)
$(
CUDA_CC_LINK
)
-o
$@
$^
$(
CFLAGS
)
$(
LIBS_PATH
)
$(
LIBS_SELECT
)
all
:
miniBUDE
run
:
miniBUDE
mpirun
--oversubscribe
-np
2 ./miniBUDE
.PHONY
:
clean all run
clean
:
rm
-f
*
.o
*
~ core miniBUDE
This diff is collapsed.
Click to expand it.
example/Performance/memBW/main.cu
0 → 100644
+
62
−
0
View file @
a2ce93f1
#include
"Vector/map_vector.hpp"
#include
"util/stat/common_statistics.hpp"
template
<
typename
vector_type
,
typename
vector_type2
>
__attribute__
((
always_inline
))
inline
__global__
void
translate_fill_prop
(
vector_type
&
vd_out
,
vector_type2
&
vd_in
)
{
auto
p
=
blockIdx
.
x
*
blockDim
.
x
+
threadIdx
.
x
;
vd_out
.
template
get
<
0
>(
p
)
=
vd_in
.
template
get
<
0
>(
p
)[
0
]
+
vd_in
.
template
get
<
0
>(
p
)[
1
];
vd_out
.
template
get
<
1
>(
p
)[
0
]
=
vd_in
.
template
get
<
0
>(
p
)[
0
];
vd_out
.
template
get
<
1
>(
p
)[
1
]
=
vd_in
.
template
get
<
0
>(
p
)[
1
];
vd_out
.
template
get
<
2
>(
p
)[
0
][
0
]
=
vd_in
.
template
get
<
0
>(
p
)[
0
];
vd_out
.
template
get
<
2
>(
p
)[
0
][
1
]
=
vd_in
.
template
get
<
0
>(
p
)[
1
];
vd_out
.
template
get
<
2
>(
p
)[
1
][
0
]
=
vd_in
.
template
get
<
0
>(
p
)[
0
]
+
vd_in
.
template
get
<
0
>(
p
)[
1
];
vd_out
.
template
get
<
2
>(
p
)[
1
][
1
]
=
vd_in
.
template
get
<
0
>(
p
)[
1
]
-
vd_in
.
template
get
<
0
>(
p
)[
0
];
vd_in
.
template
get
<
0
>(
p
)[
0
]
+=
0.01
f
;
vd_in
.
template
get
<
0
>(
p
)[
1
]
+=
0.01
f
;
}
int
main
(
int
argc
,
char
*
argv
[])
{
init_wrappers
();
openfpm
::
vector_gpu
<
aggregate
<
float
,
float
[
2
],
float
[
2
][
2
]
>>
out
;
openfpm
::
vector_gpu
<
aggregate
<
float
[
2
]
>>
in
;
int
nele
=
16777216
;
out
.
resize
(
nele
);
in
.
resize
(
nele
);
for
(
int
i
=
0
;
i
<
16777216
;
i
++
)
{
in
.
template
get
<
0
>(
i
)[
0
]
=
i
;
in
.
template
get
<
0
>(
i
)[
1
]
=
i
+
100.0
;
}
auto
ite
=
out
.
getGPUIterator
(
256
);
for
(
int
i
=
0
;
i
<
100
;
i
++
)
{
cudaDeviceSynchronize
();
timer
t
;
t
.
start
();
auto
vout
=
out
.
toKernel
();
auto
vin
=
in
.
toKernel
();
CUDA_LAUNCH
(
translate_fill_prop
,
ite
,
vout
,
vin
);
cudaDeviceSynchronize
();
t
.
stop
();
std
::
cout
<<
"Time: "
<<
t
.
getwct
()
<<
std
::
endl
;
std
::
cout
<<
"BW: "
<<
nele
*
4
*
19
/
t
.
getwct
()
*
1e-9
<<
" GB/s"
<<
std
::
endl
;
}
}
This diff is collapsed.
Click to expand it.
Preview
0%
Loading
Try again
or
attach a new file
.
Cancel
You are about to add
0
people
to the discussion. Proceed with caution.
Finish editing this message first!
Save comment
Cancel
Please
register
or
sign in
to comment