| tFixed linarr memtransfer error - sphere - GPU-based 3D discrete element method… | |
| git clone git://src.adamsgaard.dk/sphere | |
| Log | |
| Files | |
| Refs | |
| LICENSE | |
| --- | |
| commit 4bcc1be72de8c27479a909e12eb194e587ceccb7 | |
| parent 46318ad0230597131e26c3766048dde88cd38796 | |
| Author: Anders Damsgaard <[email protected]> | |
| Date: Tue, 18 Dec 2012 13:15:03 +0100 | |
| Fixed linarr memtransfer error | |
| Diffstat: | |
| M src/raytracer.cuh | 7 ++++--- | |
| 1 file changed, 4 insertions(+), 3 deletions(-) | |
| --- | |
| diff --git a/src/raytracer.cuh b/src/raytracer.cuh | |
| t@@ -415,6 +415,8 @@ __host__ void DEM::render( | |
| Float* linarr; // Linear array to use for color visualization | |
| Float* dev_linarr; // Device linear array to use for color visualization | |
| + cudaMalloc((void**)&dev_linarr, np*sizeof(Float)); | |
| + checkForCudaErrors("Error during cudaMalloc of linear array"); | |
| std::string desc; // Description of parameter visualized | |
| std::string unit; // Unit of parameter values visualized | |
| unsigned int i; | |
| t@@ -492,9 +494,8 @@ __host__ void DEM::render( | |
| // Copy linarr to dev_linarr if required | |
| if (transfer == 1) { | |
| - cudaMalloc((void**)&dev_linarr, np*sizeof(Float)); | |
| - cudaMemcpy(dev_linarr, &linarr, np*sizeof(Float), cudaMemcpyHostTo… | |
| - checkForCudaErrors("Error during cudaMalloc or cudaMemcpy of linea… | |
| + cudaMemcpy(dev_linarr, linarr, np*sizeof(Float), cudaMemcpyHostToD… | |
| + checkForCudaErrors("Error during cudaMemcpy of linear array"); | |
| } | |
| // Start raytracing kernel |