p2p_distr_mpi.hpp 41.6 KB
Newer Older
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
#ifndef P2P_DISTR_MPI_HPP
#define P2P_DISTR_MPI_HPP

#include <mpi.h>

#include <vector>
#include <memory>
#include <cassert>

#include <type_traits>
#include <omp.h>
#include <algorithm>

#include "scope_timer.hpp"
#include "particles_utils.hpp"
#include "p2p_tree.hpp"
17
#include "lock_free_bool_array.hpp"
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41

template <class partsize_t, class real_number>
class p2p_distr_mpi {
protected:
    static const int MaxNbRhs = 100;

    enum MpiTag{
        TAG_NB_PARTICLES,
        TAG_POSITION_PARTICLES,
        TAG_RESULT_PARTICLES,
    };

    struct NeighborDescriptor{
        partsize_t nbParticlesToExchange;
        int destProc;
        int nbLevelsToExchange;
        bool isRecv;

        std::unique_ptr<real_number[]> toRecvAndMerge;
        std::unique_ptr<real_number[]> toCompute;
        std::unique_ptr<real_number[]> results;
    };

    enum Action{
42
        NOTHING_TODO = 512,
43
44
45
        RECV_PARTICLES,
        COMPUTE_PARTICLES,
        RELEASE_BUFFER_PARTICLES,
46
        MERGE_PARTICLES
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
    };

    MPI_Comm current_com;

    int my_rank;
    int nb_processes;
    int nb_processes_involved;

    const std::pair<int,int> current_partition_interval;
    const int current_partition_size;
    const std::array<size_t,3> field_grid_dim;

    std::unique_ptr<int[]> partition_interval_size_per_proc;
    std::unique_ptr<int[]> partition_interval_offset_per_proc;

    std::unique_ptr<partsize_t[]> current_offset_particles_for_partition;

    std::vector<std::pair<Action,int>> whatNext;
    std::vector<MPI_Request> mpiRequests;
    std::vector<NeighborDescriptor> neigDescriptors;

    std::array<real_number,3> spatial_box_width;
    std::array<real_number,3> spatial_box_offset;

71
    const real_number cutoff_radius_compute;
72
    const int nb_cells_factor;
73
74
75
    const real_number cutoff_radius;
    std::array<long int,3> nb_cell_levels;

76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
    template <class DataType, int sizeElement>
    static void permute_copy(const partsize_t offsetIdx, const partsize_t nbElements,
                             const std::pair<long int,partsize_t> permutation[],
                             DataType data[], std::vector<unsigned char>* buffer){
        buffer->resize(nbElements*sizeof(DataType)*sizeElement);
        DataType* dataBuffer = reinterpret_cast<DataType*>(buffer->data());

        // Permute
        for(partsize_t idxPart = 0 ; idxPart < nbElements ; ++idxPart){
            const partsize_t srcData = permutation[idxPart].second;
            const partsize_t destData = idxPart;
            for(int idxVal = 0 ; idxVal < sizeElement ; ++idxVal){
                dataBuffer[destData*sizeElement + idxVal]
                        = data[srcData*sizeElement + idxVal];
            }
        }

        // Copy back
        for(partsize_t idxPart = 0 ; idxPart < nbElements ; ++idxPart){
            const partsize_t srcData = idxPart;
            const partsize_t destData = idxPart+offsetIdx;
            for(int idxVal = 0 ; idxVal < sizeElement ; ++idxVal){
                data[destData*sizeElement + idxVal]
                        = dataBuffer[srcData*sizeElement + idxVal];
            }
        }
    }

104
    static int foundGridFactor(const real_number in_cutoff_radius, const std::array<real_number,3>& in_spatial_box_width){
105
106
107
108
        int idx_factor = 1;
        while(in_cutoff_radius <= in_spatial_box_width[IDX_Z]/real_number(idx_factor+1)){
            idx_factor += 1;
        }
109
        return idx_factor;
110
111
    }

112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
public:
    ////////////////////////////////////////////////////////////////////////////

    p2p_distr_mpi(MPI_Comm in_current_com,
                     const std::pair<int,int>& in_current_partitions,
                     const std::array<size_t,3>& in_field_grid_dim,
                     const std::array<real_number,3>& in_spatial_box_width,
                     const std::array<real_number,3>& in_spatial_box_offset,
                     const real_number in_cutoff_radius)
        : current_com(in_current_com),
            my_rank(-1), nb_processes(-1),nb_processes_involved(-1),
            current_partition_interval(in_current_partitions),
            current_partition_size(current_partition_interval.second-current_partition_interval.first),
            field_grid_dim(in_field_grid_dim),
            spatial_box_width(in_spatial_box_width), spatial_box_offset(in_spatial_box_offset),
127
            cutoff_radius_compute(in_cutoff_radius),
128
129
            nb_cells_factor(foundGridFactor(in_cutoff_radius, in_spatial_box_width)),
            cutoff_radius(in_spatial_box_width[IDX_Z]/real_number(nb_cells_factor)){
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158

        AssertMpi(MPI_Comm_rank(current_com, &my_rank));
        AssertMpi(MPI_Comm_size(current_com, &nb_processes));

        partition_interval_size_per_proc.reset(new int[nb_processes]);
        AssertMpi( MPI_Allgather( const_cast<int*>(&current_partition_size), 1, MPI_INT,
                                  partition_interval_size_per_proc.get(), 1, MPI_INT,
                                  current_com) );
        assert(partition_interval_size_per_proc[my_rank] == current_partition_size);

        partition_interval_offset_per_proc.reset(new int[nb_processes+1]);
        partition_interval_offset_per_proc[0] = 0;
        for(int idxProc = 0 ; idxProc < nb_processes ; ++idxProc){
            partition_interval_offset_per_proc[idxProc+1] = partition_interval_offset_per_proc[idxProc] + partition_interval_size_per_proc[idxProc];
        }

        current_offset_particles_for_partition.reset(new partsize_t[current_partition_size+1]);

        nb_processes_involved = nb_processes;
        while(nb_processes_involved != 0 && partition_interval_size_per_proc[nb_processes_involved-1] == 0){
            nb_processes_involved -= 1;
        }
        assert(nb_processes_involved != 0);
        for(int idx_proc_involved = 0 ; idx_proc_involved < nb_processes_involved ; ++idx_proc_involved){
            assert(partition_interval_size_per_proc[idx_proc_involved] != 0);
        }

        assert(int(field_grid_dim[IDX_Z]) == partition_interval_offset_per_proc[nb_processes_involved]);

159
160
161
        nb_cell_levels[IDX_X] = nb_cells_factor;
        nb_cell_levels[IDX_Y] = nb_cells_factor;
        nb_cell_levels[IDX_Z] = nb_cells_factor;
162
163
164
165
166
167
    }

    virtual ~p2p_distr_mpi(){}

    ////////////////////////////////////////////////////////////////////////////

168
169
170
171
172
173
174
175
    int getGridFactor() const{
        return nb_cells_factor;
    }

    real_number getGridCutoff() const{
        return cutoff_radius;
    }

176
177
178
179
180
    long int get_cell_coord_x_from_index(const long int index) const{
        return index % nb_cell_levels[IDX_X];
    }

    long int get_cell_coord_y_from_index(const long int index) const{
181
        return (index % (nb_cell_levels[IDX_X]*nb_cell_levels[IDX_Y]))
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
                / nb_cell_levels[IDX_X];
    }

    long int get_cell_coord_z_from_index(const long int index) const{
        return index / (nb_cell_levels[IDX_X]*nb_cell_levels[IDX_Y]);
    }

    long int first_cell_level_proc(const int dest_proc) const{
        const real_number field_section_width_z = spatial_box_width[IDX_Z]/real_number(field_grid_dim[IDX_Z]);
        return static_cast<long int>((field_section_width_z*real_number(partition_interval_offset_per_proc[dest_proc]))/cutoff_radius);
    }

    long int last_cell_level_proc(const int dest_proc) const{
        const real_number field_section_width_z = spatial_box_width[IDX_Z]/real_number(field_grid_dim[IDX_Z]);
        return static_cast<long int>((field_section_width_z*real_number(partition_interval_offset_per_proc[dest_proc+1])
                                     - std::numeric_limits<real_number>::epsilon())/cutoff_radius);
    }

200
201
202
203
204
205
206
207
208
209
    real_number apply_pbc(real_number pos, IDXS_3D dim) const{
        while( pos < spatial_box_offset[dim] ){
            pos += spatial_box_width[dim];
        }
        while( spatial_box_width[dim]+spatial_box_offset[dim] <= pos){
            pos -= spatial_box_width[dim];
        }
        return pos;
    }

210
211
    std::array<long int,3> get_cell_coordinate(const real_number pos_x, const real_number pos_y,
                                               const real_number pos_z) const {
212
213
214
        const real_number diff_x = apply_pbc(pos_x,IDX_X) - spatial_box_offset[IDX_X];
        const real_number diff_y = apply_pbc(pos_y,IDX_Y) - spatial_box_offset[IDX_Y];
        const real_number diff_z = apply_pbc(pos_z,IDX_Z) - spatial_box_offset[IDX_Z];
215
216
217
218
219
220
221
222
223
224
225
226
227
228
        std::array<long int,3> coord;
        coord[IDX_X] = static_cast<long int>(diff_x/cutoff_radius);
        coord[IDX_Y] = static_cast<long int>(diff_y/cutoff_radius);
        coord[IDX_Z] = static_cast<long int>(diff_z/cutoff_radius);
        return coord;
    }

    long int get_cell_idx(const real_number pos_x, const real_number pos_y,
                          const real_number pos_z) const {
        std::array<long int,3> coord = get_cell_coordinate(pos_x, pos_y, pos_z);
        return ((coord[IDX_Z]*nb_cell_levels[IDX_Y])+coord[IDX_Y])*nb_cell_levels[IDX_X]+coord[IDX_X];
    }

    real_number compute_distance_r2(const real_number x1, const real_number y1, const real_number z1,
Berenger Bramas's avatar
Berenger Bramas committed
229
230
                                    const real_number x2, const real_number y2, const real_number z2,
                                    const real_number xshift_coef, const real_number yshift_coef, const real_number zshift_coef) const {
231
        real_number diff_x = std::abs(apply_pbc(x1,IDX_X)-apply_pbc(x2,IDX_X)+xshift_coef*spatial_box_width[IDX_X]);
Berenger Bramas's avatar
Berenger Bramas committed
232
        assert(diff_x <= 2*cutoff_radius);
233

234
        real_number diff_y = std::abs(apply_pbc(y1,IDX_X)-apply_pbc(y2,IDX_X)+yshift_coef*spatial_box_width[IDX_Y]);
Berenger Bramas's avatar
Berenger Bramas committed
235
        assert(diff_y <= 2*cutoff_radius);
236

237
        real_number diff_z = std::abs(apply_pbc(z1,IDX_X)-apply_pbc(z2,IDX_X)+zshift_coef*spatial_box_width[IDX_Z]);
Berenger Bramas's avatar
Berenger Bramas committed
238
        assert(diff_z <= 2*cutoff_radius);
Berenger Bramas's avatar
Debug    
Berenger Bramas committed
239
240
241

        return (diff_x*diff_x) + (diff_y*diff_y) + (diff_z*diff_z);
    }
242

243
    template <class computer_class, int size_particle_positions, int size_particle_rhs>
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
    void compute_distr(computer_class& in_computer,
                       const partsize_t current_my_nb_particles_per_partition[],
                       real_number particles_positions[],
                       real_number particles_current_rhs[],
                       partsize_t inout_index_particles[]){
        TIMEZONE("compute_distr");

        // Some processes might not be involved
        if(nb_processes_involved <= my_rank){
            return;
        }

        const long int my_top_z_cell_level = last_cell_level_proc(my_rank);
        const long int my_down_z_cell_level = first_cell_level_proc(my_rank);
        const long int my_nb_cell_levels = 1+my_top_z_cell_level-my_down_z_cell_level;

        current_offset_particles_for_partition[0] = 0;
        partsize_t myTotalNbParticles = 0;
        for(int idxPartition = 0 ; idxPartition < current_partition_size ; ++idxPartition){
            myTotalNbParticles += current_my_nb_particles_per_partition[idxPartition];
            current_offset_particles_for_partition[idxPartition+1] = current_offset_particles_for_partition[idxPartition] + current_my_nb_particles_per_partition[idxPartition];
        }

        // Compute box idx for each particle
        std::unique_ptr<long int[]> particles_coord(new long int[current_offset_particles_for_partition[current_partition_size]]);

        {
            for(int idxPartition = 0 ; idxPartition < current_partition_size ; ++idxPartition){
                #pragma omp parallel for schedule(static)
                for(partsize_t idxPart = current_offset_particles_for_partition[idxPartition] ; idxPart < current_offset_particles_for_partition[idxPartition+1] ; ++idxPart ){
                    particles_coord[idxPart] = get_cell_idx(particles_positions[(idxPart)*size_particle_positions + IDX_X],
                                                                              particles_positions[(idxPart)*size_particle_positions + IDX_Y],
                                                                              particles_positions[(idxPart)*size_particle_positions + IDX_Z]);
                    assert(my_down_z_cell_level <= get_cell_coord_z_from_index(particles_coord[idxPart]));
                    assert(get_cell_coord_z_from_index(particles_coord[idxPart]) <= my_top_z_cell_level);
                }
            }

282
            std::vector<std::pair<long int,partsize_t>> part_to_sort;
283
284
285
286
287
288
289

            // Sort each partition in cells
            for(int idxPartition = 0 ; idxPartition < current_partition_size ; ++idxPartition){
                part_to_sort.clear();

                for(partsize_t idxPart = current_offset_particles_for_partition[idxPartition] ; idxPart < current_offset_particles_for_partition[idxPartition+1] ; ++idxPart ){
                    part_to_sort.emplace_back();
290
291
                    part_to_sort.back().first = particles_coord[idxPart];
                    part_to_sort.back().second = idxPart;
292
293
                }

294
                assert(partsize_t(part_to_sort.size()) == (current_my_nb_particles_per_partition[idxPartition]));
295
296

                std::sort(part_to_sort.begin(), part_to_sort.end(),
297
298
299
                          [](const std::pair<long int,partsize_t>& p1,
                             const std::pair<long int,partsize_t>& p2){
                    return p1.first < p2.first;
300
                });
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315

                // Permute array using buffer
                std::vector<unsigned char> buffer;
                permute_copy<real_number, size_particle_positions>(current_offset_particles_for_partition[idxPartition],
                                                                   current_my_nb_particles_per_partition[idxPartition],
                                                                   part_to_sort.data(), particles_positions, &buffer);
                permute_copy<real_number, size_particle_rhs>(current_offset_particles_for_partition[idxPartition],
                                                             current_my_nb_particles_per_partition[idxPartition],
                                                             part_to_sort.data(), particles_current_rhs, &buffer);
                permute_copy<partsize_t, 1>(current_offset_particles_for_partition[idxPartition],
                                            current_my_nb_particles_per_partition[idxPartition],
                                            part_to_sort.data(), inout_index_particles, &buffer);
                permute_copy<long int, 1>(current_offset_particles_for_partition[idxPartition],
                                            current_my_nb_particles_per_partition[idxPartition],
                                            part_to_sort.data(), particles_coord.get(), &buffer);
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
            }
        }

        // Build the tree
        p2p_tree<std::vector<std::pair<partsize_t,partsize_t>>> my_tree(nb_cell_levels);

        for(int idxPartition = 0 ; idxPartition < current_partition_size ; ++idxPartition){
            long int current_cell_idx = -1;
            partsize_t current_nb_particles_in_cell = 0;
            partsize_t current_cell_offset = 0;

            for(partsize_t idx_part = current_offset_particles_for_partition[idxPartition] ;
                            idx_part != current_offset_particles_for_partition[idxPartition+1]; ++idx_part){
                if(particles_coord[idx_part] != current_cell_idx){
                    if(current_nb_particles_in_cell){
                        my_tree.getCell(current_cell_idx).emplace_back(current_cell_offset,current_nb_particles_in_cell);
                    }
                    current_cell_idx = particles_coord[idx_part];
                    current_nb_particles_in_cell = 1;
                    current_cell_offset = idx_part;
336
337
338
                }
                else{
                    current_nb_particles_in_cell += 1;
339
340
341
342
343
344
345
346
347
                }
            }
            if(current_nb_particles_in_cell){
                my_tree.getCell(current_cell_idx).emplace_back(current_cell_offset,current_nb_particles_in_cell);

            }
        }

        // Offset per cell layers
Berenger Bramas's avatar
Debug    
Berenger Bramas committed
348
        long int previous_index = 0;
349
350
351
352
        std::unique_ptr<partsize_t[]> particles_offset_layers(new partsize_t[my_nb_cell_levels+1]());
        for(int idxPartition = 0 ; idxPartition < current_partition_size ; ++idxPartition){
            for(partsize_t idx_part = current_offset_particles_for_partition[idxPartition] ;
                            idx_part != current_offset_particles_for_partition[idxPartition+1]; ++idx_part){
Berenger Bramas's avatar
Debug    
Berenger Bramas committed
353
354
355
356
357
358
                const long int part_box_z_index = get_cell_coord_z_from_index(particles_coord[idx_part]);
                assert(my_down_z_cell_level <= part_box_z_index);
                assert(part_box_z_index <= my_top_z_cell_level);
                particles_offset_layers[part_box_z_index+1-my_down_z_cell_level] += 1;
                assert(previous_index <= part_box_z_index);
                previous_index = part_box_z_index;
359
360
            }
        }
361
        for(long int idx_layer = 0 ; idx_layer < my_nb_cell_levels ; ++idx_layer){
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
            particles_offset_layers[idx_layer+1] += particles_offset_layers[idx_layer];
        }

        // Reset vectors
        assert(whatNext.size() == 0);
        assert(mpiRequests.size() == 0);
        neigDescriptors.clear();

        // Find process with at least one neighbor
        {
            int dest_proc = (my_rank+1)%nb_processes_involved;
            while(dest_proc != my_rank
                  && (my_top_z_cell_level == first_cell_level_proc(dest_proc)
                      || (my_top_z_cell_level+1)%nb_cell_levels[IDX_Z] == first_cell_level_proc(dest_proc))){
                // Find if we have to send 1 or 2 cell levels
                int nb_levels_to_send = 1;
                if(my_nb_cell_levels > 1 // I have more than one level
                        && (my_top_z_cell_level-1+2)%nb_cell_levels[IDX_Z] <= last_cell_level_proc(dest_proc)){
                    nb_levels_to_send += 1;
                }

                NeighborDescriptor descriptor;
                descriptor.destProc = dest_proc;
                descriptor.nbLevelsToExchange = nb_levels_to_send;
                descriptor.nbParticlesToExchange = particles_offset_layers[my_nb_cell_levels] - particles_offset_layers[my_nb_cell_levels-nb_levels_to_send];
                descriptor.isRecv = false;

                neigDescriptors.emplace_back(std::move(descriptor));

                dest_proc = (dest_proc+1)%nb_processes_involved;
            }

            int src_proc = (my_rank-1+nb_processes_involved)%nb_processes_involved;
            while(src_proc != my_rank
                  && (last_cell_level_proc(src_proc) == my_down_z_cell_level
                      || (last_cell_level_proc(src_proc)+1)%nb_cell_levels[IDX_Z] == my_down_z_cell_level)){
                // Find if we have to send 1 or 2 cell levels
                int nb_levels_to_recv = 1;
                if(my_nb_cell_levels > 1 // I have more than one level
                        && first_cell_level_proc(src_proc) <= (my_down_z_cell_level-1+2)%nb_cell_levels[IDX_Z]){
                    nb_levels_to_recv += 1;
                }

                NeighborDescriptor descriptor;
                descriptor.destProc = src_proc;
                descriptor.nbLevelsToExchange = nb_levels_to_recv;
                descriptor.nbParticlesToExchange = -1;
                descriptor.isRecv = true;

                neigDescriptors.emplace_back(std::move(descriptor));

                src_proc = (src_proc-1+nb_processes_involved)%nb_processes_involved;
            }
        }

        //////////////////////////////////////////////////////////////////////
        /// Exchange the number of particles in each partition
        /// Could involve only here but I do not think it will be a problem
        //////////////////////////////////////////////////////////////////////

        assert(whatNext.size() == 0);
        assert(mpiRequests.size() == 0);


        for(int idxDescr = 0 ; idxDescr < int(neigDescriptors.size()) ; ++idxDescr){
            NeighborDescriptor& descriptor = neigDescriptors[idxDescr];

            if(descriptor.isRecv == false){
                whatNext.emplace_back(std::pair<Action,int>{NOTHING_TODO, -1});
                mpiRequests.emplace_back();
                AssertMpi(MPI_Isend(const_cast<partsize_t*>(&descriptor.nbParticlesToExchange),
                                    1, particles_utils::GetMpiType(partsize_t()),
                                    descriptor.destProc, TAG_NB_PARTICLES,
                                    current_com, &mpiRequests.back()));

                if(descriptor.nbParticlesToExchange){
                    whatNext.emplace_back(std::pair<Action,int>{NOTHING_TODO, -1});
                    mpiRequests.emplace_back();
                    assert(descriptor.nbParticlesToExchange*size_particle_positions < std::numeric_limits<int>::max());
Berenger Bramas's avatar
Debug    
Berenger Bramas committed
441
                    AssertMpi(MPI_Isend(const_cast<real_number*>(&particles_positions[particles_offset_layers[my_nb_cell_levels-descriptor.nbLevelsToExchange]*size_particle_positions]),
442
443
444
445
446
447
                              int(descriptor.nbParticlesToExchange*size_particle_positions), particles_utils::GetMpiType(real_number()),
                              descriptor.destProc, TAG_POSITION_PARTICLES,
                              current_com, &mpiRequests.back()));

                    assert(descriptor.toRecvAndMerge == nullptr);
                    descriptor.toRecvAndMerge.reset(new real_number[descriptor.nbParticlesToExchange*size_particle_rhs]);
Berenger Bramas's avatar
Debug    
Berenger Bramas committed
448
                    whatNext.emplace_back(std::pair<Action,int>{MERGE_PARTICLES, idxDescr});
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
                    mpiRequests.emplace_back();
                    assert(descriptor.nbParticlesToExchange*size_particle_rhs < std::numeric_limits<int>::max());
                    AssertMpi(MPI_Irecv(descriptor.toRecvAndMerge.get(), int(descriptor.nbParticlesToExchange*size_particle_rhs),
                                        particles_utils::GetMpiType(real_number()), descriptor.destProc, TAG_RESULT_PARTICLES,
                                        current_com, &mpiRequests.back()));
                }
            }
            else{
                whatNext.emplace_back(std::pair<Action,int>{RECV_PARTICLES, idxDescr});
                mpiRequests.emplace_back();
                AssertMpi(MPI_Irecv(&descriptor.nbParticlesToExchange,
                      1, particles_utils::GetMpiType(partsize_t()), descriptor.destProc, TAG_NB_PARTICLES,
                      current_com, &mpiRequests.back()));
            }
        }

465
466
        lock_free_bool_array cells_locker(512);

467
468
469
470
471
472
        TIMEZONE_OMP_INIT_PREPARALLEL(omp_get_max_threads())
        #pragma omp parallel default(shared)
        {
            #pragma omp master
            {
                while(mpiRequests.size()){
Berenger Bramas's avatar
Berenger Bramas committed
473
                    TIMEZONE("wait-loop");
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
                    assert(mpiRequests.size() == whatNext.size());

                    int idxDone = int(mpiRequests.size());
                    {
                        TIMEZONE("wait");
                        AssertMpi(MPI_Waitany(int(mpiRequests.size()), mpiRequests.data(), &idxDone, MPI_STATUSES_IGNORE));
                    }
                    const std::pair<Action, int> releasedAction = whatNext[idxDone];
                    std::swap(mpiRequests[idxDone], mpiRequests[mpiRequests.size()-1]);
                    std::swap(whatNext[idxDone], whatNext[mpiRequests.size()-1]);
                    mpiRequests.pop_back();
                    whatNext.pop_back();

                    //////////////////////////////////////////////////////////////////////
                    /// Data to exchange particles
                    //////////////////////////////////////////////////////////////////////
                    if(releasedAction.first == RECV_PARTICLES){
Berenger Bramas's avatar
Berenger Bramas committed
491
                        TIMEZONE("post-recv-particles");
492
                        NeighborDescriptor& descriptor = neigDescriptors[releasedAction.second];
Berenger Bramas's avatar
Debug    
Berenger Bramas committed
493
                        assert(descriptor.isRecv);
494
495
496
497
498
499
500
                        const int destProc = descriptor.destProc;
                        const partsize_t NbParticlesToReceive = descriptor.nbParticlesToExchange;
                        assert(NbParticlesToReceive != -1);
                        assert(descriptor.toCompute == nullptr);

                        if(NbParticlesToReceive){
                            descriptor.toCompute.reset(new real_number[NbParticlesToReceive*size_particle_positions]);
501
                            whatNext.emplace_back(std::pair<Action,int>{COMPUTE_PARTICLES, releasedAction.second});
502
503
504
505
506
507
508
509
510
511
512
513
                            mpiRequests.emplace_back();
                            assert(NbParticlesToReceive*size_particle_positions < std::numeric_limits<int>::max());
                            AssertMpi(MPI_Irecv(descriptor.toCompute.get(), int(NbParticlesToReceive*size_particle_positions),
                                                particles_utils::GetMpiType(real_number()), destProc, TAG_POSITION_PARTICLES,
                                                current_com, &mpiRequests.back()));
                        }
                    }

                    //////////////////////////////////////////////////////////////////////
                    /// Computation
                    //////////////////////////////////////////////////////////////////////
                    if(releasedAction.first == COMPUTE_PARTICLES){
Berenger Bramas's avatar
Berenger Bramas committed
514
                        TIMEZONE("compute-particles");
515
                        NeighborDescriptor& descriptor = neigDescriptors[releasedAction.second];
Berenger Bramas's avatar
Debug    
Berenger Bramas committed
516
                        assert(descriptor.isRecv);
517
518
519
520
                        const partsize_t NbParticlesToReceive = descriptor.nbParticlesToExchange;

                        assert(descriptor.toCompute != nullptr);
                        descriptor.results.reset(new real_number[NbParticlesToReceive*size_particle_rhs]);
Berenger Bramas's avatar
Debug    
Berenger Bramas committed
521
                        in_computer.template init_result_array<size_particle_rhs>(descriptor.results.get(), NbParticlesToReceive);
522
523
524
525
526
527
528

                        // Compute
                        partsize_t idxPart = 0;
                        while(idxPart != NbParticlesToReceive){
                            const long int current_cell_idx = get_cell_idx(descriptor.toCompute[idxPart*size_particle_positions + IDX_X],
                                                                           descriptor.toCompute[idxPart*size_particle_positions + IDX_Y],
                                                                           descriptor.toCompute[idxPart*size_particle_positions + IDX_Z]);
Berenger Bramas's avatar
Debug    
Berenger Bramas committed
529
                            partsize_t nb_parts_in_cell = 1;
530
531
532
533
534
535
536
                            while(idxPart+nb_parts_in_cell != NbParticlesToReceive
                                  && current_cell_idx == get_cell_idx(descriptor.toCompute[(idxPart+nb_parts_in_cell)*size_particle_positions + IDX_X],
                                                                     descriptor.toCompute[(idxPart+nb_parts_in_cell)*size_particle_positions + IDX_Y],
                                                                     descriptor.toCompute[(idxPart+nb_parts_in_cell)*size_particle_positions + IDX_Z])){
                                nb_parts_in_cell += 1;
                            }

537
538
539
540
541
542
                            #pragma omp task default(shared) firstprivate(idxPart, nb_parts_in_cell, current_cell_idx)
                            {
                                const std::vector<std::pair<partsize_t,partsize_t>>* neighbors[27];
                                long int neighbors_indexes[27];
                                std::array<real_number,3> shift[27];
                                const int nbNeighbors = my_tree.getNeighbors(current_cell_idx, neighbors, neighbors_indexes, shift, true);
543

544
                                // with other interval
545
                                for(int idx_neighbor = 0 ; idx_neighbor < nbNeighbors ; ++idx_neighbor){
546
547
548
549
550
551
552
553
554
555
556
557
558
                                    cells_locker.lock(neighbors_indexes[idx_neighbor]);

                                    for(size_t idx_2 = 0 ; idx_2 < (*neighbors[idx_neighbor]).size() ; ++idx_2){
                                        for(partsize_t idx_p1 = 0 ; idx_p1 < nb_parts_in_cell ; ++idx_p1){
                                            for(partsize_t idx_p2 = 0 ; idx_p2 < (*neighbors[idx_neighbor])[idx_2].second ; ++idx_p2){
                                                const real_number dist_r2 = compute_distance_r2(descriptor.toCompute[(idxPart+idx_p1)*size_particle_positions + IDX_X],
                                                                                                descriptor.toCompute[(idxPart+idx_p1)*size_particle_positions + IDX_Y],
                                                                                                descriptor.toCompute[(idxPart+idx_p1)*size_particle_positions + IDX_Z],
                                                                                                particles_positions[((*neighbors[idx_neighbor])[idx_2].first+idx_p2)*size_particle_positions + IDX_X],
                                                                                                particles_positions[((*neighbors[idx_neighbor])[idx_2].first+idx_p2)*size_particle_positions + IDX_Y],
                                                                                                particles_positions[((*neighbors[idx_neighbor])[idx_2].first+idx_p2)*size_particle_positions + IDX_Z],
                                                                                                shift[idx_neighbor][IDX_X], shift[idx_neighbor][IDX_Y], shift[idx_neighbor][IDX_Z]);
                                                if(dist_r2 < cutoff_radius_compute*cutoff_radius_compute){
559
                                                    in_computer.template compute_interaction<size_particle_positions, size_particle_rhs>(
560
561
562
563
                                                                        &descriptor.toCompute[(idxPart+idx_p1)*size_particle_positions],
                                                                        &descriptor.results[(idxPart+idx_p1)*size_particle_rhs],
                                                                        &particles_positions[((*neighbors[idx_neighbor])[idx_2].first+idx_p2)*size_particle_positions],
                                                                        &particles_current_rhs[((*neighbors[idx_neighbor])[idx_2].first+idx_p2)*size_particle_rhs],
564
                                                                        dist_r2, cutoff_radius_compute, shift[idx_neighbor][IDX_X], shift[idx_neighbor][IDX_Y], shift[idx_neighbor][IDX_Z]);
565
                                                }
Berenger Bramas's avatar
Debug    
Berenger Bramas committed
566
                                            }
567
568
                                        }
                                    }
569
570

                                    cells_locker.unlock(neighbors_indexes[idx_neighbor]);
571
572
573
574
575
576
                                }
                            }

                            idxPart += nb_parts_in_cell;
                        }

577
578
                        #pragma omp taskwait

579
580
581
582
583
584
585
586
                        // Send back
                        const int destProc = descriptor.destProc;
                        whatNext.emplace_back(std::pair<Action,int>{RELEASE_BUFFER_PARTICLES, releasedAction.second});
                        mpiRequests.emplace_back();
                        assert(NbParticlesToReceive*size_particle_rhs < std::numeric_limits<int>::max());
                        AssertMpi(MPI_Isend(descriptor.results.get(), int(NbParticlesToReceive*size_particle_rhs),
                                            particles_utils::GetMpiType(real_number()), destProc, TAG_RESULT_PARTICLES,
                                            current_com, &mpiRequests.back()));
587
                        descriptor.toCompute.release();
588
589
                    }
                    //////////////////////////////////////////////////////////////////////
590
                    /// Release memory that was sent back
591
                    //////////////////////////////////////////////////////////////////////
592
593
                    if(releasedAction.first == RELEASE_BUFFER_PARTICLES){
                        NeighborDescriptor& descriptor = neigDescriptors[releasedAction.second];
594
                        assert(descriptor.results != nullptr);
Berenger Bramas's avatar
Debug    
Berenger Bramas committed
595
                        assert(descriptor.isRecv);
596
                        descriptor.results.release();
597
598
599
600
                    }
                    //////////////////////////////////////////////////////////////////////
                    /// Merge
                    //////////////////////////////////////////////////////////////////////
Berenger Bramas's avatar
Debug    
Berenger Bramas committed
601
                    if(releasedAction.first == MERGE_PARTICLES){
Berenger Bramas's avatar
Berenger Bramas committed
602
                        TIMEZONE("merge");
603
                        NeighborDescriptor& descriptor = neigDescriptors[releasedAction.second];
Berenger Bramas's avatar
Debug    
Berenger Bramas committed
604
                        assert(descriptor.isRecv == false);
605
                        assert(descriptor.toRecvAndMerge != nullptr);
Berenger Bramas's avatar
Debug    
Berenger Bramas committed
606
607
                        in_computer.template reduce_particles_rhs<size_particle_rhs>(&particles_current_rhs[particles_offset_layers[my_nb_cell_levels-descriptor.nbLevelsToExchange]*size_particle_rhs],
                                descriptor.toRecvAndMerge.get(), descriptor.nbParticlesToExchange);
608
609
610
611
612
613
614
615
616
617
618
                        descriptor.toRecvAndMerge.release();
                    }
                }
            }
        }

        assert(whatNext.size() == 0);
        assert(mpiRequests.size() == 0);

        // Compute self data
        for(const auto& iter_cell : my_tree){
Berenger Bramas's avatar
Berenger Bramas committed
619
            TIMEZONE("proceed-leaf");
620
621
            const long int currenct_cell_idx = iter_cell.first;
            const std::vector<std::pair<partsize_t,partsize_t>>* intervals_ptr = &iter_cell.second;
622

623
624
625
#pragma omp task default(shared) firstprivate(currenct_cell_idx, intervals_ptr)
            {
                const std::vector<std::pair<partsize_t,partsize_t>>& intervals = (*intervals_ptr);
626

627
                cells_locker.lock(currenct_cell_idx);
628

629
630
                for(size_t idx_1 = 0 ; idx_1 < intervals.size() ; ++idx_1){
                    // self interval
631
                    for(partsize_t idx_p1 = 0 ; idx_p1 < intervals[idx_1].second ; ++idx_p1){
632
                        for(partsize_t idx_p2 = idx_p1+1 ; idx_p2 < intervals[idx_1].second ; ++idx_p2){
633
634
635
                            const real_number dist_r2 = compute_distance_r2(particles_positions[(intervals[idx_1].first+idx_p1)*size_particle_positions + IDX_X],
                                                                            particles_positions[(intervals[idx_1].first+idx_p1)*size_particle_positions + IDX_Y],
                                                                            particles_positions[(intervals[idx_1].first+idx_p1)*size_particle_positions + IDX_Z],
636
637
638
                                                                            particles_positions[(intervals[idx_1].first+idx_p2)*size_particle_positions + IDX_X],
                                                                            particles_positions[(intervals[idx_1].first+idx_p2)*size_particle_positions + IDX_Y],
                                                                            particles_positions[(intervals[idx_1].first+idx_p2)*size_particle_positions + IDX_Z],
Berenger Bramas's avatar
Berenger Bramas committed
639
                                                                            0, 0, 0);
640
                            if(dist_r2 < cutoff_radius_compute*cutoff_radius_compute){
641
                                in_computer.template compute_interaction<size_particle_positions,size_particle_rhs>(
642
643
                                                    &particles_positions[(intervals[idx_1].first+idx_p1)*size_particle_positions],
                                                    &particles_current_rhs[(intervals[idx_1].first+idx_p1)*size_particle_rhs],
644
645
                                                    &particles_positions[(intervals[idx_1].first+idx_p2)*size_particle_positions],
                                                    &particles_current_rhs[(intervals[idx_1].first+idx_p2)*size_particle_rhs],
646
                                                    dist_r2, cutoff_radius_compute, 0, 0, 0);
647
                            }
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
                        }
                    }

                    // with other interval
                    for(size_t idx_2 = idx_1+1 ; idx_2 < intervals.size() ; ++idx_2){
                        for(partsize_t idx_p1 = 0 ; idx_p1 < intervals[idx_1].second ; ++idx_p1){
                            for(partsize_t idx_p2 = 0 ; idx_p2 < intervals[idx_2].second ; ++idx_p2){
                                const real_number dist_r2 = compute_distance_r2(particles_positions[(intervals[idx_1].first+idx_p1)*size_particle_positions + IDX_X],
                                                                                particles_positions[(intervals[idx_1].first+idx_p1)*size_particle_positions + IDX_Y],
                                                                                particles_positions[(intervals[idx_1].first+idx_p1)*size_particle_positions + IDX_Z],
                                                                                particles_positions[(intervals[idx_2].first+idx_p2)*size_particle_positions + IDX_X],
                                                                                particles_positions[(intervals[idx_2].first+idx_p2)*size_particle_positions + IDX_Y],
                                                                                particles_positions[(intervals[idx_2].first+idx_p2)*size_particle_positions + IDX_Z],
                                                                                0, 0, 0);
                                if(dist_r2 < cutoff_radius_compute*cutoff_radius_compute){
663
                                    in_computer.template compute_interaction<size_particle_positions,size_particle_rhs>(
664
665
666
667
                                                        &particles_positions[(intervals[idx_1].first+idx_p1)*size_particle_positions],
                                                        &particles_current_rhs[(intervals[idx_1].first+idx_p1)*size_particle_rhs],
                                                        &particles_positions[(intervals[idx_2].first+idx_p2)*size_particle_positions],
                                                        &particles_current_rhs[(intervals[idx_2].first+idx_p2)*size_particle_rhs],
668
                                                        dist_r2, cutoff_radius_compute, 0, 0, 0);
669
670
                                }
                            }
671
672
673
674
                        }
                    }
                }

675
676
677
678
679
680
681
                const std::vector<std::pair<partsize_t,partsize_t>>* neighbors[27];
                long int neighbors_indexes[27];
                std::array<real_number,3> shift[27];
                const int nbNeighbors = my_tree.getNeighbors(currenct_cell_idx, neighbors, neighbors_indexes, shift, false);

                for(size_t idx_1 = 0 ; idx_1 < intervals.size() ; ++idx_1){
                    // with other interval
682
                    for(int idx_neighbor = 0 ; idx_neighbor < nbNeighbors ; ++idx_neighbor){
683
684
685
686
687
688
689
690
691
692
693
694
695
696
                        if(currenct_cell_idx < neighbors_indexes[idx_neighbor]){
                            cells_locker.lock(neighbors_indexes[idx_neighbor]);

                            for(size_t idx_2 = 0 ; idx_2 < (*neighbors[idx_neighbor]).size() ; ++idx_2){
                                for(partsize_t idx_p1 = 0 ; idx_p1 < intervals[idx_1].second ; ++idx_p1){
                                    for(partsize_t idx_p2 = 0 ; idx_p2 < (*neighbors[idx_neighbor])[idx_2].second ; ++idx_p2){
                                        const real_number dist_r2 = compute_distance_r2(particles_positions[(intervals[idx_1].first+idx_p1)*size_particle_positions + IDX_X],
                                                                                        particles_positions[(intervals[idx_1].first+idx_p1)*size_particle_positions + IDX_Y],
                                                                                        particles_positions[(intervals[idx_1].first+idx_p1)*size_particle_positions + IDX_Z],
                                                                                        particles_positions[((*neighbors[idx_neighbor])[idx_2].first+idx_p2)*size_particle_positions + IDX_X],
                                                                                        particles_positions[((*neighbors[idx_neighbor])[idx_2].first+idx_p2)*size_particle_positions + IDX_Y],
                                                                                        particles_positions[((*neighbors[idx_neighbor])[idx_2].first+idx_p2)*size_particle_positions + IDX_Z],
                                                                                        shift[idx_neighbor][IDX_X], shift[idx_neighbor][IDX_Y], shift[idx_neighbor][IDX_Z]);
                                        if(dist_r2 < cutoff_radius_compute*cutoff_radius_compute){
697
                                            in_computer.template compute_interaction<size_particle_positions,size_particle_rhs>(
698
699
700
701
                                                                &particles_positions[(intervals[idx_1].first+idx_p1)*size_particle_positions],
                                                                &particles_current_rhs[(intervals[idx_1].first+idx_p1)*size_particle_rhs],
                                                                &particles_positions[((*neighbors[idx_neighbor])[idx_2].first+idx_p2)*size_particle_positions],
                                                                &particles_current_rhs[((*neighbors[idx_neighbor])[idx_2].first+idx_p2)*size_particle_rhs],
702
                                                                dist_r2, cutoff_radius_compute, shift[idx_neighbor][IDX_X], shift[idx_neighbor][IDX_Y], shift[idx_neighbor][IDX_Z]);
703
                                        }
704
                                    }
705
706
                                }
                            }
707
                            cells_locker.unlock(neighbors_indexes[idx_neighbor]);
708
709
710
                        }
                    }
                }
711
712

                cells_locker.unlock(currenct_cell_idx);
713
714
715
716
717
718
            }
        }
    }
};

#endif