1
/*============================================================================
2
* Functions related to the gathering of local arrays based on global ranks
4
*============================================================================*/
7
This file is part of Code_Saturne, a general-purpose CFD tool.
9
Copyright (C) 1998-2011 EDF S.A.
11
This program is free software; you can redistribute it and/or modify it under
12
the terms of the GNU General Public License as published by the Free Software
13
Foundation; either version 2 of the License, or (at your option) any later
16
This program is distributed in the hope that it will be useful, but WITHOUT
17
ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
18
FOR A PARTICULAR PURPOSE. See the GNU General Public License for more
21
You should have received a copy of the GNU General Public License along with
22
this program; if not, write to the Free Software Foundation, Inc., 51 Franklin
23
Street, Fifth Floor, Boston, MA 02110-1301, USA.
26
/*----------------------------------------------------------------------------*/
28
#if defined(HAVE_CONFIG_H)
29
#include "cs_config.h"
32
/*----------------------------------------------------------------------------
33
* Standard C library headers
34
*----------------------------------------------------------------------------*/
40
/*----------------------------------------------------------------------------
42
*----------------------------------------------------------------------------*/
45
#include <bft_printf.h>
47
/*----------------------------------------------------------------------------
49
*----------------------------------------------------------------------------*/
51
#include "fvm_config_defs.h"
53
#include "fvm_io_num.h"
54
#include "fvm_parall.h"
56
/*----------------------------------------------------------------------------
57
* Header for the current file
58
*----------------------------------------------------------------------------*/
60
#include "fvm_gather.h"
62
/*----------------------------------------------------------------------------*/
67
} /* Fake brace to force back Emacs auto-indentation back to column 0 */
69
#endif /* __cplusplus */
71
/*============================================================================
72
* Local structure definitions
73
*============================================================================*/
77
/*----------------------------------------------------------------------------
78
* Structure keeping track of the status of a series of fvm_gather_...()
79
* operations by slices of element global I/O number intervals.
80
*----------------------------------------------------------------------------*/
82
struct _fvm_gather_slice_t {
84
int local_rank; /* Local rank in communicator */
85
int n_ranks; /* Number of ranks in communicator */
87
/* Initial and final global element numbers for all slices */
89
fvm_gnum_t global_num_initial;
90
fvm_gnum_t global_num_final;
92
/* First and past the last global element numbers for current slice */
94
fvm_gnum_t ref_slice_size;
95
fvm_gnum_t global_num_slice_start;
96
fvm_gnum_t global_num_slice_end;
98
/* Local indexes corresponding to first element of a series with
99
global number >= global_num_slice_start initially, and
100
global number > global_num_end after a gather operation
101
(initialized to 0 before first of a series of calls) */
103
fvm_lnum_t local_index_start; /* Current start value */
104
fvm_lnum_t local_index_last; /* Last value reached, to advance */
106
/* Number of local entities (constant after initialization) */
108
fvm_lnum_t local_index_end;
110
/* Next global element number expected for each rank before and after
111
current call (rank 0 only); initialized upon the first call of
112
a series of fvm_gather_...() operations. Keeping track of this
113
allows us to avoid having to receive empty messages from
114
processors having no element in a given slice, which could
115
represent the majority of MPI calls for high processor counts */
117
fvm_gnum_t *next_global_num; /* Current values */
118
fvm_gnum_t *next_global_num_last; /* Next values, to advance */
120
/* We may only use the next global element number to avoid
121
sending and receiving empty messages when it is up to date;
122
this is not yet the case when a slice has just been initialized
123
or resized and no operations on it have been done yet */
125
_Bool use_next_global_num;
127
/* Buffers associated with MPI Datatype creation */
129
size_t recv_buf_size; /* Receive buffer size, in bytes */
130
void *recv_buf; /* Receive buffer */
132
int *blocklengths; /* Required only for indexed operations */
133
fvm_gnum_t *displacements; /* Always required for all ranks; size
134
slice_size + 1 (last position used to
135
exchange next_global_num[] values) */
139
#endif /* defined(HAVE_MPI) */
141
/*=============================================================================
142
* Private function definitions
143
*============================================================================*/
145
#if defined(HAVE_MPI)
147
/*----------------------------------------------------------------------------
148
* Check size of receive buffer for array data, resizing if necessary.
151
* this_slice <-- pointer to structure that should be destroyed
152
* n_entities <-- number of entities associated with receive
153
* stride <-- number of values per entity
154
* type size <-- sizeof(datatype)
155
*----------------------------------------------------------------------------*/
158
_slice_recv_buf_size_array(fvm_gather_slice_t * this_slice,
163
size_t size_mult = stride*type_size;
165
_Bool reallocate = false;
169
if (this_slice->recv_buf_size < this_slice->ref_slice_size*size_mult) {
170
this_slice->recv_buf_size = this_slice->ref_slice_size*size_mult;
174
/* Safety resizing (should not be necessary) */
176
while (this_slice->recv_buf_size < n_entities*size_mult) {
177
this_slice->recv_buf_size *= 2;
181
if (reallocate == true)
182
BFT_REALLOC(this_slice->recv_buf, this_slice->recv_buf_size, char);
185
/*----------------------------------------------------------------------------
186
* Check size of receive buffer for indexed data, resizing if necessary.
189
* this_slice <-- pointer to structure that should be destroyed
190
* n_values <-- number of values to receive
191
* type size <-- sizeof(datatype)
192
*----------------------------------------------------------------------------*/
195
_slice_recv_buf_size_indexed(fvm_gather_slice_t * this_slice,
199
size_t recv_size = n_values*type_size;
200
_Bool reallocate = false;
204
if ( this_slice->recv_buf_size < this_slice->ref_slice_size*type_size) {
205
this_slice->recv_buf_size = this_slice->ref_slice_size * type_size;
209
/* Safety resizing (may be necessary) */
211
while (this_slice->recv_buf_size < recv_size) {
212
this_slice->recv_buf_size *= 2;
216
if (reallocate == true)
217
BFT_REALLOC(this_slice->recv_buf, this_slice->recv_buf_size, char);
220
#if 0 && defined(DEBUG) && !defined(NDEBUG)
222
/*----------------------------------------------------------------------------
223
* Print information on a fvm_gather_slice_t structure.
226
* this_slice <-- pointer to structure that should be printed
227
*----------------------------------------------------------------------------*/
230
fvm_gather_slice_dump(fvm_gather_slice_t * this_slice)
232
if (this_slice != NULL) {
234
bft_printf("slice info:\n"
238
" global_num_initial = %llu\n"
239
" global_num_final = %llu\n"
240
" ref_slice_size = %llu\n"
241
" global_num_slice_start = %llu\n"
242
" global_num_slice_end = %llu\n"
243
" local_index_start = %ld\n"
244
" local_index_last = %ld\n"
245
" local_index_end = %ld\n"
246
" use_next_global_num = %d\n",
248
this_slice->local_rank,
250
(long)(this_slice->local_index_end),
251
(unsigned long long)(this_slice->global_num_initial),
252
(unsigned long long)(this_slice->global_num_final),
253
(unsigned long long)(this_slice->ref_slice_size),
254
(unsigned long long)(this_slice->global_num_slice_start),
255
(unsigned long long)(this_slice->global_num_slice_end),
256
(long)(this_slice->local_index_start),
257
(long)(this_slice->local_index_last),
258
(long)(this_slice->local_index_end),
259
(int)(this_slice->use_next_global_num));
261
if (this_slice->next_global_num != NULL) {
263
bft_printf(" rank | next_global_num | next_global_num_last\n");
264
for (i = 0; i < this_slice->n_ranks; i++) {
265
bft_printf(" %7d | %12llu | %12llu\n",
267
(unsigned long long)(this_slice->next_global_num[i]),
268
(unsigned long long)(this_slice->next_global_num_last[i]));
274
#endif /* 0 && defined(DEBUG) && !defined(NDEBUG) */
276
#endif /* defined(HAVE_MPI) */
278
/*=============================================================================
279
* Public function definitions
280
*============================================================================*/
282
#if defined(HAVE_MPI)
284
/*----------------------------------------------------------------------------
285
* Create a fvm_gather_slice_t structure.
288
* entity_io_num <-- I/O numbering structure associated with slice entity
289
* slice_size <-- reference slice size
290
* comm <-- associated MPI communicator
291
*----------------------------------------------------------------------------*/
294
fvm_gather_slice_create(const fvm_io_num_t *entity_io_num,
295
const fvm_gnum_t slice_size,
299
int local_rank, n_ranks;
300
fvm_gather_slice_t *this_slice;
302
/* Get local rank and size of the current MPI communicator */
303
MPI_Comm_rank(comm, &local_rank);
304
MPI_Comm_size(comm, &n_ranks);
306
/* Allocate and initialize slice structure */
308
BFT_MALLOC(this_slice, 1, fvm_gather_slice_t);
310
this_slice->local_rank = local_rank;
311
this_slice->n_ranks = n_ranks;
313
this_slice->global_num_initial = 1;
314
this_slice->global_num_final = fvm_io_num_get_global_count(entity_io_num);
316
this_slice->ref_slice_size = slice_size;
317
this_slice->global_num_slice_start = 1;
318
this_slice->global_num_slice_end = 1;
320
this_slice->local_index_end = fvm_io_num_get_local_count(entity_io_num);
322
this_slice->local_index_start = 0;
323
this_slice->local_index_last = 0;
325
/* Allocate and initialize "next expected global number" arrays */
327
if (local_rank == 0) {
329
BFT_MALLOC(this_slice->next_global_num, n_ranks, fvm_gnum_t);
330
BFT_MALLOC(this_slice->next_global_num_last, n_ranks, fvm_gnum_t);
332
for (i = 0; i < n_ranks; i++) {
333
this_slice->next_global_num[i] = 0;
334
this_slice->next_global_num_last[i] = 0;
339
this_slice->next_global_num = NULL;
340
this_slice->next_global_num_last = NULL;
343
this_slice->use_next_global_num = false;
345
/* Allocated buffers (blocklengths allocated only if needed for ranks > 0) */
347
this_slice->recv_buf_size = 0;
348
this_slice->recv_buf = NULL;
350
this_slice->blocklengths = NULL;
352
BFT_MALLOC(this_slice->displacements, slice_size + 1, fvm_gnum_t);
357
/*----------------------------------------------------------------------------
358
* Destroy a fvm_gather_slice_t structure.
361
* this_slice <-- pointer to structure that should be destroyed
365
*----------------------------------------------------------------------------*/
368
fvm_gather_slice_destroy(fvm_gather_slice_t * this_slice)
370
if (this_slice != NULL) {
372
if (this_slice->next_global_num != NULL) {
373
BFT_FREE(this_slice->next_global_num);
374
BFT_FREE(this_slice->next_global_num_last);
377
if (this_slice->recv_buf != NULL)
378
BFT_FREE(this_slice->recv_buf);
380
if (this_slice->blocklengths != NULL)
381
BFT_FREE(this_slice->blocklengths);
383
BFT_FREE(this_slice->displacements);
385
BFT_FREE(this_slice);
392
/*----------------------------------------------------------------------------
393
* Advance a fvm_gather_slice_t structure to the next start and end values.
395
* Elements within this slice will be those for whose global number
396
* is >= global_num_start and < global_num_end.
399
* this_slice <-- pointer to structure that should be advanced
400
* global_num_start --> new current global slice start number
401
* global_num_end --> new current global slice past the end number
404
* 0 if the end of the slice has not been reached before this call,
405
* 1 if we have already attained the end of the slice.
406
*----------------------------------------------------------------------------*/
409
fvm_gather_slice_advance(fvm_gather_slice_t *this_slice,
410
fvm_gnum_t *global_num_start,
411
fvm_gnum_t *global_num_end)
415
if (this_slice != NULL) {
417
if (this_slice->global_num_slice_end > this_slice->global_num_final)
420
this_slice->global_num_slice_start
421
= this_slice->global_num_slice_end;
423
this_slice->global_num_slice_end
424
= this_slice->global_num_slice_start + this_slice->ref_slice_size;
426
if ( this_slice->global_num_slice_end
427
> this_slice->global_num_final + 1)
428
this_slice->global_num_slice_end = this_slice->global_num_final + 1;
430
this_slice->local_index_start = this_slice->local_index_last;
432
if (this_slice->next_global_num != NULL) {
434
for (i = 0; i < this_slice->n_ranks; i++)
435
this_slice->next_global_num[i]
436
= this_slice->next_global_num_last[i];
439
if ( this_slice->global_num_slice_start
440
!= this_slice->global_num_initial)
441
this_slice->use_next_global_num = true;
443
*global_num_start = this_slice->global_num_slice_start;
444
*global_num_end = this_slice->global_num_slice_end;
451
/*----------------------------------------------------------------------------
452
* Reset a fvm_gather_slice_t structure to its initial state.
455
* this_slice <-- pointer to structure that should be reinitialized
456
*----------------------------------------------------------------------------*/
459
fvm_gather_slice_reinitialize(fvm_gather_slice_t *this_slice)
461
if (this_slice != NULL) {
463
this_slice->global_num_slice_start = this_slice->global_num_initial;
464
this_slice->global_num_slice_end = this_slice->global_num_initial;
466
this_slice->local_index_start = 0;
467
this_slice->local_index_last = 0;
469
if (this_slice->next_global_num != NULL) {
471
for (i = 0; i < this_slice->n_ranks; i++) {
472
this_slice->next_global_num[i] = 0;
473
this_slice->next_global_num_last[i] = 0;
477
this_slice->use_next_global_num = false;
481
/*----------------------------------------------------------------------------
482
* Limit an fvm_gather_slice_t structure's end value.
484
* This allows setting a lower global_num_end value than that previously
485
* defined (which may be necessary when buffer size limits require it).
488
* this_slice <-- pointer to structure that should be advanced
489
* global_num_end --> new current global slice past the end number
490
*----------------------------------------------------------------------------*/
493
fvm_gather_slice_limit(fvm_gather_slice_t *this_slice,
494
fvm_gnum_t *global_num_end)
496
if (*global_num_end != this_slice->global_num_slice_end) {
498
if (*global_num_end > this_slice->global_num_final)
499
*global_num_end = this_slice->global_num_final;
501
this_slice->global_num_slice_end = *global_num_end;
503
/* If slice size is changed, the next_global_num[] array of rank 0
504
may not be up to date in certain cases, so do not use it */
506
this_slice->use_next_global_num = false;
511
/*----------------------------------------------------------------------------
512
* Build a slice index (0 to n-1 numbering) on rank 0 from local index arrays.
514
* This is done by computing the local block lengths from the local
515
* index, gathering those lengths to rank 0, and rebuilding a 0 to n-1
516
* numbered slice index on rank 0.
518
* This function is intended to be used within a loop on subsets of the global
519
* lengths array (so as to enable writing to file or sending to an
520
* external process without requiring the full array to reside in the process
521
* directly handling I/O's memory). As such, it avoids allocating its own
522
* working arrays (buffers), so that they may be allocated outside the loop
523
* and reused for each call (avoiding the overhead associated with memory
526
* All or most elements in a given portion may belong to a same process rank
527
* (depending on mesh numbering and domain splitting). To account for
528
* this, for each process rank, the slice_index[] arrays must be large
529
* enough to contain (slice_size * stride) values, even though most processes
533
* local_index <-- local index array
534
* slice_index --> global slice index section for elements
535
* slice global_num_start to global_num_end
536
* (output for rank 0, working array only for others)
537
* element_io_num <-- I/O numbering structure associated with elements
538
* comm <-- MPI communicator for structures considered
539
* this_slice <-> structure for management of slice status
540
*----------------------------------------------------------------------------*/
543
fvm_gather_slice_index(const fvm_lnum_t local_index[],
544
fvm_gnum_t slice_index[],
545
const fvm_io_num_t *element_io_num,
547
fvm_gather_slice_t *this_slice)
550
int n_local_entities, n_distant_entities;
551
fvm_lnum_t local_index_start, local_index_stop;
553
/* MPI related variables */
556
const int local_rank = this_slice->local_rank;
557
const int n_ranks = this_slice->n_ranks;
558
fvm_gnum_t *const displacements = this_slice->displacements;
560
/* Local number of elements */
561
const fvm_lnum_t local_index_end = this_slice->local_index_end;
563
/* Global numbering */
564
const fvm_gnum_t global_num_start = this_slice->global_num_slice_start;
565
const fvm_gnum_t global_num_end = this_slice->global_num_slice_end;
566
const fvm_gnum_t *entity_global_num
567
= fvm_io_num_get_global_num(element_io_num);
569
/* Initialize blocklengths and displacements */
571
local_index_start = this_slice->local_index_start;
573
assert( local_index_start >= local_index_end
574
|| entity_global_num[local_index_start] >= global_num_start);
576
for (i = 0, j = local_index_start;
577
i < local_index_end && j < local_index_end
578
&& entity_global_num[j] < global_num_end;
580
displacements[i] = entity_global_num[j] - global_num_start;
583
n_local_entities = i;
584
local_index_stop = local_index_start + n_local_entities;
585
this_slice->local_index_last = local_index_stop;
587
if (local_index_stop < local_index_end)
588
displacements[n_local_entities] = entity_global_num[j];
590
displacements[n_local_entities] = this_slice->global_num_final + 1;
594
For rank 0, set "final" values directly;
595
values are lengths and not indexes at this stage.
598
if (local_rank == 0) {
599
for (i = 0, j = local_index_start;
600
i < n_local_entities;
602
slice_index[displacements[i]] = local_index[j+1] - local_index[j];
606
int n_local_values = n_local_entities;
607
for (i = 0, j = local_index_start;
610
slice_index[i] = local_index[j+1] - local_index[j];
614
/* Gather lengths information from ranks > 0 */
616
if (local_rank == 0) {
618
for (distant_rank = 1; distant_rank < n_ranks; distant_rank++) {
620
/* Get index from distant rank */
622
if ( this_slice->next_global_num[distant_rank] < global_num_end
623
|| this_slice->use_next_global_num == false) {
625
MPI_Send(&distant_rank, 1, MPI_INT, distant_rank, FVM_MPI_TAG, comm);
626
MPI_Recv(&n_distant_entities, 1, MPI_INT,
627
distant_rank, FVM_MPI_TAG, comm, &status);
629
MPI_Recv(displacements, n_distant_entities, FVM_MPI_GNUM,
630
distant_rank, FVM_MPI_TAG, comm, &status);
632
n_distant_entities -= 1;
633
this_slice->next_global_num_last[distant_rank]
634
= displacements[n_distant_entities];
636
if (n_distant_entities > 0) {
638
fvm_gnum_t *recv_buf = NULL;
640
_slice_recv_buf_size_array(this_slice, n_distant_entities,
641
1, sizeof(fvm_gnum_t));
643
recv_buf = this_slice->recv_buf;
645
MPI_Recv(this_slice->recv_buf, n_distant_entities,
646
FVM_MPI_GNUM, distant_rank, FVM_MPI_TAG, comm, &status);
648
for (i = 0 ; i < n_distant_entities ; i++)
649
slice_index[displacements[i]] = recv_buf[i];
657
/* Transform lengths to global index (0 to n-1)*/
662
int n_entities_max = (int)(global_num_end - global_num_start);
663
for (i = 0; i < n_entities_max; i++) {
664
l_cur = slice_index[i];
665
slice_index[i] = l_sum;
668
slice_index[n_entities_max] = l_sum;
674
if ( n_local_entities > 0
675
|| this_slice->use_next_global_num == false) {
677
/* Send local index to rank zero */
681
MPI_Recv(&buf_val, 1, MPI_INT, 0, FVM_MPI_TAG, comm, &status);
683
buf_val = n_local_entities + 1;
684
MPI_Send(&buf_val, 1, MPI_INT, 0, FVM_MPI_TAG, comm);
686
MPI_Send(displacements, n_local_entities + 1,
687
FVM_MPI_GNUM, 0, FVM_MPI_TAG, comm);
689
/* Send local portion of lengths array to rank zero */
691
if (n_local_entities > 0)
692
MPI_Send(slice_index, (int)n_local_entities,
693
FVM_MPI_GNUM, 0, FVM_MPI_TAG, comm);
699
#if 0 && defined(DEBUG) && !defined(NDEBUG)
705
/*----------------------------------------------------------------------------
706
* Recompute maximum value of global_num_end and slice connectivity size for
707
* an indexed connectivity slice.
709
* Given an initial global connectivity buffer size associated with the
710
* slice (global_connect_s_size), this function verifies that the connectivity
711
* associated with the slice from global_num_start to global_num_end may fit
712
* in this buffer. If this is not the case, global_num_end is reduced to the
713
* largest value such that the associated indexed connectivity or values may
714
* fit in the indicated buffer size.
716
* In any case, slice size will neither be increased above the current
717
* slice size, nor be reduced to less than
718
* min(n_g_elements, n_elements_s_min) if initially larger than this.
719
* If necessary, global_connect_s_size is increased so that this minimal
720
* slice may fit in a buffer of this same size.
723
* n_elements_s_min <-- minimum number of elements per slice desired
724
* global_num_end --> new current global slice past the end number
725
* global_connect_s_size <-> pointer to global connectivity slice size
726
* comm <-- associated MPI communicator
727
* slice_index <-- index of blocks corresponding to a given
728
* element in the global_connect_s array
729
* (required for rank 0 only)
730
* this_slice <-> structure for management of slice status
731
*----------------------------------------------------------------------------*/
734
fvm_gather_resize_indexed_slice(const fvm_gnum_t n_elements_s_min,
735
fvm_gnum_t *global_num_end,
736
fvm_gnum_t *global_connect_s_size,
738
const fvm_gnum_t slice_index[],
739
fvm_gather_slice_t *this_slice)
744
fvm_gnum_t global_num_start = this_slice->global_num_slice_start;
746
const int local_rank = this_slice->local_rank;
748
*global_num_end = this_slice->global_num_slice_end;
750
/* Recompute maximum value of global_num_end for this slice */
752
if (local_rank == 0) {
754
for (i_s = 1; i_s < *global_num_end - global_num_start + 1; i_s++) {
755
if (slice_index[i_s] > *global_connect_s_size) {
756
*global_num_end = global_num_start + i_s - 1;
761
/* If possible, size to at least n_elements_s_min (unless
762
reference slice size is smaller or already at end of slice) */
764
if (*global_num_end - global_num_start < n_elements_s_min) {
766
*global_num_end = global_num_start + n_elements_s_min;
768
if (*global_num_end - global_num_start > this_slice->ref_slice_size)
769
*global_num_end = global_num_start + this_slice->ref_slice_size;
771
if (*global_num_end > this_slice->global_num_final + 1)
772
*global_num_end = this_slice->global_num_final + 1;
774
/* Limit to initial global_num_end */
776
if (*global_num_end > this_slice->global_num_slice_end)
777
*global_num_end = this_slice->global_num_slice_end;
779
*global_connect_s_size =
780
FVM_MAX(*global_connect_s_size,
781
slice_index[*global_num_end - global_num_start]);
785
buf[0] = *global_num_end, buf[1] = *global_connect_s_size;
789
MPI_Bcast(buf, 2, FVM_MPI_GNUM, 0, comm);
791
/* Modify (reduce) slice size limit if necessary */
793
fvm_gather_slice_limit(this_slice, &(buf[0]));
795
/* Set output values */
797
*global_num_end = buf[0];
798
*global_connect_s_size = buf[1];
800
#if 0 && defined(DEBUG) && !defined(NDEBUG)
806
/*----------------------------------------------------------------------------
807
* Gather a given portion of an array to rank 0.
809
* This function is intended to be used within a loop on subsets of the global
810
* array (so as to enable writing to file or sending to an external process
811
* without requiring the full array to reside in the process directly
812
* handling I/O's memory). As such, it avoids allocating its own working arrays
813
* (buffers), so that they may be allocated outside the loop and reused for
814
* each call (avoiding the overhead associated with memory allocation).
816
* All or most elements in a given portion may belong to a same process rank
817
* (depending on mesh numbering and domain splitting). To account for
818
* this, for each process rank, the global_array_s[] array must be large
819
* enough to contain (slice_size * stride) values, even though most processes
823
* local_array <-- local array (size n_local_elements * stride)
824
* global_array_s --> global array section for elements
825
* slice global_num_start to global_num_end
826
* (output for rank 0, working array only for others)
827
* datatype <-- MPI datatype of each value
828
* stride <-- number of (interlaced) values per element
829
* element_io_num <-- I/O numbering structure associated with elements
830
* comm <-- MPI communicator for structures considered
831
* this_slice <-> structure for management of slice status
832
*----------------------------------------------------------------------------*/
835
fvm_gather_array(const void *local_array,
836
void *global_array_s,
837
MPI_Datatype datatype,
839
const fvm_io_num_t *element_io_num,
841
fvm_gather_slice_t *this_slice)
843
int n_local_entities, n_distant_entities;
846
fvm_lnum_t local_index_start, local_index_stop;
848
const char *local_array_val = local_array;
849
char *global_array_s_val = global_array_s;
851
/* MPI related variables */
853
int size, distant_rank;
854
const int local_rank = this_slice->local_rank;
855
const int n_ranks = this_slice->n_ranks;
856
fvm_gnum_t *const displacements = this_slice->displacements;
858
/* Local number of elements */
859
const fvm_lnum_t local_index_end = this_slice->local_index_end;
861
/* Global numbering */
862
const fvm_gnum_t global_num_start = this_slice->global_num_slice_start;
863
const fvm_gnum_t global_num_end = this_slice->global_num_slice_end;
864
const fvm_gnum_t *entity_global_num
865
= fvm_io_num_get_global_num(element_io_num);
867
/* Get info on the current MPI communicator */
869
MPI_Type_size(datatype, &size);
871
/* Initialize displacements */
873
local_index_start = this_slice->local_index_start;
875
assert( local_index_start >= local_index_end
876
|| entity_global_num[local_index_start] >= global_num_start);
878
/* Displacements should be expressed in bytes */
880
size_mult = size * stride;
882
for (i = 0, j = local_index_start;
883
( i < (size_t)local_index_end
884
&& j < (size_t)local_index_end
885
&& entity_global_num[j] < global_num_end);
887
displacements[i] = (entity_global_num[j] - global_num_start) * size_mult;
890
n_local_entities = i;
891
local_index_stop = local_index_start + n_local_entities;
892
this_slice->local_index_last = local_index_stop;
894
if (local_index_stop < local_index_end)
895
displacements[n_local_entities] = entity_global_num[j];
897
displacements[n_local_entities] = this_slice->global_num_final + 1;
900
Prepare send buffer (we use a copy to ensure constedness of input)
901
For rank 0, set final values directly.
904
if (local_rank == 0) {
905
for (i = 0, j = (size_t)local_index_start;
906
i < (size_t)n_local_entities;
908
for (k = 0; k < size_mult; k++) {
909
global_array_s_val[displacements[i] + k]
910
= local_array_val[j*size_mult + k];
915
memcpy(global_array_s_val,
916
local_array_val+(local_index_start*size_mult),
917
n_local_entities*size_mult);
919
/* Gather connectivity information from ranks > 0 */
921
if (local_rank == 0) {
923
for (distant_rank = 1; distant_rank < n_ranks; distant_rank++) {
925
/* Get index from distant rank */
927
if ( this_slice->next_global_num[distant_rank] < global_num_end
928
|| this_slice->use_next_global_num == false) {
930
MPI_Send(&distant_rank, 1, MPI_INT, distant_rank, FVM_MPI_TAG, comm);
931
MPI_Recv(&n_distant_entities, 1, MPI_INT,
932
distant_rank, FVM_MPI_TAG, comm, &status);
934
MPI_Recv(displacements, n_distant_entities, FVM_MPI_GNUM,
935
distant_rank, FVM_MPI_TAG, comm, &status);
937
n_distant_entities -= 1;
938
this_slice->next_global_num_last[distant_rank]
939
= displacements[n_distant_entities];
941
if (n_distant_entities > 0) {
943
char *_recv_buf = NULL;
945
_slice_recv_buf_size_array(this_slice, n_distant_entities,
948
_recv_buf = this_slice->recv_buf;
950
MPI_Recv(this_slice->recv_buf, n_distant_entities*stride, datatype,
951
distant_rank, FVM_MPI_TAG, comm, &status);
953
for (i = 0 ; i < (size_t)n_distant_entities ; i++) {
954
for (j = 0 ; j < size_mult ; j++) {
955
global_array_s_val[displacements[i] + j]
956
= _recv_buf[i*size_mult + j];
969
if ( n_local_entities > 0
970
|| this_slice->use_next_global_num ==false) {
972
/* Send local index to rank zero */
976
MPI_Recv(&buf_val, 1, MPI_INT, 0, FVM_MPI_TAG, comm, &status);
978
buf_val = n_local_entities + 1;
979
MPI_Send(&buf_val, 1, MPI_INT, 0, FVM_MPI_TAG, comm);
981
MPI_Send(displacements, n_local_entities + 1, FVM_MPI_GNUM, 0,
984
/* Send local portion of connectivity array to rank zero */
986
if (n_local_entities > 0)
987
MPI_Send(global_array_s, (int)(n_local_entities * stride),
988
datatype, 0, FVM_MPI_TAG, comm);
994
#if 0 && defined(DEBUG) && !defined(NDEBUG)
1000
/*----------------------------------------------------------------------------
1001
* Gather a given portion of an indexed array of to rank 0.
1003
* A slice_index[] array indicating the index (0 to n-1) of blocks in
1004
* the slice is required for rank 0. This implies that the block sizes in
1005
* the slice have already been gathered through the use of
1006
* fvm_gather_slice_index() or some similar method, and used to build this
1009
* This function is intended to be used within a loop on subsets of the global
1010
* lengths array (so as to enable writing to file or sending to an
1011
* external process without requiring the full array to reside in the process
1012
* directly handling I/O's memory). As such, it avoids allocating its own
1013
* working arrays (buffers), so that they may be allocated outside the loop
1014
* and reused for each call (avoiding the overhead associated with memory
1017
* All or most elements in a given portion may belong to a same process rank
1018
* (depending on mesh numbering and domain splitting). To account for
1019
* this, for each process rank, the global_lengths_s[] arrays must be large
1020
* enough to contain (slice_index[current_slice_size] - 1) values, even
1021
* though most processes will require less.
1022
* Use fvm_gather_resize_indexed_slice() to adjust current_slice_size.
1025
* local_array <-- local array
1026
* (size: local_index[n_local_elements] * stride)
1027
* global_array_s --> global array section for elements
1028
* slice global_num_start to global_num_end
1029
* (output for rank 0, working array only for others)
1030
* datatype <-- MPI datatype of each value
1031
* local_index <-- local index array
1032
* element_io_num <-- I/O numbering structure associated with elements
1033
* comm <-- MPI communicator for structures considered
1034
* slice_index <-- index of blocks corresponding to a given
1035
* element in the global_numbers_s array
1036
* (required for rank 0 only)
1037
* this_slice <-> structure for management of slice status
1038
*----------------------------------------------------------------------------*/
1041
fvm_gather_indexed(const void *local_array,
1042
void *global_array_s,
1043
const MPI_Datatype datatype,
1044
const fvm_lnum_t local_index[],
1045
const fvm_io_num_t *element_io_num,
1047
const fvm_gnum_t slice_index[],
1048
fvm_gather_slice_t *this_slice)
1051
int n_local_entities, n_distant_entities;
1052
int n_values_send = 0;
1053
fvm_lnum_t local_index_start, local_index_stop;
1055
const char *local_array_val = local_array;
1056
char *global_array_s_val = global_array_s;
1058
/* MPI related variables */
1060
int size, distant_rank;
1061
const int local_rank = this_slice->local_rank;
1062
const int n_ranks = this_slice->n_ranks;
1063
int *blocklengths = this_slice->blocklengths;
1064
fvm_gnum_t *const displacements = this_slice->displacements;
1066
/* Local number of elements */
1067
const fvm_lnum_t local_index_end = this_slice->local_index_end;
1069
/* Global numbering */
1070
const fvm_gnum_t global_num_start = this_slice->global_num_slice_start;
1071
const fvm_gnum_t global_num_end = this_slice->global_num_slice_end;
1072
const fvm_gnum_t *entity_global_num
1073
= fvm_io_num_get_global_num(element_io_num);
1075
MPI_Type_size(datatype, &size);
1077
/* Initialize blocklengths and displacements; this operation
1078
requires a blocklengths[] array for all ranks (contrary
1079
to others), so it updates the slice structure accordingly */
1081
if (blocklengths == NULL) {
1082
BFT_MALLOC(this_slice->blocklengths, this_slice->ref_slice_size, int);
1083
blocklengths = this_slice->blocklengths;
1086
local_index_start = this_slice->local_index_start;
1088
assert( local_index_start >= local_index_end
1089
|| entity_global_num[local_index_start] >= global_num_start);
1091
/* Displacements are first used to transfer the global slice index position
1092
of a given entity, and will be set to the index value later on rank 0 */
1094
for (i = 0, j = local_index_start;
1095
j < local_index_end && entity_global_num[j] < global_num_end;
1097
displacements[i] = entity_global_num[j] - global_num_start;
1099
n_local_entities = i;
1100
local_index_stop = local_index_start + n_local_entities;
1101
this_slice->local_index_last = local_index_stop;
1103
if (local_index_stop < local_index_end)
1104
displacements[n_local_entities] = entity_global_num[j];
1106
displacements[n_local_entities] = this_slice->global_num_final + 1;
1109
Prepare send buffer:
1110
For rank 0, set final values directly; for others, set blocklengths
1111
to be sent to rank 0.
1114
if (local_rank == 0) {
1115
for (i = 0, j = local_index_start;
1116
i < n_local_entities;
1118
const size_t displacement = slice_index[displacements[i]] * size;
1119
for (k = local_index[j]*(int)size, l = 0;
1120
k < (local_index[j+1])*(int)size;
1122
global_array_s_val[displacement + l] = local_array_val[k];
1126
n_values_send = ( local_index[local_index_start + n_local_entities]
1127
- local_index[local_index_start]);
1128
memcpy(global_array_s_val,
1129
local_array_val + ((local_index[local_index_start]) * size),
1130
n_values_send * size);
1131
for (i = 0, j = local_index_start;
1132
i < n_local_entities;
1134
blocklengths[i] = (local_index[j+1] - local_index[j]);
1138
/* Gather numbers information from ranks > 0 */
1140
if (local_rank == 0) {
1142
for (distant_rank = 1; distant_rank < n_ranks; distant_rank++) {
1144
size_t recv_size = 0;
1146
/* Get index from distant rank */
1148
if ( this_slice->next_global_num[distant_rank] < global_num_end
1149
|| this_slice->use_next_global_num == false) {
1151
MPI_Send(&distant_rank, 1, MPI_INT, distant_rank, FVM_MPI_TAG, comm);
1152
MPI_Recv(&n_distant_entities, 1, MPI_INT,
1153
distant_rank, FVM_MPI_TAG, comm, &status);
1155
/* Get index from distant rank */
1157
MPI_Recv(displacements, n_distant_entities, FVM_MPI_GNUM,
1158
distant_rank, FVM_MPI_TAG, comm, &status);
1160
n_distant_entities -= 1;
1161
this_slice->next_global_num_last[distant_rank]
1162
= displacements[n_distant_entities];
1164
for (i = 0; i < n_distant_entities; i++) {
1165
j = displacements[i];
1166
blocklengths[i] = (slice_index[j+1] - slice_index[j])*size;
1167
displacements[i] = slice_index[j] * size;
1168
recv_size += blocklengths[i];
1171
if (n_distant_entities > 0) {
1174
char *_recv_buf = NULL;
1176
_slice_recv_buf_size_indexed(this_slice,
1180
MPI_Recv(this_slice->recv_buf, (int)recv_size, datatype,
1181
distant_rank, FVM_MPI_TAG, comm, &status);
1183
_recv_buf = this_slice->recv_buf;
1185
for (i = 0, recv_id = 0 ; i < n_distant_entities ; i++) {
1186
for (j = 0 ; j < blocklengths[i] ; j++) {
1187
global_array_s_val[displacements[i] + j]
1188
= _recv_buf[recv_id++];
1201
if ( n_local_entities > 0
1202
|| this_slice->use_next_global_num == false) {
1204
/* Send local index to rank zero */
1208
MPI_Recv(&buf_val, 1, MPI_INT, 0, FVM_MPI_TAG, comm, &status);
1210
buf_val = n_local_entities + 1;
1211
MPI_Send(&buf_val, 1, MPI_INT, 0, FVM_MPI_TAG, comm);
1213
MPI_Send(displacements, n_local_entities + 1, FVM_MPI_GNUM,
1214
0, FVM_MPI_TAG, comm);
1216
/* Send local portion of numbers array to rank zero */
1218
if (n_local_entities > 0)
1219
MPI_Send(global_array_s, (int)n_values_send,
1220
datatype, 0, FVM_MPI_TAG, comm);
1226
#if 0 && defined(DEBUG) && !defined(NDEBUG)
1232
/*----------------------------------------------------------------------------
1233
* Gather a given portion of a strided (i.e. regular) connectivity array
1234
* to rank 0. Connectivity values are converted from local to global values
1235
* (both with 1 to n type numbering).
1237
* This function is intended to be used within a loop on subsets of the global
1238
* connectivity array (so as to enable writing to file or sending to an
1239
* external process without requiring the full array to reside in the process
1240
* directly handling I/O's memory). As such, it avoids allocating its own
1241
* working arrays (buffers), so that they may be allocated outside the loop
1242
* and reused for each call (avoiding the overhead associated with memory
1245
* All or most elements in a given portion may belong to a same process rank
1246
* (depending on mesh numbering and domain splitting). To account for
1247
* this, for each process rank, the global_connect_s[] array must be large
1248
* enough to contain (slice_size * stride) values, even though most processes
1249
* will require less.
1252
* local_connect <-- local connectivity array (1 to n numbering)
1253
* global_connect_s --> global connectivity array section for elements
1254
* slice global_num_start to global_num_end
1255
* (output for rank 0, working array only for others)
1256
* stride <-- number of connected entities (i.e. vertices in
1257
* a nodal connectivity) per element
1258
* connected_io_num <-- I/O numbering structure associated with "connected"
1259
* entities (i.e. vertices in a nodal connectivity)
1260
* element_io_num <-- I/O numbering structure associated with elements
1261
* comm <-- MPI communicator for structures considered
1262
* this_slice <-> structure for management of slice status
1263
*----------------------------------------------------------------------------*/
1266
fvm_gather_strided_connect(const fvm_lnum_t local_connect[],
1267
fvm_gnum_t global_connect_s[],
1269
const fvm_io_num_t *connected_io_num,
1270
const fvm_io_num_t *element_io_num,
1272
fvm_gather_slice_t *this_slice)
1275
int n_local_entities, n_distant_entities;
1276
fvm_lnum_t local_index_start, local_index_stop;
1278
/* MPI related variables */
1281
const int local_rank = this_slice->local_rank;
1282
const int n_ranks = this_slice->n_ranks;
1283
fvm_gnum_t *const displacements = this_slice->displacements;
1285
/* Local number of elements */
1286
const fvm_lnum_t local_index_end = this_slice->local_index_end;
1288
/* Global numbering */
1289
const fvm_gnum_t global_num_start = this_slice->global_num_slice_start;
1290
const fvm_gnum_t global_num_end = this_slice->global_num_slice_end;
1291
const fvm_gnum_t *connected_global_num
1292
= fvm_io_num_get_global_num(connected_io_num);
1293
const fvm_gnum_t *entity_global_num
1294
= fvm_io_num_get_global_num(element_io_num);
1296
/* Initialize displacements */
1298
local_index_start = this_slice->local_index_start;
1300
assert( local_index_start >= local_index_end
1301
|| entity_global_num[local_index_start] >= global_num_start);
1303
for (i = 0, j = local_index_start;
1304
i < local_index_end && j < local_index_end
1305
&& entity_global_num[j] < global_num_end;
1307
displacements[i] = (entity_global_num[j] - global_num_start) * stride;
1310
n_local_entities = i;
1311
local_index_stop = local_index_start + n_local_entities;
1312
this_slice->local_index_last = local_index_stop;
1314
if (local_index_stop < local_index_end)
1315
displacements[n_local_entities] = entity_global_num[j];
1317
displacements[n_local_entities] = this_slice->global_num_final + 1;
1320
Prepare send buffer:
1321
replace local connected entity numbers with their global counterparts.
1322
For rank 0, set final values directly.
1325
if (local_rank == 0) {
1326
for (i = 0, j = local_index_start;
1327
i < n_local_entities;
1329
for (k = 0; k < stride; k++)
1330
global_connect_s[displacements[i] + k]
1331
= connected_global_num[local_connect[j*stride + k] - 1];
1335
int n_local_values = n_local_entities * stride;
1336
for (i = 0, j = (fvm_gnum_t)(local_index_start * stride);
1339
global_connect_s[i] = connected_global_num[local_connect[j] - 1];
1343
/* Gather connectivity information from ranks > 0 */
1345
if (local_rank == 0) {
1347
for (distant_rank = 1; distant_rank < n_ranks; distant_rank++) {
1349
/* Get index from distant rank */
1351
if ( this_slice->next_global_num[distant_rank] < global_num_end
1352
|| this_slice->use_next_global_num == false) {
1354
MPI_Send(&distant_rank, 1, MPI_INT, distant_rank, FVM_MPI_TAG, comm);
1355
MPI_Recv(&n_distant_entities, 1, MPI_INT,
1356
distant_rank, FVM_MPI_TAG, comm, &status);
1358
MPI_Recv(displacements, n_distant_entities, FVM_MPI_GNUM,
1359
distant_rank, FVM_MPI_TAG, comm, &status);
1361
n_distant_entities -= 1;
1362
this_slice->next_global_num_last[distant_rank]
1363
= displacements[n_distant_entities];
1365
if (n_distant_entities > 0) {
1367
fvm_gnum_t *_recv_buf;
1369
_slice_recv_buf_size_array(this_slice, n_distant_entities,
1370
stride, sizeof(fvm_gnum_t));
1372
_recv_buf = this_slice->recv_buf;
1374
MPI_Recv(this_slice->recv_buf, (int)(n_distant_entities*stride),
1375
FVM_MPI_GNUM, distant_rank, FVM_MPI_TAG, comm, &status);
1377
for (i = 0 ; i < n_distant_entities ; i++) {
1378
for (j = 0 ; j < stride ; j++)
1379
global_connect_s[displacements[i] + j]
1380
= _recv_buf[i*stride + j];
1392
if ( n_local_entities > 0
1393
|| this_slice->use_next_global_num == false) {
1395
/* Send local index to rank zero */
1399
MPI_Recv(&buf_val, 1, MPI_INT, 0, FVM_MPI_TAG, comm, &status);
1401
buf_val = n_local_entities + 1;
1402
MPI_Send(&buf_val, 1, MPI_INT, 0, FVM_MPI_TAG, comm);
1404
MPI_Send(displacements, n_local_entities + 1, FVM_MPI_GNUM, 0,
1407
/* Send local portion of connectivity array to rank zero */
1409
if (n_local_entities > 0)
1410
MPI_Send(global_connect_s, (int)(n_local_entities * stride),
1411
FVM_MPI_GNUM, 0, FVM_MPI_TAG, comm);
1417
#if 0 && defined(DEBUG) && !defined(NDEBUG)
1423
/*----------------------------------------------------------------------------
1424
* Gather a given portion of an indexed array of numbers to rank 0.
1425
* If the connected_io_num argument is non-NULL, these numbers
1426
* are assumed to represent connectivity values, and are converted from
1427
* local to global values (both with 1 to n type numbering).
1428
* Otherwise, they are considered to represent any other type of positive
1429
* integer (such as the number of vertices for each of a polyhedron's faces).
1431
* A slice_index[] array indicating the index (0 to n-1) of blocks in
1432
* the slice is required for rank 0. This implies that the block sizes in
1433
* the slice have already been gathered through the use of
1434
* fvm_gather_slice_index() or some similar method, and used to build this
1437
* This function is intended to be used within a loop on subsets of the global
1438
* lengths array (so as to enable writing to file or sending to an
1439
* external process without requiring the full array to reside in the process
1440
* directly handling I/O's memory). As such, it avoids allocating its own
1441
* working arrays (buffers), so that they may be allocated outside the loop
1442
* and reused for each call (avoiding the overhead associated with memory
1445
* All or most elements in a given portion may belong to a same process rank
1446
* (depending on mesh numbering and domain splitting). To account for
1447
* this, for each process rank, the global_lengths_s[] arrays must be large
1448
* enough to contain (slice_index[current_slice_size] - 1) values, even
1449
* though most processes will require less.
1450
* Use fvm_gather_resize_indexed_slice() to adjust current_slice_size.
1453
* local_index <-- local index array
1454
* local_numbers <-- local numbers array
1455
* global_numbers_s --> global numbers array section for elements
1456
* slice global_num_start to global_num_end
1457
* (output for rank 0, working array only for others)
1458
* connected_io_num <-- I/O numbering structure associated with "connected"
1459
* entities (i.e. vertices in a nodal connectivity)
1460
* element_io_num <-- I/O numbering structure associated with elements
1461
* comm <-- MPI communicator for structures considered
1462
* slice_index <-- index of blocks corresponding to a given
1463
* element in the global_numbers_s array
1464
* (required for rank 0 only)
1465
* this_slice <-> structure for management of slice status
1466
*----------------------------------------------------------------------------*/
1469
fvm_gather_indexed_numbers(const fvm_lnum_t local_index[],
1470
const fvm_lnum_t local_numbers[],
1471
fvm_gnum_t global_numbers_s[],
1472
const fvm_io_num_t *connected_io_num,
1473
const fvm_io_num_t *element_io_num,
1475
const fvm_gnum_t slice_index[],
1476
fvm_gather_slice_t *this_slice)
1479
int n_local_entities, n_distant_entities;
1480
int n_values_send = 0;
1481
fvm_lnum_t local_index_start, local_index_stop;
1483
/* MPI related variables */
1486
const int local_rank = this_slice->local_rank;
1487
const int n_ranks = this_slice->n_ranks;
1488
int *blocklengths = this_slice->blocklengths;
1489
fvm_gnum_t *const displacements = this_slice->displacements;
1491
/* Local number of elements */
1492
const fvm_lnum_t local_index_end = this_slice->local_index_end;
1494
/* Global numbering */
1495
const fvm_gnum_t global_num_start = this_slice->global_num_slice_start;
1496
const fvm_gnum_t global_num_end = this_slice->global_num_slice_end;
1497
const fvm_gnum_t *connected_global_num = NULL;
1498
const fvm_gnum_t *entity_global_num
1499
= fvm_io_num_get_global_num(element_io_num);
1501
if (connected_io_num != NULL)
1502
connected_global_num = fvm_io_num_get_global_num(connected_io_num);
1504
/* Initialize blocklengths and displacements; this operation
1505
requires a blocklengths[] array for all ranks (contrary
1506
to others), so it updates the slice structure accordingly */
1508
if (blocklengths == NULL) {
1509
BFT_MALLOC(this_slice->blocklengths, this_slice->ref_slice_size, int);
1510
blocklengths = this_slice->blocklengths;
1513
local_index_start = this_slice->local_index_start;
1515
assert( local_index_start >= local_index_end
1516
|| entity_global_num[local_index_start] >= global_num_start);
1518
/* Displacements are first used to transfer the global slice index position
1519
of a given entity, and will be set to the index value later on rank 0 */
1521
for (i = 0, j = local_index_start;
1522
i < local_index_end && j < local_index_end
1523
&& entity_global_num[j] < global_num_end;
1525
displacements[i] = entity_global_num[j] - global_num_start;
1527
n_local_entities = i;
1528
local_index_stop = local_index_start + n_local_entities;
1529
this_slice->local_index_last = local_index_stop;
1531
if (local_index_stop < local_index_end)
1532
displacements[n_local_entities] = entity_global_num[j];
1534
displacements[n_local_entities] = this_slice->global_num_final + 1;
1537
Prepare send buffer:
1538
For rank 0, set final values directly; for others, set blocklengths
1539
to be sent to rank 0.
1542
if (connected_io_num == NULL) {
1544
if (local_rank == 0) {
1545
for (i = 0, j = local_index_start;
1546
i < n_local_entities;
1548
displacements[i] = slice_index[displacements[i]];
1549
for (k = local_index[j], l = 0; k < local_index[j+1]; k++, l++)
1550
global_numbers_s[displacements[i] + l] = local_numbers[k];
1555
for (i = 0, j = local_index_start;
1556
i < n_local_entities;
1558
blocklengths[i] = local_index[j+1] - local_index[j];
1559
for (k = local_index[j]; k < local_index[j+1]; k++)
1560
global_numbers_s[l++] = local_numbers[k];
1568
if (local_rank == 0) {
1569
for (i = 0, j = local_index_start;
1570
i < n_local_entities;
1572
displacements[i] = slice_index[displacements[i]];
1573
for (k = local_index[j], l = 0; k < local_index[j+1]; k++, l++)
1574
global_numbers_s[displacements[i] + l]
1575
= connected_global_num[local_numbers[k] - 1];
1580
for (i = 0, j = local_index_start;
1581
i < n_local_entities;
1583
blocklengths[i] = local_index[j+1] - local_index[j];
1584
for (k = local_index[j]; k < local_index[j+1]; k++)
1585
global_numbers_s[l++] = connected_global_num[local_numbers[k] - 1];
1592
/* Gather numbers information from ranks > 0 */
1594
if (local_rank == 0) {
1596
for (distant_rank = 1; distant_rank < n_ranks; distant_rank++) {
1598
/* Get index from distant rank */
1600
if ( this_slice->next_global_num[distant_rank] < global_num_end
1601
|| this_slice->use_next_global_num == false) {
1603
size_t recv_size = 0;
1605
MPI_Send(&distant_rank, 1, MPI_INT, distant_rank, FVM_MPI_TAG, comm);
1606
MPI_Recv(&n_distant_entities, 1, MPI_INT,
1607
distant_rank, FVM_MPI_TAG, comm, &status);
1609
/* Get index from distant rank */
1611
MPI_Recv(displacements, n_distant_entities, FVM_MPI_GNUM,
1612
distant_rank, FVM_MPI_TAG, comm, &status);
1614
n_distant_entities -= 1;
1615
this_slice->next_global_num_last[distant_rank]
1616
= displacements[n_distant_entities];
1618
for (i = 0; i < n_distant_entities; i++) {
1619
j = displacements[i];
1620
blocklengths[i] = slice_index[j+1] - slice_index[j];
1621
displacements[i] = slice_index[j];
1622
recv_size += blocklengths[i];
1625
if (n_distant_entities > 0) {
1628
fvm_gnum_t *_recv_buf = NULL;
1630
_slice_recv_buf_size_indexed(this_slice,
1632
sizeof(fvm_gnum_t));
1634
MPI_Recv(this_slice->recv_buf, recv_size, FVM_MPI_GNUM,
1635
distant_rank, FVM_MPI_TAG, comm, &status);
1637
_recv_buf = this_slice->recv_buf;
1639
for (i = 0, recv_id = 0 ; i < n_distant_entities ; i++) {
1640
for (j = 0 ; j < blocklengths[i] ; j++)
1641
global_numbers_s[displacements[i] + j]
1642
= _recv_buf[recv_id++];
1654
if ( n_local_entities > 0
1655
|| this_slice->use_next_global_num == false) {
1657
/* Send local index to rank zero */
1661
MPI_Recv(&buf_val, 1, MPI_INT, 0, FVM_MPI_TAG, comm, &status);
1663
buf_val = n_local_entities + 1;
1664
MPI_Send(&buf_val, 1, MPI_INT, 0, FVM_MPI_TAG, comm);
1666
MPI_Send(displacements, n_local_entities + 1, FVM_MPI_GNUM,
1667
0, FVM_MPI_TAG, comm);
1669
/* Send local portion of numbers array to rank zero */
1671
if (n_local_entities > 0)
1672
MPI_Send(global_numbers_s, (int)n_values_send,
1673
FVM_MPI_GNUM, 0, FVM_MPI_TAG, comm);
1679
#if 0 && defined(DEBUG) && !defined(NDEBUG)
1685
#endif /* defined(HAVE_MPI) */
1687
/*----------------------------------------------------------------------------*/
1691
#endif /* __cplusplus */