Merge branch 'release-4-6'
[alexxy/gromacs.git] / src / gromacs / legacyheaders / thread_mpi / atomic / xlc_ppc.h
1 /*
2    This source code file is part of thread_mpi.
3    Written by Sander Pronk, Erik Lindahl, and possibly others.
4
5    Copyright (c) 2009, Sander Pronk, Erik Lindahl.
6    All rights reserved.
7
8    Redistribution and use in source and binary forms, with or without
9    modification, are permitted provided that the following conditions are met:
10    1) Redistributions of source code must retain the above copyright
11    notice, this list of conditions and the following disclaimer.
12    2) Redistributions in binary form must reproduce the above copyright
13    notice, this list of conditions and the following disclaimer in the
14    documentation and/or other materials provided with the distribution.
15    3) Neither the name of the copyright holders nor the
16    names of its contributors may be used to endorse or promote products
17    derived from this software without specific prior written permission.
18
19    THIS SOFTWARE IS PROVIDED BY US ''AS IS'' AND ANY
20    EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE IMPLIED
21    WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE ARE
22    DISCLAIMED. IN NO EVENT SHALL WE BE LIABLE FOR ANY
23    DIRECT, INDIRECT, INCIDENTAL, SPECIAL, EXEMPLARY, OR CONSEQUENTIAL DAMAGES
24    (INCLUDING, BUT NOT LIMITED TO, PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES;
25    LOSS OF USE, DATA, OR PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND
26    ON ANY THEORY OF LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT
27    (INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
28    SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
29
30    If you want to redistribute modifications, please consider that
31    scientific software is very special. Version control is crucial -
32    bugs must be traceable. We will be happy to consider code for
33    inclusion in the official distribution, but derived work should not
34    be called official thread_mpi. Details are found in the README & COPYING
35    files.
36  */
37
38 /* PowerPC using xlC inline assembly.
39  * Recent versions of xlC (>=7.0) _partially_ support GCC inline assembly
40  * if you use the option -qasm=gcc but we have had to hack things a bit, in
41  * particular when it comes to clobbered variables. Since this implementation
42  * _could_ be buggy, we have separated it from the known-to-be-working gcc
43  * one above.
44  *
45  * For now, we just disable the inline keyword if we're compiling C code:
46  */
47 #if 1
48 #if (!defined(__cplusplus)) && (!defined(inline))
49 #define inline_defined_in_atomic 1
50 #define inline
51 #endif
52 #endif
53
54
55 /* IBM xlC compiler */
56 #ifdef __cplusplus
57 #include <builtins.h>
58 #endif
59
60
61 #define TMPI_XLC_INTRINSICS
62
63 /* ppc has many memory synchronization instructions */
64 /*#define tMPI_Atomic_memory_barrier() { __fence(); __sync(); __fence();}*/
65 /*#define tMPI_Atomic_memory_barrier() __isync();*/
66 /*#define tMPI_Atomic_memory_barrier() __lwsync();*/
67
68 /* for normal memory, this should be enough: */
69 #define tMPI_Atomic_memory_barrier() { __fence(); __eieio(); __fence(); }
70 #define tMPI_Atomic_memory_barrier_acq() { __eieio(); __fence(); }
71 #define tMPI_Atomic_memory_barrier_rel() { __fence(); __eieio(); }
72 #define TMPI_HAVE_ACQ_REL_BARRIERS
73
74 /*#define tMPI_Atomic_memory_barrier() __eieio();*/
75
76
77 typedef struct tMPI_Atomic
78 {
79     volatile int value __attribute__ ((aligned(64)));
80 }
81 tMPI_Atomic_t;
82
83
84 typedef struct tMPI_Atomic_ptr
85 {
86     volatile char* volatile* value __attribute__ ((aligned(64)));  /*!< Volatile, to avoid compiler aliasing */
87 }
88 tMPI_Atomic_ptr_t;
89
90
91 typedef struct tMPI_Spinlock
92 {
93     volatile int lock __attribute__ ((aligned(64)));
94 }
95 tMPI_Spinlock_t;
96 #define TMPI_ATOMIC_HAVE_NATIVE_SPINLOCK
97
98
99
100
101 #define tMPI_Atomic_get(a)   (int)((a)->value)
102 #define tMPI_Atomic_set(a, i)  (((a)->value) = (i))
103 #define tMPI_Atomic_ptr_get(a)   ((a)->value)
104 #define tMPI_Atomic_ptr_set(a, i)  (((a)->value) = (i))
105
106 #define TMPI_SPINLOCK_INITIALIZER   { 0 }
107
108
109 static inline int tMPI_Atomic_cas(tMPI_Atomic_t *a, int oldval, int newval)
110 {
111 #ifdef TMPI_XLC_INTRINSICS
112     int ret;
113
114     __fence(); /* this one needs to be here to avoid ptr. aliasing issues */
115     __eieio();
116     ret = (__compare_and_swap(&(a->value), &oldval, newval));
117     __isync();
118     __fence(); /* and this one needs to be here to avoid aliasing issues */
119     return ret;
120 #else
121     int prev;
122     __asm__ __volatile__ ("1:    lwarx   %0,0,%2 \n"
123                           "\t cmpw    0,%0,%3 \n"
124                           "\t bne     2f \n"
125                           "\t stwcx.  %4,0,%2 \n"
126                           "\t bne-    1b \n"
127                           "\t sync \n"
128                           "2: \n"
129                           : "=&r" (prev), "=m" (a->value)
130                           : "r" (&a->value), "r" (oldval), "r" (newval),
131                           "m" (a->value));
132
133     return prev == oldval;
134 #endif
135 }
136
137
138 static inline int tMPI_Atomic_ptr_cas(tMPI_Atomic_ptr_t *a, void* oldval,
139                                       void* newval)
140 {
141     int ret;
142     volatile char* volatile* oldv = oldval;
143     volatile char* volatile* newv = newval;
144
145     __fence(); /* this one needs to be here to avoid ptr. aliasing issues */
146     __eieio();
147 #if (!defined (__LP64__) ) && (!defined(__powerpc64__) )
148     ret = __compare_and_swap((int *)&(a->value), (int*)&oldv, (int)newv);
149 #else
150     ret = __compare_and_swaplp((long *)&(a->value), (long*)&oldv, (long)newv);
151 #endif
152     __isync();
153     __fence();
154
155     return ret;
156 }
157
158
159
160
161 static inline int tMPI_Atomic_add_return(tMPI_Atomic_t *a, int i)
162 {
163 #ifdef TMPI_XLC_INTRINSICS
164     int oldval, newval;
165
166     do
167     {
168         __fence();
169         __eieio(); /* these memory barriers are neccesary */
170         oldval = tMPI_Atomic_get(a);
171         newval = oldval + i;
172     }
173     /*while(!__compare_and_swap( &(a->value), &oldval, newval));*/
174     while (__check_lock_mp( (int*)&(a->value), oldval, newval));
175
176     /*__isync();*/
177     __fence();
178
179     return newval;
180 #else
181     int t;
182
183     __asm__ __volatile__("1:     lwarx   %0,0,%2 \n"
184                          "\t add     %0,%1,%0 \n"
185                          "\t stwcx.  %0,0,%2 \n"
186                          "\t bne-    1b \n"
187                          "\t isync \n"
188                          : "=&r" (t)
189                          : "r" (i), "r" (&a->value) );
190     return t;
191 #endif
192 }
193 #define TMPI_ATOMIC_HAVE_NATIVE_ADD_RETURN
194
195
196
197 static inline int tMPI_Atomic_fetch_add(tMPI_Atomic_t *a, int i)
198 {
199 #ifdef TMPI_XLC_INTRINSICS
200     int oldval, newval;
201
202     do
203     {
204         __fence();
205         __eieio(); /* these memory barriers are neccesary */
206         oldval = tMPI_Atomic_get(a);
207         newval = oldval + i;
208     }
209     /*while(__check_lock_mp((const int*)&(a->value), oldval, newval));*/
210     while (__check_lock_mp( (int*)&(a->value), oldval, newval));
211     /*while(!__compare_and_swap( &(a->value), &oldval, newval));*/
212     /*__isync();*/
213     __fence();
214
215     return oldval;
216 #else
217     int t;
218
219     __asm__ __volatile__("\t eieio\n"
220                          "1:     lwarx   %0,0,%2 \n"
221                          "\t add     %0,%1,%0 \n"
222                          "\t stwcx.  %0,0,%2 \n"
223                          "\t bne-    1b \n"
224                          "\t isync \n"
225                          : "=&r" (t)
226                          : "r" (i), "r" (&a->value));
227
228     return (t - i);
229 #endif
230 }
231 #define TMPI_ATOMIC_HAVE_NATIVE_FETCH_ADD
232
233
234 static inline void tMPI_Spinlock_init(tMPI_Spinlock_t *x)
235 {
236     __fence();
237     __clear_lock_mp((const int*)x, 0);
238     __fence();
239 }
240
241
242 static inline void tMPI_Spinlock_lock(tMPI_Spinlock_t *x)
243 {
244     __fence();
245     do
246     {
247     }
248     while (__check_lock_mp((int*)&(x->lock), 0, 1));
249     tMPI_Atomic_memory_barrier_acq();
250 }
251
252
253 static inline int tMPI_Spinlock_trylock(tMPI_Spinlock_t *x)
254 {
255     int ret;
256     /* Return 0 if we got the lock */
257     __fence();
258     ret = __check_lock_mp((int*)&(x->lock), 0, 1);
259     tMPI_Atomic_memory_barrier_acq();
260     return ret;
261 }
262
263
264 static inline void tMPI_Spinlock_unlock(tMPI_Spinlock_t *x)
265 {
266     tMPI_Atomic_memory_barrier_rel();
267     __clear_lock_mp((int*)&(x->lock), 0);
268 }
269
270
271 static inline int tMPI_Spinlock_islocked(const tMPI_Spinlock_t *x)
272 {
273     int ret;
274     __fence();
275     ret = ((x->lock) != 0);
276     tMPI_Atomic_memory_barrier_acq();
277     return ret;
278 }
279
280
281 static inline void tMPI_Spinlock_wait(tMPI_Spinlock_t *x)
282 {
283     do
284     {
285     }
286     while (tMPI_Spinlock_islocked(x));
287 }