Lattice Land Core Library
Loading...
Searching...
No Matches
fixpoint.hpp
Go to the documentation of this file.
1// Copyright 2022 Pierre Talbot
2
3#ifndef FIXPOINT_HPP
4#define FIXPOINT_HPP
5
6#include "logic/logic.hpp"
8#include "battery/memory.hpp"
9
10#ifdef __CUDACC__
11 #include <cooperative_groups.h>
12#endif
13
14namespace lala {
15
16/** A simple form of sequential fixpoint computation based on Kleene fixpoint.
17 * At each iteration, the refinement operations \f$ f_1, \ldots, f_n \f$ are simply composed by functional composition \f$ f = f_n \circ \ldots \circ f_1 \f$.
18 * This strategy basically corresponds to the Gauss-Seidel iteration method. */
20public:
21 CUDA void barrier() {}
22
23 template <class A>
24 CUDA void iterate(A& a, local::BInc& has_changed) {
25 size_t n = a.num_refinements();
26 for(size_t i = 0; i < n; ++i) {
27 a.refine(i, has_changed);
28 }
29 }
30
31 template <class A>
32 CUDA size_t fixpoint(A& a, local::BInc& has_changed) {
33 size_t iterations = 0;
34 local::BInc changed(true);
35 while(changed && !a.is_top()) {
36 changed.dtell_bot();
37 iterate(a, changed);
38 has_changed.tell(changed);
39 iterations++;
40 }
41 return iterations;
42 }
43
44 template <class A>
45 CUDA local::BInc fixpoint(A& a) {
46 local::BInc has_changed(false);
47 fixpoint(a, has_changed);
48 return has_changed;
49 }
50};
51
52#ifdef __CUDACC__
53
54/** A simple form of fixpoint computation based on Kleene fixpoint.
55 * At each iteration, the refinement operations \f$ f_1, \ldots, f_n \f$ are composed by parallel composition \f$ f = f_1 \| \ldots \| f_n \f$ meaning they are executed in parallel by different threads.
56 * This is called an asynchronous iteration and it is due to (Cousot, Asynchronous iterative methods for solving a fixed point system of monotone equations in a complete lattice, 1977).
57 * The underlying lattice on which we iterate must provide two methods:
58 * - `a.refine(int, BInc&)`: call the ith refinement functions and set `has_changed` to `true` if `a` has changed.
59 * - `a.num_refinements()`: return the number of refinement functions.
60 * \tparam Group is a CUDA cooperative group class.
61 * \tparam Memory is an atomic memory, that must be compatible with the cooperative group chosen (e.g., don't use atomic_memory_block if the group contains multiple blocks). */
62template <class Group, class Memory, class Allocator>
63class AsynchronousIterationGPU {
64public:
65 using memory_type = Memory;
66 using allocator_type = Allocator;
67 using group_type = Group;
68private:
69 using atomic_binc = BInc<memory_type>;
70 battery::vector<atomic_binc, allocator_type> changed;
71 battery::vector<atomic_binc, allocator_type> is_top;
72 Group group;
73
74 CUDA void assert_cuda_arch() {
75 printf("AsynchronousIterationGPU must be used on the GPU device only.\n");
76 assert(0);
77 }
78
79 CUDA void reset() {
80 changed[0].tell_top();
81 changed[1].dtell_bot();
82 changed[2].dtell_bot();
83 for(int i = 0; i < is_top.size(); ++i) {
84 is_top[i].dtell_bot();
85 }
86 }
87
88public:
89 CUDA AsynchronousIterationGPU(const Group& group, const allocator_type& alloc = allocator_type()):
90 group(group), changed(3, alloc), is_top(3, alloc)
91 {}
92
93 CUDA void barrier() {
94 #ifndef __CUDA_ARCH__
95 assert_cuda_arch();
96 #else
97 cooperative_groups::sync(group);
98 #endif
99 }
100
101 template <class A, class M>
102 CUDA void iterate(A& a, BInc<M>& has_changed) {
103 #ifndef __CUDA_ARCH__
104 assert_cuda_arch();
105 #else
106 size_t n = a.num_refinements();
107 for (size_t t = group.thread_rank(); t < n; t += group.num_threads()) {
108 a.refine(t, has_changed);
109 if((t-group.thread_rank()) + group.num_threads() < n) __syncwarp();
110 }
111 #endif
112 }
113
114 template <class A, class M>
115 CUDA size_t fixpoint(A& a, BInc<M>& has_changed, volatile bool* stop) {
116 #ifndef __CUDA_ARCH__
117 assert_cuda_arch();
118 return 0;
119 #else
120 reset();
121 barrier();
122 size_t i;
123 for(i = 1; changed[(i-1)%3] && !is_top[(i-1)%3]; ++i) {
124 iterate(a, changed[i%3]);
125 changed[(i+1)%3].dtell_bot(); // reinitialize changed for the next iteration.
126 is_top[i%3].tell(a.is_top());
127 is_top[i%3].tell(local::BInc{*stop});
128 barrier();
129 }
130 // It changes if we performed several iteration, or if the first iteration changed the abstract domain.
131 has_changed.tell(changed[1]);
132 has_changed.tell(changed[2]);
133 return i - 1;
134 #endif
135 }
136
137 template <class A, class M>
138 CUDA size_t fixpoint(A& a, BInc<M>& has_changed) {
139 bool stop = false;
140 return fixpoint(a, has_changed, &stop);
141 }
142
143 template <class A>
144 CUDA local::BInc fixpoint(A& a) {
145 local::BInc has_changed(false);
146 fixpoint(a, has_changed);
147 return has_changed;
148 }
149};
150
151template <class Allocator>
152using BlockAsynchronousIterationGPU = AsynchronousIterationGPU<cooperative_groups::thread_block, battery::atomic_memory_block, Allocator>;
153
154template <class Allocator>
155using GridAsynchronousIterationGPU = AsynchronousIterationGPU<cooperative_groups::grid_group, battery::atomic_memory_grid, Allocator>;
156
157#endif
158
159} // namespace lala
160
161#endif
Definition fixpoint.hpp:19
CUDA void iterate(A &a, local::BInc &has_changed)
Definition fixpoint.hpp:24
CUDA local::BInc fixpoint(A &a)
Definition fixpoint.hpp:45
CUDA void barrier()
Definition fixpoint.hpp:21
CUDA size_t fixpoint(A &a, local::BInc &has_changed)
Definition fixpoint.hpp:32
Definition primitive_upset.hpp:118
CUDA constexpr local::BInc is_top() const
Definition primitive_upset.hpp:224
CUDA constexpr this_type & dtell_bot()
Definition primitive_upset.hpp:259
CUDA constexpr this_type & tell(const this_type2< M1 > &other, BInc< M2 > &has_changed)
Definition primitive_upset.hpp:239
Definition abstract_deps.hpp:14