3b4f985db3e34f0f2186bb24016332814c5d618d
[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 is not a bug, but means a volatile pointer
87        to a volatile value. This is needed for older versions of
88        xlc. */
89     volatile char* volatile value __attribute__ ((aligned(64)));  /*!< Volatile, to avoid compiler aliasing */
90 }
91 tMPI_Atomic_ptr_t;
92
93
94 typedef struct tMPI_Spinlock
95 {
96     volatile int lock __attribute__ ((aligned(64)));
97 }
98 tMPI_Spinlock_t;
99 #define TMPI_ATOMIC_HAVE_NATIVE_SPINLOCK
100
101
102
103
104 #define tMPI_Atomic_get(a)   (int)((a)->value)
105 #define tMPI_Atomic_set(a, i)  (((a)->value) = (i))
106 #define tMPI_Atomic_ptr_get(a)   ((a)->value)
107 #define tMPI_Atomic_ptr_set(a, i)  (((a)->value) = (i))
108
109 #define TMPI_SPINLOCK_INITIALIZER   { 0 }
110
111
112 static inline int tMPI_Atomic_cas(tMPI_Atomic_t *a, int oldval, int newval)
113 {
114 #ifdef TMPI_XLC_INTRINSICS
115     int ret;
116
117     __fence(); /* this one needs to be here to avoid ptr. aliasing issues */
118     __eieio();
119     ret = (__compare_and_swap(&(a->value), &oldval, newval));
120     __isync();
121     __fence(); /* and this one needs to be here to avoid aliasing issues */
122     return ret;
123 #else
124     int prev;
125     __asm__ __volatile__ ("1:    lwarx   %0,0,%2 \n"
126                           "\t cmpw    0,%0,%3 \n"
127                           "\t bne     2f \n"
128                           "\t stwcx.  %4,0,%2 \n"
129                           "\t bne-    1b \n"
130                           "\t sync \n"
131                           "2: \n"
132                           : "=&r" (prev), "=m" (a->value)
133                           : "r" (&a->value), "r" (oldval), "r" (newval),
134                           "m" (a->value));
135
136     return prev == oldval;
137 #endif
138 }
139
140
141 static inline int tMPI_Atomic_ptr_cas(tMPI_Atomic_ptr_t *a, void* oldval,
142                                       void* newval)
143 {
144     int ret;
145     volatile char* volatile oldv = (char*)oldval;
146     volatile char* volatile newv = (char*)newval;
147
148     __fence(); /* this one needs to be here to avoid ptr. aliasing issues */
149     __eieio();
150 #if (!defined (__LP64__) ) && (!defined(__powerpc64__) )
151     ret = __compare_and_swap((int *)&(a->value), (int*)&oldv, (int)newv);
152 #else
153     ret = __compare_and_swaplp((long *)&(a->value), (long*)&oldv, (long)newv);
154 #endif
155     __isync();
156     __fence();
157
158     return ret;
159 }
160
161
162
163
164 static inline int tMPI_Atomic_add_return(tMPI_Atomic_t *a, int i)
165 {
166 #ifdef TMPI_XLC_INTRINSICS
167     int oldval, newval;
168
169     do
170     {
171         __fence();
172         __eieio(); /* these memory barriers are neccesary */
173         oldval = tMPI_Atomic_get(a);
174         newval = oldval + i;
175     }
176     /*while(!__compare_and_swap( &(a->value), &oldval, newval));*/
177     while (__check_lock_mp( (int*)&(a->value), oldval, newval));
178
179     /*__isync();*/
180     __fence();
181
182     return newval;
183 #else
184     int t;
185
186     __asm__ __volatile__("1:     lwarx   %0,0,%2 \n"
187                          "\t add     %0,%1,%0 \n"
188                          "\t stwcx.  %0,0,%2 \n"
189                          "\t bne-    1b \n"
190                          "\t isync \n"
191                          : "=&r" (t)
192                          : "r" (i), "r" (&a->value) );
193     return t;
194 #endif
195 }
196 #define TMPI_ATOMIC_HAVE_NATIVE_ADD_RETURN
197
198
199
200 static inline int tMPI_Atomic_fetch_add(tMPI_Atomic_t *a, int i)
201 {
202 #ifdef TMPI_XLC_INTRINSICS
203     int oldval, newval;
204
205     do
206     {
207         __fence();
208         __eieio(); /* these memory barriers are neccesary */
209         oldval = tMPI_Atomic_get(a);
210         newval = oldval + i;
211     }
212     /*while(__check_lock_mp((const int*)&(a->value), oldval, newval));*/
213     while (__check_lock_mp( (int*)&(a->value), oldval, newval));
214     /*while(!__compare_and_swap( &(a->value), &oldval, newval));*/
215     /*__isync();*/
216     __fence();
217
218     return oldval;
219 #else
220     int t;
221
222     __asm__ __volatile__("\t eieio\n"
223                          "1:     lwarx   %0,0,%2 \n"
224                          "\t add     %0,%1,%0 \n"
225                          "\t stwcx.  %0,0,%2 \n"
226                          "\t bne-    1b \n"
227                          "\t isync \n"
228                          : "=&r" (t)
229                          : "r" (i), "r" (&a->value));
230
231     return (t - i);
232 #endif
233 }
234 #define TMPI_ATOMIC_HAVE_NATIVE_FETCH_ADD
235
236
237 static inline void tMPI_Spinlock_init(tMPI_Spinlock_t *x)
238 {
239     __fence();
240     __clear_lock_mp((const int*)x, 0);
241     __fence();
242 }
243
244
245 static inline void tMPI_Spinlock_lock(tMPI_Spinlock_t *x)
246 {
247     __fence();
248     do
249     {
250     }
251     while (__check_lock_mp((int*)&(x->lock), 0, 1));
252     tMPI_Atomic_memory_barrier_acq();
253 }
254
255
256 static inline int tMPI_Spinlock_trylock(tMPI_Spinlock_t *x)
257 {
258     int ret;
259     /* Return 0 if we got the lock */
260     __fence();
261     ret = __check_lock_mp((int*)&(x->lock), 0, 1);
262     tMPI_Atomic_memory_barrier_acq();
263     return ret;
264 }
265
266
267 static inline void tMPI_Spinlock_unlock(tMPI_Spinlock_t *x)
268 {
269     tMPI_Atomic_memory_barrier_rel();
270     __clear_lock_mp((int*)&(x->lock), 0);
271 }
272
273
274 static inline int tMPI_Spinlock_islocked(const tMPI_Spinlock_t *x)
275 {
276     int ret;
277     __fence();
278     ret = ((x->lock) != 0);
279     tMPI_Atomic_memory_barrier_acq();
280     return ret;
281 }
282
283
284 static inline void tMPI_Spinlock_wait(tMPI_Spinlock_t *x)
285 {
286     do
287     {
288     }
289     while (tMPI_Spinlock_islocked(x));
290 }