openCARP
Doxygen code documentation for the open cardiac electrophysiology simulator openCARP
SF_abstract_vector.h
Go to the documentation of this file.
1 // ----------------------------------------------------------------------------
2 // openCARP is an open cardiac electrophysiology simulator.
3 //
4 // Copyright (C) 2020 openCARP project
5 //
6 // This program is licensed under the openCARP Academic Public License (APL)
7 // v1.0: You can use and redistribute it and/or modify it in non-commercial
8 // academic environments under the terms of APL as published by the openCARP
9 // project v1.0, or (at your option) any later version. Commercial use requires
10 // a commercial license (info@opencarp.org).
11 //
12 // This program is distributed without any warranty; see the openCARP APL for
13 // more details.
14 //
15 // You should have received a copy of the openCARP APL along with this program
16 // and can find it online: http://www.opencarp.org/license
17 // ----------------------------------------------------------------------------
18 
19 #ifndef _SF_ABSTRACT_VECTOR_H
20 #define _SF_ABSTRACT_VECTOR_H
21 
22 #include <mpi.h>
23 
24 #include <petscvec.h> // TODO: For scattering, VecScatter and friends
25 #include <petscis.h>
26 #include <petscsys.h> // TODO: For PETSC_COMM_WORLD
27 
28 #include "SF_container.h"
29 #include "SF_globals.h"
30 #include "SF_parallel_layout.h"
31 #include "SF_parallel_utils.h"
32 
33 #include "hashmap.hpp"
34 
35 namespace SF {
36 
37 // Forward declare the scattering class
38 class scattering;
39 
53 template<class T, class S>
55 
56  public:
57 
60 
62  int dpn = 0;
64 
66  virtual ~abstract_vector() = default;
67 
75  virtual void init(const meshdata<mesh_int_t, mesh_real_t>& imesh,
76  int idpn,
77  ltype inp_layout) = 0;
78 
84  virtual void init(const abstract_vector<T, S>& vec) = 0;
85 
94  virtual void init(int igsize, int ilsize, int idpn = 1, ltype ilayout = unset) = 0;
95 
105  inline std::tuple<int, int> init_common(const meshdata<mesh_int_t, mesh_real_t>& imesh,
106  int idpn,
107  ltype inp_layout)
108  {
109  mesh = &imesh;
110  dpn = idpn;
111  int N = 0, n = 0;
112 
113  switch(inp_layout) {
114  case algebraic:
115  N = mesh->pl.num_global_idx()*dpn;
116  n = mesh->pl.num_algebraic_idx()*dpn;
117  layout = algebraic;
118  break;
119 
120  case nodal:
121  N = mesh->l_numpts*dpn;
122  MPI_Allreduce(MPI_IN_PLACE, &N, 1, MPI_INT, MPI_SUM, mesh->comm);
123  n = mesh->l_numpts*dpn;
124  layout = nodal;
125  break;
126 
127  case elemwise:
128  N = mesh->g_numelem;
129  n = mesh->l_numelem;
130  layout = elemwise;
131  break;
132 
133  case ipwise:
134  N = mesh->g_numelem * dpn; //* EP_MAX_IPOINTS;
135  n = mesh->l_numelem * dpn; //* EP_MAX_IPOINTS;
136  layout = ipwise;
137  break;
138 
139  default: break;
140  }
141 
142  return std::tuple<int, int>(N, n);
143  }
144 
152  virtual void set(const vector<T>& idx, const vector<S>& vals, const bool additive = false) = 0;
153 
160  virtual void set(const vector<T>& idx, const S val) = 0;
161 
167  virtual void set(const S val) = 0;
168 
177  virtual void set(const T idx, const S val) = 0;
178 
186  virtual void get(const vector<T> & idx, S *out) = 0;
187 
195  virtual S get(const T idx) = 0;
196 
202  virtual void operator*=(const S sca) = 0;
203 
209  virtual void operator /= (const S sca) = 0;
210 
216  virtual void operator*=(const abstract_vector<T, S>& vec) = 0;
217 
224  virtual void add_scaled(const abstract_vector<T, S>& vec, S k) = 0;
225 
231  virtual void operator += (const abstract_vector<T,S> & vec) = 0;
232 
238  virtual void operator-=(const abstract_vector<T, S>& vec) = 0;
239 
245  virtual void operator+=(S c) = 0;
246 
252  virtual void operator=(const vector<S>& rhs) = 0;
253 
259  virtual void operator=(const abstract_vector<T, S>& rhs) = 0;
260 
268  virtual void shallow_copy(const abstract_vector<T, S>& v) = 0;
269 
277  virtual void deep_copy(const abstract_vector<T, S>& v) = 0;
278 
291  virtual void overshadow(const abstract_vector<T, S>& sub, bool member, int offset, int sz, bool share) = 0;
292 
298  virtual T lsize() const = 0;
299 
305  virtual T gsize() const = 0;
306 
313  virtual void get_ownership_range(T& start, T& stop) const = 0;
314 
322  virtual S* ptr() = 0;
323 
332  virtual const S* const_ptr() const = 0;
333 
341  virtual void release_ptr(S*& p) = 0;
342 
350  virtual void const_release_ptr(const S*& p) const = 0;
351 
357  virtual S mag() const = 0;
358 
364  virtual S sum() const = 0;
365 
371  virtual S min() const = 0;
372 
380  virtual S dot(const abstract_vector<T, S>& v) const = 0;
381 
387  virtual bool is_init() const = 0;
388 
394  virtual std::string to_string() const = 0;
395 
401  virtual bool equals(const abstract_vector<T,S> & rhs) const = 0;
402 
406  virtual void finish_assembly() = 0;
407 
415  virtual void forward(abstract_vector<T, S>& out, scattering &sc, bool add = false) = 0;
416 
424  virtual void backward(abstract_vector<T,S> & out, scattering &sc, bool add = false) = 0;
425 
432  virtual void apply_scattering(scattering& sc, bool fwd) = 0;
433 
442  inline size_t write_ascii(const char* file, bool write_header)
443  {
444  int size, rank;
445  MPI_Comm comm = this->mesh != NULL ? this->mesh->comm : PETSC_COMM_WORLD;
446  MPI_Comm_rank(comm, &rank);
447  MPI_Comm_size(comm, &size);
448 
449  int glb_size = this->gsize();
450  int loc_size = this->lsize();
451 
452  int err = 0;
453  FILE* fd = NULL;
454  long int nwr = 0;
455 
456  if(rank == 0) {
457  fd = fopen(file, "w");
458  if(!fd) err = 1;
459  }
460 
461  MPI_Allreduce(MPI_IN_PLACE, &err, 1, MPI_INT, MPI_MAX, comm);
462  if(err) {
463  treat_file_open_error(file, __func__, errno, false, rank);
464  return nwr;
465  }
466 
467  S* p = this->ptr();
468 
469  if(rank == 0) {
470  if(write_header)
471  fprintf(fd, "%d\n", glb_size);
472 
473  for(int i=0; i<loc_size/this->dpn; i++) {
474  for(int j=0; j<this->dpn; j++)
475  nwr += fprintf(fd, "%f ", p[i*this->dpn+j]);
476  nwr += fprintf(fd, "\n");
477  }
478 
479  vector<S> wbuff;
480  for(int pid=1; pid < size; pid++)
481  {
482  int rsize;
483  MPI_Status stat;
484 
485  MPI_Recv(&rsize, 1, MPI_INT, pid, SF_MPITAG, comm, &stat);
486  wbuff.resize(rsize);
487  MPI_Recv(wbuff.data(), rsize*sizeof(S), MPI_BYTE, pid, SF_MPITAG, comm, &stat);
488 
489  for(int i=0; i<rsize/this->dpn; i++) {
490  for(int j=0; j<this->dpn; j++)
491  nwr += fprintf(fd, "%f ", wbuff[i*this->dpn+j]);
492  nwr += fprintf(fd, "\n");
493  }
494  }
495  fclose(fd);
496  }
497  else {
498  MPI_Send(&loc_size, 1, MPI_INT, 0, SF_MPITAG, comm);
499  MPI_Send(p, loc_size*sizeof(S), MPI_BYTE, 0, SF_MPITAG, comm);
500  }
501 
502  this->release_ptr(p);
503  MPI_Bcast(&nwr, 1, MPI_LONG, 0, comm);
504 
505  return nwr;
506  }
507 
515  template<typename V>
516  inline size_t write_binary(FILE* fd)
517  {
518  int size, rank;
519  MPI_Comm comm = this->mesh != NULL ? this->mesh->comm : PETSC_COMM_WORLD;
520  MPI_Comm_rank(comm, &rank); MPI_Comm_size(comm, &size);
521 
522  long int loc_size = this->lsize();
523  S* p = this->ptr();
524 
525  vector<V> buff(loc_size);
526  for(long int i=0; i<loc_size; i++) buff[i] = p[i];
527 
528  long int nwr = root_write(fd, buff, comm);
529 
530  this->release_ptr(p);
531  return nwr;
532  }
533 
535  template<typename V>
536  inline size_t write_binary(std::string file)
537  {
538  size_t nwr = 0;
539  MPI_Comm comm = this->mesh != NULL ? this->mesh->comm : PETSC_COMM_WORLD;
540  int rank; MPI_Comm_rank(comm, &rank);
541 
542  FILE* fd = NULL;
543  int error = 0;
544 
545  if(rank == 0) {
546  fd = fopen(file.c_str(), "w");
547  if(fd == NULL)
548  error++;
549  }
550 
551  MPI_Allreduce(MPI_IN_PLACE, &error, 1, MPI_INT, MPI_SUM, comm);
552 
553  if(error == 0) {
554  nwr = this->write_binary<V>(fd);
555  fclose(fd);
556  }
557  else {
558  treat_file_open_error(file.c_str(), __func__, errno, false, rank);
559  }
560 
561  return nwr;
562  }
563 
564  template<typename V>
565  inline size_t read_binary(FILE* fd)
566  {
567  MPI_Comm comm = this->mesh != NULL ? this->mesh->comm : PETSC_COMM_WORLD;
568 
569  size_t loc_size = this->lsize();
570  S* p = this->ptr();
571  vector<V> buff(loc_size);
572 
573  size_t nrd = root_read(fd, buff, comm);
574 
575  for(size_t i=0; i<loc_size; i++)
576  p[i] = buff[i];
577 
578  this->release_ptr(p);
579 
580  return nrd;
581  }
582 
583  template<typename V>
584  inline size_t read_binary(std::string file)
585  {
586  MPI_Comm comm = mesh != NULL ? mesh->comm : PETSC_COMM_WORLD;
587 
588  size_t nrd = 0;
589  int rank; MPI_Comm_rank(comm, &rank);
590 
591  FILE* fd = NULL;
592  int error = 0;
593 
594  if(rank == 0) {
595  fd = fopen(file.c_str(), "r");
596  if(fd == NULL)
597  error++;
598  }
599 
600  MPI_Allreduce(MPI_IN_PLACE, &error, 1, MPI_INT, MPI_SUM, comm);
601 
602  if(error == 0) {
603  nrd = read_binary<V>(fd);
604  fclose(fd);
605  }
606  else {
607  treat_file_open_error(file.c_str(), __func__, errno, false, rank);
608  }
609 
610  return nrd;
611  }
612 
613  inline size_t read_ascii(FILE* fd)
614  {
615  MPI_Comm comm = this->mesh != NULL ? this->mesh->comm : PETSC_COMM_WORLD;
616 
617  size_t loc_size = this->lsize();
618  S* p = this->ptr();
619 
620  size_t nrd = root_read_ascii(fd, p, loc_size, comm, false);
621  return nrd;
622  }
623 
624  inline size_t read_ascii(std::string file)
625  {
626  MPI_Comm comm = this->mesh != NULL ? this->mesh->comm : PETSC_COMM_WORLD;
627  int rank; MPI_Comm_rank(comm, &rank);
628 
629  size_t nrd = 0;
630 
631  FILE* fd = NULL;
632  int error = 0;
633 
634  if(rank == 0) {
635  fd = fopen(file.c_str(), "r");
636  if(fd == NULL)
637  error++;
638  }
639 
640  MPI_Allreduce(MPI_IN_PLACE, &error, 1, MPI_INT, MPI_SUM, comm);
641 
642  if(error == 0) {
643  nrd = read_ascii(fd);
644  if(fd) fclose(fd);
645  }
646  else {
647  treat_file_open_error(file.c_str(), __func__, errno, false, rank);
648  }
649 
650  return nrd;
651  }
652 };
653 
654 
657 template<class T, class S>
658 inline bool is_init(const abstract_vector<T,S>* v)
659 {
660  return (v && v->is_init());
661 }
662 
663 
664 // TODO: This is a type of Vector and can most likely be integrated into the
665 // `abstract_vector` type in some way. Everything here and also in
666 // `SF_parallel_layout.h` probably need to be abstracted away, or made less
667 // PETSc specific?
676 {
677  public:
678  Vec b_buff;
679  VecScatter vec_sc;
680 
682 
684  scattering() : b_buff(NULL), vec_sc(NULL)
685  {}
686 
689  {
690  if(b_buff) VecDestroy(&b_buff);
691  if(vec_sc) VecScatterDestroy(&vec_sc);
692  }
699  template<class T, class S>
700  inline void forward(abstract_vector<T,S> & in, abstract_vector<T,S> & out, bool add = false)
701  {
702  in.forward(out, *this, add);
703  }
704 
711  template<class T, class S>
712  inline void backward(abstract_vector<T,S> & in, abstract_vector<T,S> & out, bool add = false)
713  {
714  in.backward(out, *this, add);
715  }
716 
723  template<class T, class S>
724  inline void operator()(abstract_vector<T,S> & v, bool fwd)
725  {
726  v.apply_scattering(*this, fwd);
727  }
728 };
729 
738 template<class T>
740 {
741  private:
744 
745  public:
758  inline
760  const vector<T> & nbr_a,
761  const vector<T> & nbr_b,
762  const size_t gsize_a,
763  const size_t gsize_b,
764  const short dpn)
765  {
766  assert(_registry.count(spec) == 0);
767 
768  scattering* sc = new scattering();
769  _registry[spec] = sc;
770 
771  IS is_a, is_b;
772  int err = 0;
773  vector<SF_int> idx_a(nbr_a.size() * dpn);
774  vector<SF_int> idx_b(nbr_b.size() * dpn);
775 
776  // we copy indexing into new containers because of dpn and also
777  // because of the implicit typecast between T and SF_int
778  for(size_t i=0; i<nbr_a.size(); i++)
779  for(short j=0; j<dpn; j++) idx_a[i*dpn+j] = nbr_a[i]*dpn+j;
780 
781  for(size_t i=0; i<nbr_b.size(); i++)
782  for(short j=0; j<dpn; j++) idx_b[i*dpn+j] = nbr_b[i]*dpn+j;
783 
784  err += ISCreateGeneral(PETSC_COMM_WORLD, idx_a.size(), idx_a.data(), PETSC_COPY_VALUES, &is_a);
785  err += ISCreateGeneral(PETSC_COMM_WORLD, idx_b.size(), idx_b.data(), PETSC_COPY_VALUES, &is_b);
786 
787  err += ISSetPermutation(is_b);
788 
789  Vec a;
790  err += VecCreateMPI(PETSC_COMM_WORLD, idx_a.size(), gsize_a*dpn, &a);
791  err += VecCreateMPI(PETSC_COMM_WORLD, idx_b.size(), gsize_b*dpn, &sc->b_buff);
792 
793  err += VecScatterCreate(a, is_a, sc->b_buff, is_b, &sc->vec_sc);
794 
795  err += VecDestroy(&a);
796  err += ISDestroy (&is_a);
797  err += ISDestroy (&is_b);
798 
799  if(err)
800  fprintf(stderr, "%s : an error has occurred! Scattering spec: Mesh ID %d, permutation ID %d, numbering %d, dpn %d!\n",
801  __func__, int(spec.v1), int(spec.v2), int(spec.v3), int(spec.v4));
802 
803  sc->idx_a = idx_a;
804  sc->idx_b = idx_b;
805 
806  return sc;
807  }
808 
809 
815  inline
817  const vector<T> & layout_a,
818  const vector<T> & layout_b,
819  const vector<T> & idx_a,
820  const vector<T> & idx_b,
821  const int rank,
822  const int dpn)
823  {
824  // scattering spec must not exist yet
825  assert(_registry.count(spec) == 0);
826 
827  // the local index sets have to be of equal size.
828  assert(idx_a.size() == idx_b.size());
829 
830  scattering* sc = new scattering();
831  _registry[spec] = sc;
832 
833  IS is_a, is_b;
834  Vec a;
835  int err = 0;
836  vector<SF_int> sidx_a(idx_a.size() * dpn);
837  vector<SF_int> sidx_b(idx_b.size() * dpn);
838 
839  T lsize_a = layout_a[rank+1] - layout_a[rank];
840  T lsize_b = layout_b[rank+1] - layout_b[rank];
841  T gsize_a = layout_a[layout_a.size()-1];
842  T gsize_b = layout_b[layout_a.size()-1];
843 
844  // we copy indexing into new containers because of dpn and also
845  // because of the implicit typecast between T and SF_int
846  for(size_t i=0; i<idx_a.size(); i++)
847  for(short j=0; j<dpn; j++) sidx_a[i*dpn+j] = idx_a[i]*dpn+j;
848 
849  for(size_t i=0; i<idx_b.size(); i++)
850  for(short j=0; j<dpn; j++) sidx_b[i*dpn+j] = idx_b[i]*dpn+j;
851 
852  err += ISCreateGeneral(PETSC_COMM_WORLD, sidx_a.size(), sidx_a.data(), PETSC_COPY_VALUES, &is_a);
853  err += ISCreateGeneral(PETSC_COMM_WORLD, sidx_b.size(), sidx_b.data(), PETSC_COPY_VALUES, &is_b);
854 
855  err += VecCreateMPI(PETSC_COMM_WORLD, lsize_a*dpn, gsize_a*dpn, &a);
856  err += VecCreateMPI(PETSC_COMM_WORLD, lsize_b*dpn, gsize_b*dpn, &sc->b_buff);
857 
858  err += VecScatterCreate(a, is_a, sc->b_buff, is_b, &sc->vec_sc);
859 
860  err += VecDestroy(&a);
861  err += ISDestroy (&is_a);
862  err += ISDestroy (&is_b);
863 
864  if(err)
865  fprintf(stderr, "%s : an error has occurred! Scattering spec: To ID %d, From ID %d, numbering %d, dpn %d!\n",
866  __func__, int(spec.v1), int(spec.v2), int(spec.v3), int(spec.v4));
867 
868  sc->idx_a = idx_a;
869  sc->idx_b = idx_b;
870 
871  return sc;
872  }
873 
885  {
886  scattering* ret = NULL;
887 
888  // scattering must be present
889  if(_registry.count(spec))
890  ret = _registry[spec];
891 
892  return ret;
893  }
894 
898  inline void free_scatterings()
899  {
900  typename hashmap::unordered_map<quadruple<int>,scattering*>::iterator it;
901  for(it = _registry.begin(); it != _registry.end(); ++it) {
902  if(it->second) {
903  delete it->second;
904  it->second = NULL;
905  }
906  }
907  }
908 };
909 
910 
911 } // namespace SF
912 
913 
914 #endif // _SF_ABSTRACT_VECTOR_H
Basic containers.
#define SF_MPITAG
the MPI tag when communicating
Definition: SF_globals.h:29
Classes and algorithms related to the layout of distributed meshes.
Simple utility functions for parallel data.
size_t read_ascii(FILE *fd)
size_t read_binary(std::string file)
virtual void get(const vector< T > &idx, S *out)=0
virtual void apply_scattering(scattering &sc, bool fwd)=0
virtual S mag() const =0
virtual S * ptr()=0
virtual void operator-=(const abstract_vector< T, S > &vec)=0
virtual void release_ptr(S *&p)=0
virtual std::string to_string() const =0
virtual S min() const =0
virtual void operator*=(const abstract_vector< T, S > &vec)=0
virtual void init(int igsize, int ilsize, int idpn=1, ltype ilayout=unset)=0
virtual S sum() const =0
virtual T gsize() const =0
size_t write_ascii(const char *file, bool write_header)
std::tuple< int, int > init_common(const meshdata< mesh_int_t, mesh_real_t > &imesh, int idpn, ltype inp_layout)
size_t read_binary(FILE *fd)
virtual void operator/=(const S sca)=0
virtual void forward(abstract_vector< T, S > &out, scattering &sc, bool add=false)=0
virtual void set(const vector< T > &idx, const S val)=0
virtual bool is_init() const =0
ltype layout
used vector layout (nodal, algebraic, unset)
virtual void set(const vector< T > &idx, const vector< S > &vals, const bool additive=false)=0
virtual bool equals(const abstract_vector< T, S > &rhs) const =0
virtual void const_release_ptr(const S *&p) const =0
size_t read_ascii(std::string file)
virtual void operator=(const abstract_vector< T, S > &rhs)=0
virtual void deep_copy(const abstract_vector< T, S > &v)=0
size_t write_binary(std::string file)
write binary. Open file descriptor myself.
virtual void get_ownership_range(T &start, T &stop) const =0
virtual S get(const T idx)=0
virtual void shallow_copy(const abstract_vector< T, S > &v)=0
virtual void add_scaled(const abstract_vector< T, S > &vec, S k)=0
virtual void overshadow(const abstract_vector< T, S > &sub, bool member, int offset, int sz, bool share)=0
virtual void set(const S val)=0
virtual ~abstract_vector()=default
virtual void operator+=(const abstract_vector< T, S > &vec)=0
virtual void operator=(const vector< S > &rhs)=0
virtual void set(const T idx, const S val)=0
size_t write_binary(FILE *fd)
Write a vector to HD in binary. File descriptor is already set up.
virtual void backward(abstract_vector< T, S > &out, scattering &sc, bool add=false)=0
virtual T lsize() const =0
virtual void init(const abstract_vector< T, S > &vec)=0
int dpn
d.o.f. per mesh vertex.
virtual const S * const_ptr() const =0
virtual void operator*=(const S sca)=0
virtual void finish_assembly()=0
virtual void init(const meshdata< mesh_int_t, mesh_real_t > &imesh, int idpn, ltype inp_layout)=0
const meshdata< mesh_int_t, mesh_real_t > * mesh
the connected mesh
virtual void operator+=(S c)=0
virtual S dot(const abstract_vector< T, S > &v) const =0
overlapping_layout< T > pl
nodal parallel layout
Definition: SF_container.h:419
size_t l_numelem
local number of elements
Definition: SF_container.h:388
size_t l_numpts
local number of points
Definition: SF_container.h:390
size_t g_numelem
global number of elements
Definition: SF_container.h:387
MPI_Comm comm
the parallel mesh is defined on a MPI world
Definition: SF_container.h:393
The scatterer registry class.
scattering * register_scattering(const quadruple< int > spec, const vector< T > &layout_a, const vector< T > &layout_b, const vector< T > &idx_a, const vector< T > &idx_b, const int rank, const int dpn)
Register a scattering.
scattering * register_permutation(const quadruple< int > spec, const vector< T > &nbr_a, const vector< T > &nbr_b, const size_t gsize_a, const size_t gsize_b, const short dpn)
Register a permutation scattering.
scattering * get_scattering(const quadruple< int > spec)
Access an previously registered scattering.
void free_scatterings()
Free the registered scatterings.
Container for a PETSc VecScatter.
vector< SF_int > idx_a
~scattering()
Destructor.
vector< SF_int > idx_b
scattering()
Constructor.
Vec b_buff
A buffer vector which also defines the parallel layout of the "b" side.
void operator()(abstract_vector< T, S > &v, bool fwd)
Apply the scattering on a data vector.
void forward(abstract_vector< T, S > &in, abstract_vector< T, S > &out, bool add=false)
Forward scattering.
VecScatter vec_sc
The scatterer.
void backward(abstract_vector< T, S > &in, abstract_vector< T, S > &out, bool add=false)
Backward scattering.
A vector storing arbitrary data.
Definition: SF_vector.h:43
size_t size() const
The current size of the vector.
Definition: SF_vector.h:104
void resize(size_t n)
Resize a vector.
Definition: SF_vector.h:209
T * data()
Pointer to the vector's start.
Definition: SF_vector.h:91
hm_int count(const K &key) const
Check if key exists.
Definition: hashmap.hpp:579
Classes similar to unordered_set and unordered_map, but with better performance.
Definition: dense_mat.hpp:34
void treat_file_open_error(const char *file, const char *caller, const int errnum, const bool do_exit, int rank)
treat a file open error by displaying the errnum string interpretation and the caller
Definition: SF_io_base.h:83
bool is_init(const abstract_vector< T, S > *v)
size_t root_read(FILE *fd, vector< V > &vec, MPI_Comm comm)
Read binary data into a vector.
size_t root_read_ascii(FILE *fd, vector< V > &vec, MPI_Comm comm, bool int_data)
Read binary data into a vector.
size_t root_write(FILE *fd, const vector< V > &vec, MPI_Comm comm)
Write vector data binary to disk.