Skip to content
GitLab
Projects
Groups
Snippets
/
Help
Help
Support
Community forum
Keyboard shortcuts
?
Submit feedback
Contribute to GitLab
Sign in
Toggle navigation
Menu
Open sidebar
Florian Weik
waLBerla
Commits
020743a1
Commit
020743a1
authored
Oct 24, 2018
by
Martin Bauer
Browse files
Precompute x and f allocation size of GPUField
parent
319909f0
Changes
2
Hide whitespace changes
Inline
Side-by-side
src/cuda/GPUField.h
View file @
020743a1
...
...
@@ -146,6 +146,9 @@ namespace cuda {
uint_t
ySize_
;
uint_t
zSize_
;
uint_t
fSize_
;
uint_t
xAllocSize_
;
uint_t
fAllocSize_
;
Layout
layout_
;
bool
usePitchedMem_
;
};
...
...
src/cuda/GPUField.impl.h
View file @
020743a1
...
...
@@ -58,6 +58,21 @@ GPUField<T>::GPUField( uint_t _xSize, uint_t _ySize, uint_t _zSize, uint_t _fSiz
pitchedPtr_
=
make_cudaPitchedPtr
(
NULL
,
extent
.
width
,
extent
.
width
,
extent
.
height
);
WALBERLA_CUDA_CHECK
(
cudaMalloc
(
&
pitchedPtr_
.
ptr
,
extent
.
width
*
extent
.
height
*
extent
.
depth
)
);
}
// allocation size is stored in pitched pointer
// pitched pointer stores the amount of padded region in bytes
// but we keep track of the size in #elements
WALBERLA_ASSERT_EQUAL
(
pitchedPtr_
.
pitch
%
sizeof
(
T
),
0
);
if
(
layout_
==
field
::
fzyx
)
{
xAllocSize_
=
pitchedPtr_
.
pitch
/
sizeof
(
T
);
fAllocSize_
=
fSize_
;
}
else
{
fAllocSize_
=
pitchedPtr_
.
pitch
/
sizeof
(
T
);
xAllocSize_
=
xSize_
+
2
*
nrOfGhostLayers_
;
}
}
...
...
@@ -242,15 +257,7 @@ bool GPUField<T>::operator==( const GPUField & o ) const
template
<
typename
T
>
uint_t
GPUField
<
T
>::
xAllocSize
()
const
{
if
(
layout_
==
field
::
fzyx
)
{
// allocation size is stored in pitched pointer
// pitched pointer stores the amount of padded region in bytes
// but we have to return the size in #elements
WALBERLA_ASSERT_EQUAL
(
pitchedPtr_
.
pitch
%
sizeof
(
T
),
0
);
return
pitchedPtr_
.
pitch
/
sizeof
(
T
);
}
return
xSize_
+
2
*
nrOfGhostLayers_
;
return
xAllocSize_
;
}
template
<
typename
T
>
...
...
@@ -268,12 +275,7 @@ uint_t GPUField<T>::zAllocSize() const
template
<
typename
T
>
uint_t
GPUField
<
T
>::
fAllocSize
()
const
{
if
(
layout_
==
field
::
zyxf
)
{
WALBERLA_ASSERT_EQUAL
(
pitchedPtr_
.
pitch
%
sizeof
(
T
),
0
);
return
pitchedPtr_
.
pitch
/
sizeof
(
T
);
}
return
fSize_
;
return
fAllocSize_
;
}
template
<
typename
T
>
...
...
Write
Preview
Supports
Markdown
0%
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!
Cancel
Please
register
or
sign in
to comment