This file is indexed.

/usr/include/dolfin/common/MPI.h is in libdolfin-dev 2017.2.0.post0-2.

This file is owned by root:root, with mode 0o644.

The actual contents of the file can be viewed below.

  1
  2
  3
  4
  5
  6
  7
  8
  9
 10
 11
 12
 13
 14
 15
 16
 17
 18
 19
 20
 21
 22
 23
 24
 25
 26
 27
 28
 29
 30
 31
 32
 33
 34
 35
 36
 37
 38
 39
 40
 41
 42
 43
 44
 45
 46
 47
 48
 49
 50
 51
 52
 53
 54
 55
 56
 57
 58
 59
 60
 61
 62
 63
 64
 65
 66
 67
 68
 69
 70
 71
 72
 73
 74
 75
 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
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
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
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
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
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
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
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
747
748
749
750
751
752
753
754
755
756
757
758
759
760
761
762
763
764
765
766
767
768
769
770
771
772
773
774
775
776
777
778
779
780
781
782
783
784
785
786
787
788
789
790
791
792
793
794
795
796
797
798
799
800
801
802
803
804
805
806
807
808
809
810
811
812
813
814
815
816
817
818
819
820
821
822
823
824
825
826
827
828
829
830
831
832
833
834
835
836
837
838
839
840
841
842
843
844
845
846
847
848
849
850
851
852
853
854
855
856
857
858
859
860
861
862
863
864
865
866
867
868
869
870
871
872
873
874
875
876
877
878
879
880
881
882
883
884
885
886
887
888
889
890
// Copyright (C) 2007-2014 Magnus Vikstrøm and Garth N. Wells
//
// This file is part of DOLFIN.
//
// DOLFIN is free software: you can redistribute it and/or modify
// it under the terms of the GNU Lesser General Public License as published by
// the Free Software Foundation, either version 3 of the License, or
// (at your option) any later version.
//
// DOLFIN is distributed in the hope that it will be useful,
// but WITHOUT ANY WARRANTY; without even the implied warranty of
// MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
// GNU Lesser General Public License for more details.
//
// You should have received a copy of the GNU Lesser General Public License
// along with DOLFIN. If not, see <http://www.gnu.org/licenses/>.
//
// Modified by Ola Skavhaug 2008-2009
// Modified by Anders Logg 2008-2011
// Modified by Niclas Jansson 2009
// Modified by Joachim B Haga 2012
// Modified by Martin Sandve Alnes 2014

#ifndef __MPI_DOLFIN_WRAPPER_H
#define __MPI_DOLFIN_WRAPPER_H

#include <cstdint>
#include <iostream>
#include <numeric>
#include <type_traits>
#include <utility>
#include <vector>

#ifdef HAS_MPI
#define MPICH_IGNORE_CXX_SEEK 1
#include <mpi.h>
#endif

#include <dolfin/log/log.h>
#include <dolfin/log/Table.h>

#ifndef HAS_MPI
typedef int MPI_Comm;
#define MPI_COMM_WORLD 2
#define MPI_COMM_SELF 1
#define MPI_COMM_NULL 0
#endif

namespace dolfin
{

  #ifdef HAS_MPI
  class MPIInfo
  {
  public:
    MPIInfo();
    ~MPIInfo();
    MPI_Info& operator*();
  private:

    MPI_Info info;

  };
  #endif

  /// This class provides utility functions for easy communication
  /// with MPI and handles cases when DOLFIN is not configured with
  /// MPI.

  class MPI
  {
  public:

    // Create a duplicate MPI communicator and manage lifetime of the
    // communicator
    class Comm
    {
    public:

      /// Duplicate communicator and wrap duplicate
      Comm(MPI_Comm comm);

      // The disabling of the below is turned off because the SWIG
      // docstring generator fails on it.

      // Disable default constructor
      //Comm() = default;

      // Disable copy constructor
      //Comm(const Comm& comm) = delete;

      // Disable assignment operator
      //void operator=(Comm const &comm) = delete;

      /// Destructor (frees wrapped communicator)
      ~Comm();

      /// Free (destroy) communicator. Calls function 'MPI_Comm_free'.
      void free();

      /// Duplicate communivator, and free any previously created
      /// communicator
      void reset(MPI_Comm comm);

      /// Return process rank for the communicator
      unsigned int rank() const;

      /// Return size of the group (number of processes) associated
      /// with the communicator. This function will also intialise MPI
      /// if it hasn't already been intialised.
      unsigned int size() const;

      /// Set a barrier (synchronization point)
      void barrier() const;

      /// Return the underlying MPI_Comm object
      MPI_Comm comm() const;

    private:

      // MPI communicator
      MPI_Comm _comm;
    };

    /// Return process rank for the communicator
    static unsigned int rank(MPI_Comm comm);

    /// Return size of the group (number of processes) associated with
    /// the communicator
    static unsigned int size(MPI_Comm comm);

    /// Determine whether we should broadcast (based on current
    /// parallel policy)
    static bool is_broadcaster(MPI_Comm comm);

    /// Determine whether we should receive (based on current parallel
    /// policy)
    static bool is_receiver(MPI_Comm comm);

    /// Set a barrier (synchronization point)
    static void barrier(MPI_Comm comm);

    /// Send in_values[p0] to process p0 and receive values from
    /// process p1 in out_values[p1]
    template<typename T>
      static void all_to_all(MPI_Comm comm,
                             std::vector<std::vector<T>>& in_values,
                             std::vector<std::vector<T>>& out_values);


    /// Send in_values[p0] to process p0 and receive values from
    /// all processes in out_values
    template<typename T>
      static void all_to_all(MPI_Comm comm,
                             std::vector<std::vector<T>>& in_values,
                             std::vector<T>& out_values);

    /// Broadcast vector of value from broadcaster to all processes
    template<typename T>
      static void broadcast(MPI_Comm comm, std::vector<T>& value,
                            unsigned int broadcaster=0);

    /// Broadcast single primitive from broadcaster to all processes
    template<typename T>
      static void broadcast(MPI_Comm comm, T& value,
                            unsigned int broadcaster=0);

    /// Scatter vector in_values[i] to process i
    template<typename T>
      static void scatter(MPI_Comm comm,
                          const std::vector<std::vector<T>>& in_values,
                          std::vector<T>& out_value,
                          unsigned int sending_process=0);

    /// Scatter primitive in_values[i] to process i
    template<typename T>
      static void scatter(MPI_Comm comm,
                          const std::vector<T>& in_values,
                          T& out_value, unsigned int sending_process=0);

    /// Gather values on one process
    template<typename T>
      static void gather(MPI_Comm comm, const std::vector<T>& in_values,
                         std::vector<T>& out_values,
                         unsigned int receiving_process=0);

    /// Gather strings on one process
    static void gather(MPI_Comm comm, const std::string& in_values,
                       std::vector<std::string>& out_values,
                       unsigned int receiving_process=0);

    /// Gather values from all processes. Same data count from each
    /// process (wrapper for MPI_Allgather)
    template<typename T>
      static void all_gather(MPI_Comm comm,
                             const std::vector<T>& in_values,
                             std::vector<T>& out_values);

    /// Gather values from each process (variable count per process)
    template<typename T>
      static void all_gather(MPI_Comm comm,
                             const std::vector<T>& in_values,
                             std::vector<std::vector<T>>& out_values);

    /// Gather values, one primitive from each process (MPI_Allgather)
    template<typename T>
      static void all_gather(MPI_Comm comm, const T in_value,
                             std::vector<T>& out_values);

    /// Gather values, one primitive from each process (MPI_Allgather).
    /// Specialization for std::string
    static void all_gather(MPI_Comm comm, const std::string& in_values,
                           std::vector<std::string>& out_values);

    /// Return global max value
    template<typename T> static T max(MPI_Comm comm, const T& value);

    /// Return global min value
    template<typename T> static T min(MPI_Comm comm, const T& value);

    /// Sum values and return sum
    template<typename T> static T sum(MPI_Comm comm, const T& value);

    /// Return average across comm; implemented only for T == Table
    template<typename T> static T avg(MPI_Comm comm, const T& value);

    /// All reduce
    template<typename T, typename X> static
      T all_reduce(MPI_Comm comm, const T& value, X op);

    /// Find global offset (index) (wrapper for MPI_(Ex)Scan with
    /// MPI_SUM as reduction op)
    static std::size_t global_offset(MPI_Comm comm,
                                     std::size_t range, bool exclusive);

    /// Send-receive data between processes (blocking)
    template<typename T>
      static void send_recv(MPI_Comm comm,
                            const std::vector<T>& send_value,
                            unsigned int dest, int send_tag,
                            std::vector<T>& recv_value,
                            unsigned int source, int recv_tag);

    /// Send-receive data between processes
    template<typename T>
      static void send_recv(MPI_Comm comm,
                            const std::vector<T>& send_value, unsigned int dest,
                            std::vector<T>& recv_value, unsigned int source);

    /// Return local range for local process, splitting [0, N - 1] into
    /// size() portions of almost equal size
    static std::pair<std::int64_t, std::int64_t>
      local_range(MPI_Comm comm, std::int64_t N);

    /// Return local range for given process, splitting [0, N - 1] into
    /// size() portions of almost equal size
    static std::pair<std::int64_t, std::int64_t>
      local_range(MPI_Comm comm, int process, std::int64_t N);

    /// Return local range for given process, splitting [0, N - 1] into
    /// size() portions of almost equal size
    static std::pair<std::int64_t, std::int64_t>
      compute_local_range(int process, std::int64_t N, int size);

    /// Return which process owns index (inverse of local_range)
    static unsigned int index_owner(MPI_Comm comm,
                                    std::size_t index, std::size_t N);

    #ifdef HAS_MPI
    /// Return average reduction operation; recognized by
    /// all_reduce(MPI_Comm, Table&, MPI_Op)
    static MPI_Op MPI_AVG();
    #endif

  private:

    #ifndef HAS_MPI
    static void error_no_mpi(const char *where)
    {
      dolfin_error("MPI.h",
                   where,
                   "DOLFIN has been configured without MPI support");
    }
    #endif

    #ifdef HAS_MPI
    // Return MPI data type
    template<typename T>
      struct dependent_false : std::false_type {};
    template<typename T> static MPI_Datatype mpi_type()
    {
      static_assert(dependent_false<T>::value, "Unknown MPI type");
      dolfin_error("MPI.h",
                   "perform MPI operation",
                   "MPI data type unknown");
      return MPI_CHAR;
    }
    #endif

    #ifdef HAS_MPI
    // Maps some MPI_Op values to string
    static std::map<MPI_Op, std::string> operation_map;
    #endif

  };

  #ifdef HAS_MPI
  // Specialisations for MPI_Datatypes
  template<> inline MPI_Datatype MPI::mpi_type<float>() { return MPI_FLOAT; }
  template<> inline MPI_Datatype MPI::mpi_type<double>() { return MPI_DOUBLE; }
  template<> inline MPI_Datatype MPI::mpi_type<short int>() { return MPI_SHORT; }
  template<> inline MPI_Datatype MPI::mpi_type<int>() { return MPI_INT; }
  template<> inline MPI_Datatype MPI::mpi_type<long int>() { return MPI_LONG; }
  template<> inline MPI_Datatype MPI::mpi_type<unsigned int>() { return MPI_UNSIGNED; }
  template<> inline MPI_Datatype MPI::mpi_type<unsigned long int>() { return MPI_UNSIGNED_LONG; }
  template<> inline MPI_Datatype MPI::mpi_type<long long>() { return MPI_LONG_LONG; }
  #endif
  //---------------------------------------------------------------------------
  template<typename T>
    void dolfin::MPI::broadcast(MPI_Comm comm, std::vector<T>& value,
                                unsigned int broadcaster)
  {
    #ifdef HAS_MPI
    // Broadcast cast size
    std::size_t bsize = value.size();
    MPI_Bcast(&bsize, 1, mpi_type<std::size_t>(), broadcaster, comm);

    // Broadcast
    value.resize(bsize);
    MPI_Bcast(const_cast<T*>(value.data()), bsize, mpi_type<T>(),
              broadcaster, comm);
    #endif
  }
  //---------------------------------------------------------------------------
  template<typename T>
    void dolfin::MPI::broadcast(MPI_Comm comm, T& value,
                                unsigned int broadcaster)
  {
    #ifdef HAS_MPI
    MPI_Bcast(&value, 1, mpi_type<T>(), broadcaster, comm);
    #endif
  }
  //---------------------------------------------------------------------------
  template<typename T>
    void dolfin::MPI::all_to_all(MPI_Comm comm,
                                 std::vector<std::vector<T>>& in_values,
                                 std::vector<std::vector<T>>& out_values)
  {
    #ifdef HAS_MPI
    const std::size_t comm_size = MPI::size(comm);

    // Data size per destination
    dolfin_assert(in_values.size() == comm_size);
    std::vector<int> data_size_send(comm_size);
    std::vector<int> data_offset_send(comm_size + 1, 0);
    for (std::size_t p = 0; p < comm_size; ++p)
    {
      data_size_send[p] = in_values[p].size();
      data_offset_send[p + 1] = data_offset_send[p] + data_size_send[p];
    }

    // Get received data sizes
    std::vector<int> data_size_recv(comm_size);
    MPI_Alltoall(data_size_send.data(), 1, mpi_type<int>(),
                 data_size_recv.data(), 1, mpi_type<int>(),
                 comm);

    // Pack data and build receive offset
    std::vector<int> data_offset_recv(comm_size + 1 ,0);
    std::vector<T> data_send(data_offset_send[comm_size]);
    for (std::size_t p = 0; p < comm_size; ++p)
    {
      data_offset_recv[p + 1] = data_offset_recv[p] + data_size_recv[p];
      std::copy(in_values[p].begin(), in_values[p].end(),
                data_send.begin() + data_offset_send[p]);
    }

    // Send/receive data
    std::vector<T> data_recv(data_offset_recv[comm_size]);
    MPI_Alltoallv(data_send.data(), data_size_send.data(),
                  data_offset_send.data(), mpi_type<T>(),
                  data_recv.data(), data_size_recv.data(),
                  data_offset_recv.data(), mpi_type<T>(), comm);

      // Repack data
    out_values.resize(comm_size);
    for (std::size_t p = 0; p < comm_size; ++p)
    {
      out_values[p].resize(data_size_recv[p]);
      std::copy(data_recv.begin() + data_offset_recv[p],
                data_recv.begin() + data_offset_recv[p + 1],
                out_values[p].begin());
    }
    #else
    dolfin_assert(in_values.size() == 1);
    out_values = in_values;
    #endif
  }
  //---------------------------------------------------------------------------
  template<typename T>
    void dolfin::MPI::all_to_all(MPI_Comm comm,
                                 std::vector<std::vector<T>>& in_values,
                                 std::vector<T>& out_values)
  {
    #ifdef HAS_MPI
    const std::size_t comm_size = MPI::size(comm);

    // Data size per destination
    dolfin_assert(in_values.size() == comm_size);
    std::vector<int> data_size_send(comm_size);
    std::vector<int> data_offset_send(comm_size + 1, 0);
    for (std::size_t p = 0; p < comm_size; ++p)
    {
      data_size_send[p] = in_values[p].size();
      data_offset_send[p + 1] = data_offset_send[p] + data_size_send[p];
    }

    // Get received data sizes
    std::vector<int> data_size_recv(comm_size);
    MPI_Alltoall(data_size_send.data(), 1, mpi_type<int>(),
                 data_size_recv.data(), 1, mpi_type<int>(),
                 comm);

    // Pack data and build receive offset
    std::vector<int> data_offset_recv(comm_size + 1 ,0);
    std::vector<T> data_send(data_offset_send[comm_size]);
    for (std::size_t p = 0; p < comm_size; ++p)
    {
      data_offset_recv[p + 1] = data_offset_recv[p] + data_size_recv[p];
      std::copy(in_values[p].begin(), in_values[p].end(),
                data_send.begin() + data_offset_send[p]);
    }

    // Send/receive data
    out_values.resize(data_offset_recv[comm_size]);
    MPI_Alltoallv(data_send.data(), data_size_send.data(),
                  data_offset_send.data(), mpi_type<T>(),
                  out_values.data(), data_size_recv.data(),
                  data_offset_recv.data(), mpi_type<T>(), comm);

    #else
    dolfin_assert(in_values.size() == 1);
    out_values = in_values[0];
    #endif
  }
  //---------------------------------------------------------------------------
#ifndef DOXYGEN_IGNORE
  template<> inline
    void dolfin::MPI::all_to_all(MPI_Comm comm,
                                 std::vector<std::vector<bool>>& in_values,
                                 std::vector<std::vector<bool>>& out_values)
  {
    #ifdef HAS_MPI
    // Copy to short int
    std::vector<std::vector<short int>> send(in_values.size());
    for (std::size_t i = 0; i < in_values.size(); ++i)
      send[i].assign(in_values[i].begin(), in_values[i].end());

    // Communicate data
    std::vector<std::vector<short int>> recv;
    all_to_all(comm, send, recv);

    // Copy back to bool
    out_values.resize(recv.size());
    for (std::size_t i = 0; i < recv.size(); ++i)
      out_values[i].assign(recv[i].begin(), recv[i].end());
    #else
    dolfin_assert(in_values.size() == 1);
    out_values = in_values;
    #endif
  }

  template<> inline
    void dolfin::MPI::all_to_all(MPI_Comm comm,
                                 std::vector<std::vector<bool>>& in_values,
                                 std::vector<bool>& out_values)
  {
    #ifdef HAS_MPI
    // Copy to short int
    std::vector<std::vector<short int>> send(in_values.size());
    for (std::size_t i = 0; i < in_values.size(); ++i)
      send[i].assign(in_values[i].begin(), in_values[i].end());

    // Communicate data
    std::vector<short int> recv;
    all_to_all(comm, send, recv);

    // Copy back to bool
    out_values.assign(recv.begin(), recv.end());
    #else
    dolfin_assert(in_values.size() == 1);
    out_values = in_values[0];
    #endif
  }

#endif
  //---------------------------------------------------------------------------
  template<typename T>
    void dolfin::MPI::scatter(MPI_Comm comm,
                              const std::vector<std::vector<T>>& in_values,
                              std::vector<T>& out_value,
                              unsigned int sending_process)
  {
    #ifdef HAS_MPI

    // Scatter number of values to each process
    const std::size_t comm_size = MPI::size(comm);
    std::vector<int> all_num_values;
    if (MPI::rank(comm) == sending_process)
    {
      dolfin_assert(in_values.size() == comm_size);
      all_num_values.resize(comm_size);
      for (std::size_t i = 0; i < comm_size; ++i)
        all_num_values[i] = in_values[i].size();
    }
    int my_num_values = 0;
    scatter(comm, all_num_values, my_num_values, sending_process);

    // Prepare send buffer and offsets
    std::vector<T> sendbuf;
    std::vector<int> offsets;
    if (MPI::rank(comm) == sending_process)
    {
      // Build offsets
      offsets.resize(comm_size + 1, 0);
      for (std::size_t i = 1; i <= comm_size; ++i)
        offsets[i] = offsets[i - 1] + all_num_values[i - 1];

      // Allocate send buffer and fill
      const std::size_t n = std::accumulate(all_num_values.begin(),
                                            all_num_values.end(), 0);
      sendbuf.resize(n);
      for (std::size_t p = 0; p < in_values.size(); ++p)
      {
        std::copy(in_values[p].begin(), in_values[p].end(),
                  sendbuf.begin() + offsets[p]);
        }
    }

    // Scatter
    out_value.resize(my_num_values);
    MPI_Scatterv(const_cast<T*>(sendbuf.data()), all_num_values.data(),
                 offsets.data(), mpi_type<T>(),
                 out_value.data(), my_num_values,
                 mpi_type<T>(), sending_process, comm);
    #else
    dolfin_assert(sending_process == 0);
    dolfin_assert(in_values.size() == 1);
    out_value = in_values[0];
    #endif
  }
  //---------------------------------------------------------------------------
#ifndef DOXYGEN_IGNORE
  template<> inline
    void dolfin::MPI::scatter(MPI_Comm comm,
                              const std::vector<std::vector<bool>>& in_values,
                              std::vector<bool>& out_value,
                              unsigned int sending_process)
  {
    #ifdef HAS_MPI
    // Copy data
    std::vector<std::vector<short int>> in(in_values.size());
    for (std::size_t i = 0; i < in_values.size(); ++i)
      in[i] = std::vector<short int>(in_values[i].begin(), in_values[i].end());

    std::vector<short int> out;
    scatter(comm, in, out, sending_process);

    out_value.resize(out.size());
    std::copy(out.begin(), out.end(), out_value.begin());
    #else
    dolfin_assert(sending_process == 0);
    dolfin_assert(in_values.size() == 1);
    out_value = in_values[0];
    #endif
  }
#endif
  //---------------------------------------------------------------------------
  template<typename T>
    void dolfin::MPI::scatter(MPI_Comm comm,
                              const std::vector<T>& in_values,
                              T& out_value, unsigned int sending_process)
  {
    #ifdef HAS_MPI
    if (MPI::rank(comm) == sending_process)
      dolfin_assert(in_values.size() == MPI::size(comm));

    MPI_Scatter(const_cast<T*>(in_values.data()), 1, mpi_type<T>(),
                 &out_value, 1, mpi_type<T>(), sending_process, comm);
    #else
    dolfin_assert(sending_process == 0);
    dolfin_assert(in_values.size() == 1);
    out_value = in_values[0];
    #endif
  }
  //---------------------------------------------------------------------------
  template<typename T>
  void dolfin::MPI::gather(MPI_Comm comm,
                           const std::vector<T>& in_values,
                           std::vector<T>& out_values,
                           unsigned int receiving_process)
  {
    #ifdef HAS_MPI
    const std::size_t comm_size = MPI::size(comm);

    // Get data size on each process
    std::vector<int> pcounts(comm_size);
    const int local_size = in_values.size();
    MPI_Gather(const_cast<int*>(&local_size), 1, mpi_type<int>(),
               pcounts.data(), 1, mpi_type<int>(),
               receiving_process, comm);

    // Build offsets
    std::vector<int> offsets(comm_size + 1, 0);
    for (std::size_t i = 1; i <= comm_size; ++i)
      offsets[i] = offsets[i - 1] + pcounts[i - 1];

    const std::size_t n = std::accumulate(pcounts.begin(), pcounts.end(), 0);
    out_values.resize(n);
    MPI_Gatherv(const_cast<T*>(in_values.data()), in_values.size(),
                mpi_type<T>(),
                out_values.data(), pcounts.data(), offsets.data(),
                mpi_type<T>(), receiving_process, comm);
    #else
    out_values = in_values;
    #endif
  }
  //---------------------------------------------------------------------------
  inline void dolfin::MPI::gather(MPI_Comm comm,
                                  const std::string& in_values,
                                  std::vector<std::string>& out_values,
                                  unsigned int receiving_process)
  {
    #ifdef HAS_MPI
    const std::size_t comm_size = MPI::size(comm);

    // Get data size on each process
    std::vector<int> pcounts(comm_size);
    int local_size = in_values.size();
    MPI_Gather(&local_size, 1, MPI_INT,
               pcounts.data(), 1,MPI_INT,
               receiving_process, comm);

    // Build offsets
    std::vector<int> offsets(comm_size + 1, 0);
    for (std::size_t i = 1; i <= comm_size; ++i)
      offsets[i] = offsets[i - 1] + pcounts[i - 1];

    // Gather
    const std::size_t n = std::accumulate(pcounts.begin(), pcounts.end(), 0);
    std::vector<char> _out(n);
    MPI_Gatherv(const_cast<char*>(in_values.data()), in_values.size(),
                MPI_CHAR,
                _out.data(), pcounts.data(), offsets.data(),
                MPI_CHAR, receiving_process, comm);

    // Rebuild
    out_values.resize(comm_size);
    for (std::size_t p = 0; p < comm_size; ++p)
    {
      out_values[p] = std::string(_out.begin() + offsets[p],
                                  _out.begin() + offsets[p + 1]);
    }
    #else
    out_values.clear();
    out_values.push_back(in_values);
    #endif
  }
  //---------------------------------------------------------------------------
  inline void dolfin::MPI::all_gather(MPI_Comm comm,
                                      const std::string& in_values,
                                      std::vector<std::string>& out_values)
  {
    #ifdef HAS_MPI
    const std::size_t comm_size = MPI::size(comm);

    // Get data size on each process
    std::vector<int> pcounts(comm_size);
    int local_size = in_values.size();
    MPI_Allgather(&local_size, 1, MPI_INT, pcounts.data(), 1, MPI_INT, comm);

    // Build offsets
    std::vector<int> offsets(comm_size + 1, 0);
    for (std::size_t i = 1; i <= comm_size; ++i)
      offsets[i] = offsets[i - 1] + pcounts[i - 1];

    // Gather
    const std::size_t n = std::accumulate(pcounts.begin(), pcounts.end(), 0);
    std::vector<char> _out(n);
    MPI_Allgatherv(const_cast<char*>(in_values.data()), in_values.size(),
                   MPI_CHAR, _out.data(), pcounts.data(), offsets.data(),
                   MPI_CHAR, comm);

    // Rebuild
    out_values.resize(comm_size);
    for (std::size_t p = 0; p < comm_size; ++p)
    {
      out_values[p] = std::string(_out.begin() + offsets[p],
                                  _out.begin() + offsets[p + 1]);
    }
    #else
    out_values.clear();
    out_values.push_back(in_values);
    #endif
  }
  //-------------------------------------------------------------------------
  template<typename T>
    void dolfin::MPI::all_gather(MPI_Comm comm,
                                 const std::vector<T>& in_values,
                                 std::vector<T>& out_values)
  {
    #ifdef HAS_MPI
    out_values.resize(in_values.size()*MPI::size(comm));
    MPI_Allgather(const_cast<T*>(in_values.data()), in_values.size(),
                  mpi_type<T>(),
                  out_values.data(), in_values.size(), mpi_type<T>(),
                  comm);
#else
    out_values = in_values;
    #endif
  }
  //---------------------------------------------------------------------------
  template<typename T>
    void dolfin::MPI::all_gather(MPI_Comm comm,
                                 const std::vector<T>& in_values,
                                 std::vector<std::vector<T>>& out_values)
  {
    #ifdef HAS_MPI
    const std::size_t comm_size = MPI::size(comm);

    // Get data size on each process
    std::vector<int> pcounts;
    const int local_size = in_values.size();
    MPI::all_gather(comm, local_size, pcounts);
    dolfin_assert(pcounts.size() == comm_size);

    // Build offsets
    std::vector<int> offsets(comm_size + 1, 0);
    for (std::size_t i = 1; i <= comm_size; ++i)
      offsets[i] = offsets[i - 1] + pcounts[i - 1];

      // Gather data
    const std::size_t n = std::accumulate(pcounts.begin(), pcounts.end(), 0);
    std::vector<T> recvbuf(n);
    MPI_Allgatherv(const_cast<T*>(in_values.data()), in_values.size(),
                   mpi_type<T>(),
                   recvbuf.data(), pcounts.data(), offsets.data(),
                   mpi_type<T>(), comm);

    // Repack data
    out_values.resize(comm_size);
    for (std::size_t p = 0; p < comm_size; ++p)
    {
      out_values[p].resize(pcounts[p]);
      for (int i = 0; i < pcounts[p]; ++i)
        out_values[p][i] = recvbuf[offsets[p] + i];
    }
    #else
    out_values.clear();
    out_values.push_back(in_values);
    #endif
  }
  //---------------------------------------------------------------------------
  template<typename T>
    void dolfin::MPI::all_gather(MPI_Comm comm, const T in_value,
                                 std::vector<T>& out_values)
  {
    #ifdef HAS_MPI
    out_values.resize(MPI::size(comm));
    MPI_Allgather(const_cast<T*>(&in_value), 1, mpi_type<T>(),
                  out_values.data(), 1, mpi_type<T>(), comm);
    #else
    out_values.clear();
    out_values.push_back(in_value);
    #endif
  }
  //---------------------------------------------------------------------------
  template<typename T, typename X>
    T dolfin::MPI::all_reduce(MPI_Comm comm, const T& value, X op)
  {
    #ifdef HAS_MPI
    T out;
    MPI_Allreduce(const_cast<T*>(&value), &out, 1, mpi_type<T>(), op, comm);
    return out;
    #else
    return value;
    #endif
  }
  //---------------------------------------------------------------------------
  template<typename T> T dolfin::MPI::max(MPI_Comm comm, const T& value)
  {
    #ifdef HAS_MPI
    // Enforce cast to MPI_Op; this is needed because template dispatch may
    // not recognize this is possible, e.g. C-enum to unsigned int in SGI MPT
    MPI_Op op = static_cast<MPI_Op>(MPI_MAX);
    return all_reduce(comm, value, op);
    #else
    return value;
    #endif
  }
  //---------------------------------------------------------------------------
  template<typename T> T dolfin::MPI::min(MPI_Comm comm, const T& value)
  {
    #ifdef HAS_MPI
    // Enforce cast to MPI_Op; this is needed because template dispatch may
    // not recognize this is possible, e.g. C-enum to unsigned int in SGI MPT
    MPI_Op op = static_cast<MPI_Op>(MPI_MIN);
    return all_reduce(comm, value, op);
    #else
    return value;
    #endif
  }
  //---------------------------------------------------------------------------
  template<typename T> T dolfin::MPI::sum(MPI_Comm comm, const T& value)
  {
    #ifdef HAS_MPI
    // Enforce cast to MPI_Op; this is needed because template dispatch may
    // not recognize this is possible, e.g. C-enum to unsigned int in SGI MPT
    MPI_Op op = static_cast<MPI_Op>(MPI_SUM);
    return all_reduce(comm, value, op);
    #else
    return value;
    #endif
  }
  //---------------------------------------------------------------------------
  template<typename T> T dolfin::MPI::avg(MPI_Comm comm, const T& value)
  {
    #ifdef HAS_MPI
    dolfin_error("MPI.h",
                 "perform average reduction",
                 "Not implemented for this type");
    #else
    return value;
    #endif
  }
  //---------------------------------------------------------------------------
  template<typename T>
    void dolfin::MPI::send_recv(MPI_Comm comm,
                                const std::vector<T>& send_value,
                                unsigned int dest, int send_tag,
                                std::vector<T>& recv_value,
                                unsigned int source, int recv_tag)
  {
    #ifdef HAS_MPI
    std::size_t send_size = send_value.size();
    std::size_t recv_size = 0;
    MPI_Status mpi_status;
    MPI_Sendrecv(&send_size, 1, mpi_type<std::size_t>(), dest, send_tag,
                 &recv_size, 1, mpi_type<std::size_t>(), source, recv_tag,
                 comm, &mpi_status);

    recv_value.resize(recv_size);
    MPI_Sendrecv(const_cast<T*>(send_value.data()), send_value.size(),
                 mpi_type<T>(), dest, send_tag,
                 recv_value.data(), recv_size, mpi_type<T>(),
                 source, recv_tag, comm, &mpi_status);
    #else
      dolfin_error("MPI.h",
      "call MPI::send_recv",
      "DOLFIN has been configured without MPI support");
    #endif
  }
  //---------------------------------------------------------------------------
  template<typename T>
    void dolfin::MPI::send_recv(MPI_Comm comm,
                                const std::vector<T>& send_value,
                                unsigned int dest,
                                std::vector<T>& recv_value,
                                unsigned int source)
  {
    MPI::send_recv(comm, send_value, dest, 0, recv_value, source, 0);
  }
  //---------------------------------------------------------------------------
  // Specialization for dolfin::log::Table class
  // NOTE: This function is not trully "all_reduce", it reduces to rank 0
  //       and returns zero Table on other ranks.
  #ifdef HAS_MPI
  template<>
    Table dolfin::MPI::all_reduce(MPI_Comm, const Table&, MPI_Op);
  #endif
  //---------------------------------------------------------------------------
  // Specialization for dolfin::log::Table class
  #ifdef HAS_MPI
  template<>
    Table dolfin::MPI::avg(MPI_Comm, const Table&);
  #endif
  //---------------------------------------------------------------------------
}

#endif