commit 252717d95f1f36372b545b7253375b9b12b8f809
parent 28338e1177010cccfa121818b520098202ad39e1
Author: Vincent Forest <vincent.forest@meso-star.com>
Date: Mon, 11 May 2015 16:12:26 +0200
Make the library compliant with the CL compiler
Diffstat:
5 files changed, 46 insertions(+), 33 deletions(-)
diff --git a/cmake/CMakeLists.txt b/cmake/CMakeLists.txt
@@ -44,8 +44,6 @@ find_package(RSys 0.1.1 REQUIRED)
find_package(StarSP REQUIRED)
find_package(Star3D)
find_package(OpenMP)
-set(CMAKE_MODULE_PATH ${CMAKE_MODULE_PATH} ${RCMAKE_SOURCE_DIR})
-include(rcmake)
if(NOT OPENMP_FOUND)
message(STATUS "No OpenMP support: muti-threading is disabled")
@@ -53,6 +51,12 @@ endif()
include_directories(${RSys_INCLUDE_DIR} ${StarSP_INCLUDE_DIR})
+
+set(CMAKE_MODULE_PATH ${CMAKE_MODULE_PATH} ${RCMAKE_SOURCE_DIR})
+include(rcmake)
+include(rcmake_runtime)
+rcmake_append_runtime_dirs(_runtime_dirs RSys StarSP Star3D)
+
################################################################################
# Configure and define targets
################################################################################
@@ -83,9 +87,11 @@ set_target_properties(smc PROPERTIES
SOVERSION ${VERSION_MAJOR})
if(OPENMP_FOUND)
- set_target_properties(smc PROPERTIES
- COMPILE_FLAGS ${OpenMP_C_FLAGS}
- LINK_FLAGS ${OpenMP_C_FLAGS})
+ set_target_properties(smc PROPERTIES COMPILE_FLAGS ${OpenMP_C_FLAGS})
+
+ if(CMAKE_COMPILER_IS_GNUCC)
+ set_target_properties(smc PROPERTIES LINK_FLAGS ${OpenMP_C_FLAGS})
+ endif()
endif()
target_link_libraries(smc RSys StarSP)
@@ -109,10 +115,15 @@ if(NOT NO_TEST)
target_link_libraries(${_name} ${_lib})
endforeach(_lib)
add_test(${_name} ${_name})
+ rcmake_set_test_runtime_dirs(${_name} _runtime_dirs)
endfunction(new_test)
+ if(CMAKE_COMPILER_IS_GNUCC)
+ set(MATH_LIB m)
+ endif()
+
new_test(test_smc_device)
- new_test(test_smc_solve m)
+ new_test(test_smc_solve ${MATH_LIB})
if(NOT TARGET Star3D)
message(STATUS
diff --git a/src/smc.h b/src/smc.h
@@ -70,7 +70,7 @@ struct smc_type {
void (*add)(void* result, const void* op0, const void* op1);
void (*sub)(void* result, const void* op0, const void* op1);
void (*mul)(void* result, const void* op0, const void* op1);
- void (*divi)(void* result, const void* op0, const unsigned long op1);
+ void (*divi)(void* result, const void* op0, const size_t op1);
void (*sqrt)(void* result, const void* value);
};
@@ -82,17 +82,17 @@ struct smc_estimator_status {
void* E; /* Expected value */
void* V; /* Variance */
void* SE; /* Standard error, i.e. sqrt(V / N) */
- unsigned long N; /* Samples count */
+ size_t N; /* Samples count */
};
struct smc_integrator {
void (*integrand)(void* value, struct ssp_rng* rng, void* ctx);
const struct smc_type* type;
- unsigned long max_steps; /* Maximum # of experiments */
+ size_t max_steps; /* Maximum # of experiments */
};
static const struct smc_integrator SMC_INTEGRATOR_NULL =
-{ NULL, NULL, ULONG_MAX };
+{ NULL, NULL, SIZE_MAX };
/* Forward declaration of opaque types */
struct smc_device; /* Entry point of the library */
diff --git a/src/smc_device.c b/src/smc_device.c
@@ -54,7 +54,7 @@ device_release(ref_T* ref)
ssp_rng_ref_put(dev->rngs[i]);
}
sa_release(dev->rngs);
- MEM_FREE(dev->allocator, dev);
+ MEM_RM(dev->allocator, dev);
}
/*******************************************************************************
diff --git a/src/smc_integrator.c b/src/smc_integrator.c
@@ -44,7 +44,7 @@ struct smc_estimator {
struct smc_type type;
void* value;
void* square_value;
- unsigned long nsamples;
+ size_t nsamples;
struct smc_estimator_status status;
@@ -75,7 +75,8 @@ estimator_create
ref_init(&estimator->ref);
#define TYPE_CREATE(Dst) { \
- if(!((Dst) = type->create(dev->allocator))) { \
+ (Dst) = type->create(dev->allocator); \
+ if(!(Dst)) { \
res = RES_MEM_ERR; \
goto error; \
} \
@@ -122,16 +123,16 @@ estimator_release(ref_T* ref)
estimator = CONTAINER_OF(ref, struct smc_estimator, ref);
dev = estimator->dev;
if(estimator->value)
- MEM_FREE(dev->allocator, estimator->value);
+ MEM_RM(dev->allocator, estimator->value);
if(estimator->square_value)
- MEM_FREE(dev->allocator, estimator->square_value);
+ MEM_RM(dev->allocator, estimator->square_value);
if(estimator->status.E)
- MEM_FREE(dev->allocator, estimator->status.E);
+ MEM_RM(dev->allocator, estimator->status.E);
if(estimator->status.V)
- MEM_FREE(dev->allocator, estimator->status.V);
+ MEM_RM(dev->allocator, estimator->status.V);
if(estimator->status.SE)
- MEM_FREE(dev->allocator, estimator->status.SE);
- MEM_FREE(dev->allocator, estimator);
+ MEM_RM(dev->allocator, estimator->status.SE);
+ MEM_RM(dev->allocator, estimator);
SMC(device_ref_put(dev));
}
@@ -146,8 +147,8 @@ smc_solve
struct smc_estimator** out_estimator)
{
struct smc_estimator* estimator = NULL;
- unsigned long i;
- size_t nthreads = 0;
+ int64_t i;
+ int64_t nthreads = 0;
char* raw = NULL;
void** vals = NULL;
void** accums = NULL;
@@ -172,7 +173,8 @@ smc_solve
accums_sqr = (void**) (raw + 2 * sizeof(void*) * nthreads);
nsamples = (size_t*)(raw + 3 * sizeof(void*) * nthreads);
#define TYPE_CREATE(Dst) { \
- if(!((Dst) = integrator->type->create(dev->allocator))) { \
+ (Dst) = integrator->type->create(dev->allocator); \
+ if(!(Dst)) { \
res = RES_MEM_ERR; \
goto error; \
} \
@@ -192,7 +194,7 @@ smc_solve
/* Parallel evaluation of the simulation */
#pragma omp parallel for schedule(static)
- for(i = 0; i < integrator->max_steps; ++i) {
+ for(i = 0; i < (int64_t)integrator->max_steps; ++i) {
const int ithread = omp_get_thread_num();
integrator->integrand(vals[ithread], dev->rngs[ithread], ctx);
integrator->type->add(accums[ithread], accums[ithread], vals[ithread]);
@@ -216,7 +218,7 @@ exit:
if(accums[i]) integrator->type->destroy(dev->allocator, accums[i]);
if(accums_sqr[i]) integrator->type->destroy(dev->allocator, accums_sqr[i]);
}
- MEM_FREE(dev->allocator, raw);
+ MEM_RM(dev->allocator, raw);
}
if(out_estimator) *out_estimator = estimator;
return res;
@@ -238,8 +240,8 @@ smc_solve_N
struct smc_estimator* estimators[])
{
res_T res = RES_OK;
- size_t i;
- size_t nthreads;
+ int64_t i;
+ int64_t nthreads = 0;
void** vals = NULL;
if(!estimators) {
@@ -254,7 +256,7 @@ smc_solve_N
}
/* Create the estimators */
- FOR_EACH(i, 0, count ) {
+ FOR_EACH(i, 0, (int64_t)count) {
res = estimator_create(dev, integrator->type, estimators + i);
if(res != RES_OK) goto error;
}
@@ -272,8 +274,8 @@ smc_solve_N
/* Parallel estimation of N simulations */
#pragma omp parallel for schedule(static)
- for(i = 0; i < count; ++i) {
- unsigned long istep;
+ for(i = 0; i < (int64_t)count; ++i) {
+ size_t istep;
const int ithread = omp_get_thread_num();
FOR_EACH(istep, 0, integrator->max_steps) {
integrator->integrand
@@ -292,12 +294,12 @@ exit:
FOR_EACH(i, 0, nthreads) {
if(vals[i]) integrator->type->destroy(dev->allocator, vals[i]);
}
- MEM_FREE(dev->allocator, vals);
+ MEM_RM(dev->allocator, vals);
}
return res;
error:
if(estimators) {
- FOR_EACH(i, 0, count) {
+ FOR_EACH(i, 0, (int64_t)count) {
if(estimators[i]) {
SMC(estimator_ref_put(estimators[i]));
estimators[i] = NULL;
diff --git a/src/smc_type_real.h b/src/smc_type_real.h
@@ -50,7 +50,7 @@ static void
SMC_REAL_FUNC__(destroy)(struct mem_allocator* allocator, void* data)
{
ASSERT(data);
- MEM_FREE(allocator, data);
+ MEM_RM(allocator, data);
}
static void
@@ -92,7 +92,7 @@ SMC_REAL_FUNC__(mul)(void* result, const void* op0, const void* op1)
}
static void
-SMC_REAL_FUNC__(divi)(void* result, const void* op0, const unsigned long op1)
+SMC_REAL_FUNC__(divi)(void* result, const void* op0, const size_t op1)
{
ASSERT(result && op0 && op1);
*(SMC_TYPE_REAL*)result =