biod.pnpi.spb.ru
/
alexxy
/
gromacs.git
/ blobdiff
commit
grep
author
committer
pickaxe
?
search:
re
summary
|
shortlog
|
log
|
commit
|
commitdiff
|
tree
raw
|
inline
| side by side
Avoid using function calls in OpenMP directives
[alexxy/gromacs.git]
/
src
/
gromacs
/
mdlib
/
nbnxn_kernels
/
simd_4xn
/
nbnxn_kernel_simd_4xn.c
diff --git
a/src/gromacs/mdlib/nbnxn_kernels/simd_4xn/nbnxn_kernel_simd_4xn.c
b/src/gromacs/mdlib/nbnxn_kernels/simd_4xn/nbnxn_kernel_simd_4xn.c
index 0684149bad85f093838e42c941c8593d7aa05567..5974873a070147da2fd1723c8909d00449eb3dfe 100644
(file)
--- a/
src/gromacs/mdlib/nbnxn_kernels/simd_4xn/nbnxn_kernel_simd_4xn.c
+++ b/
src/gromacs/mdlib/nbnxn_kernels/simd_4xn/nbnxn_kernel_simd_4xn.c
@@
-272,6
+272,7
@@
nbnxn_kernel_simd_4xn(nbnxn_pairlist_set_t gmx_unused *nbl_list,
nbnxn_pairlist_t **nbl;
int coulkt, vdwkt = 0;
int nb;
nbnxn_pairlist_t **nbl;
int coulkt, vdwkt = 0;
int nb;
+ int nthreads gmx_unused;
nnbl = nbl_list->nnbl;
nbl = nbl_list->nbl;
nnbl = nbl_list->nnbl;
nbl = nbl_list->nbl;
@@
-343,7
+344,8
@@
nbnxn_kernel_simd_4xn(nbnxn_pairlist_set_t gmx_unused *nbl_list,
gmx_incons("Unsupported VdW interaction type");
}
gmx_incons("Unsupported VdW interaction type");
}
-#pragma omp parallel for schedule(static) num_threads(gmx_omp_nthreads_get(emntNonbonded))
+ nthreads = gmx_omp_nthreads_get(emntNonbonded);
+#pragma omp parallel for schedule(static) num_threads(nthreads)
for (nb = 0; nb < nnbl; nb++)
{
nbnxn_atomdata_output_t *out;
for (nb = 0; nb < nnbl; nb++)
{
nbnxn_atomdata_output_t *out;