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