Flow123d  JS_before_hm-978-gf0793cd
distribution.cc
Go to the documentation of this file.
1 /*!
2  *
3  * Copyright (C) 2015 Technical University of Liberec. All rights reserved.
4  *
5  * This program is free software; you can redistribute it and/or modify it under
6  * the terms of the GNU General Public License version 3 as published by the
7  * Free Software Foundation. (http://www.gnu.org/licenses/gpl-3.0.en.html)
8  *
9  * This program is distributed in the hope that it will be useful, but WITHOUT
10  * ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
11  * FOR A PARTICULAR PURPOSE. See the GNU General Public License for more details.
12  *
13  *
14  * @file distribution.cc
15  * @ingroup system
16  * @brief Objects and functions for mesh partitioning.
17  */
18 
19 #include "system/global_defs.h"
20 #include "system/system.hh"
21 #include "distribution.hh"
22 
23 /**
24  * create a Distribution from local sizes (dim = np )
25  * (collective context)
26  */
27 Distribution::Distribution(const unsigned int size, MPI_Comm comm)
28 :communicator(comm),
29 lsizes(NULL)
30 {
33 // TODO: zavest odchytavani vyjimek a pouzivat new a delete
34 
35  // communicate global sizes array
36  starts= new unsigned int [np()+1];
37  unsigned int lsize=size; // since size is const
39  // count starts
40  starts[0]=0;
41  for( unsigned int i=1 ; i<=np(); i++) starts[i]+=starts[i-1];
42 }
43 
44 /**
45  * create a Distribution from local sizes (dim = np )
46  * (local context)
47  */
48 Distribution::Distribution(const unsigned int * const sizes, MPI_Comm comm)
49 :communicator(comm),
50  lsizes(NULL)
51 {
54 // TODO: zavest odchytavani vyjimek a pouzivat new a delete
55  starts= new unsigned int [np()+1];
56  starts[0]=0;
57  for(unsigned int i=0 ; i<np(); i++) starts[i+1]=starts[i]+sizes[i];
58 }
59 
60 /**
61  * constructor from existing PETSC vector
62  * (collective context)
63  */
64 Distribution::Distribution(const Vec &petsc_vector)
65 :communicator(PETSC_COMM_WORLD),
66  lsizes(NULL)
67 {
70 
71  const PetscInt *petsc_starts;
72  chkerr(VecGetOwnershipRanges(petsc_vector,&petsc_starts));
73 
74  starts= new unsigned int [np()+1];
75  for(unsigned int i=0 ; i<=np(); i++) starts[i]=petsc_starts[i];
76 }
77 
78 /**
79  * construct from given global size
80  * (collective context)
81  */
82 Distribution::Distribution(const DistributionType &type, unsigned int global_size, MPI_Comm comm)
83 :communicator(comm),
84  lsizes(NULL)
85 {
88  OLD_ASSERT( num_of_procs > 0, "MPI size is not positive, possibly broken MPI communicator.\n");
89 
90  if (type.type_ == Block) {
91  unsigned int reminder, per_proc;
92 
93  reminder=global_size % np(); per_proc=global_size / np();
94  // set perproc rows to each proc, but for first "reminder" procs set one row more
95  starts= new unsigned int [np()+1];
96  starts[0]=0;
97  for(unsigned int i=0; i<np(); i++)
98  starts[i+1]=starts[i]+per_proc+(i<reminder?1:0);
99 
100  } else if (type.type_ == Localized) {
101 
102  starts= new unsigned int [np()+1];
103  starts[0]=0;
104  for(unsigned int i=1; i<=np(); i++) starts[i]=global_size;
105  }
106  else {
107  OLD_ASSERT( 0 , "Cyclic distribution is not yet implemented.\n");
108  }
109  }
110 /**
111  * copy constructor
112  *
113  */
115 : communicator(distr.communicator)
116 {
118  my_proc=distr.my_proc;
119  starts= new unsigned int [np()+1];
120  memcpy(starts,distr.starts,(np()+1) * sizeof(unsigned int));
121  lsizes=NULL;
122 }
123 
124 
125 /**
126  * find the proc to which belongs index "idx" in the distribution
127  * use simple linear search, better binary search could be implemented
128  * (local context)
129  */
130 unsigned int Distribution::get_proc(unsigned int idx) const
131 {
132  OLD_ASSERT( starts,"Distribution is not initialized.\n");
133  OLD_ASSERT(idx < size(), "Index %d greater then distribution size %d.\n", idx, size());
134 
135  for(unsigned int i=0; i<np(); i++) {
136  if (is_on_proc(idx,i)) return (i);
137  }
138  OLD_ASSERT( 0 , "Can not find owner of index %d. \n", idx);
139  return (-1);
140 }
141 
142 const unsigned int * Distribution::get_lsizes_array()
143 {
144  if ( lsizes == NULL ) {
145  lsizes= new unsigned int [np()];
146  for(unsigned int i=0;i<np();i++) lsizes[i]=lsize(i);
147  }
148 
149  return lsizes;
150 }
151 
152 
153 
154 const unsigned int * Distribution::get_starts_array() const {
155  return starts;
156 }
157 
158 
159 
160 void Distribution::view(std::ostream &stream) const
161 {
162  stream << "[" <<myp() << "]" << "Distribution size: " << size() << " lsize: " << lsize() << " offset: " << begin() << " mpi_size: " << np() << endl;
163  for(unsigned int i=0; i<np();++i)
164  stream << "[" <<myp() << "]" << "proc: " << i << " offset: " << begin(i) << " lsize: " << lsize(i) << endl;
165 }
166 
167 /**
168  * Destructor.
169  */
171 {
172 // TODO: zavest odchytavani vyjimek a pouzivat new a delete
173  delete [] starts;
174  if (lsizes) delete [] lsizes;
175 }
176 
unsigned int size() const
get global size
unsigned int get_proc(unsigned int idx) const
get processor of the given index
int my_proc
my proc number
unsigned int * starts
starts[i] index of the first index on the proc i; starts[n_procs]=size of whole array ...
int MPI_Comm
Definition: mpi.h:141
unsigned int begin() const
Distribution(const unsigned int size, MPI_Comm comm)
Definition: distribution.cc:27
void chkerr(unsigned int ierr)
Replacement of new/delete operator in the spirit of xmalloc.
Definition: system.hh:148
int num_of_procs
number of procs
unsigned int lsize() const
const unsigned int * get_lsizes_array()
get local sizes array
#define OLD_ASSERT(...)
Definition: global_defs.h:131
void view(std::ostream &stream) const
distribution view
Global macros to enhance readability and debugging, general constants.
const MPI_Comm communicator
communicator
const unsigned int * get_starts_array() const
get local starts array
#define MPI_Allgather(sendbuf, sendcount, sendtype, recvbuf, recvcount, recvtype, comm)
Definition: mpi.h:580
unsigned int np() const
get num of processors
#define MPI_Comm_size
Definition: mpi.h:235
#define MPI_Comm_rank
Definition: mpi.h:236
unsigned int myp() const
get my processor
unsigned int * lsizes
local sizes
Support classes for parallel programing.
#define MPI_INT
Definition: mpi.h:160
bool is_on_proc(unsigned int idx, unsigned int proc) const