Flow123d  jenkins-Flow123d-windows-release-multijob-285
distribution.cc
Go to the documentation of this file.
1 /*!
2  *
3  * Copyright (C) 2007 Technical University of Liberec. All rights reserved.
4  *
5  * Please make a following refer to Flow123d on your project site if you use the program for any purpose,
6  * especially for academic research:
7  * Flow123d, Research Centre: Advanced Remedial Technologies, Technical University of Liberec, Czech Republic
8  *
9  * This program is free software; you can redistribute it and/or modify it under the terms
10  * of the GNU General Public License version 3 as published by the Free Software Foundation.
11  *
12  * This program is distributed in the hope that it will be useful, but WITHOUT ANY WARRANTY;
13  * without even the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.
14  * See the GNU General Public License for more details.
15  *
16  * You should have received a copy of the GNU General Public License along with this program; if not,
17  * write to the Free Software Foundation, Inc., 59 Temple Place - Suite 330, Boston, MA 021110-1307, USA.
18  *
19  *
20  * $Id$
21  * $Revision$
22  * $LastChangedBy$
23  * $LastChangedDate$
24  *
25  * @file
26  * @ingroup system
27  * @brief Objects and functions for mesh partitioning.
28  *
29  */
30 
31 #include "system/global_defs.h"
32 #include "system/system.hh"
33 #include "distribution.hh"
34 
35 /**
36  * create a Distribution from local sizes (dim = np )
37  * (collective context)
38  */
39 Distribution::Distribution(const unsigned int size, MPI_Comm comm)
40 :communicator(comm),
41 lsizes(NULL)
42 {
43  int ierr;
45  ASSERT( ! ierr , "Can not get MPI rank.\n" );
47  ASSERT( ! ierr , "Can not get MPI size.\n" );
48 // TODO: zavest odchytavani vyjimek a pouzivat new a delete
49 
50  // communicate global sizes array
51  starts=(unsigned int *)xmalloc((np()+1)*sizeof(unsigned int));
52  unsigned int lsize=size; // since size is const
54  // count starts
55  starts[0]=0;
56  for( unsigned int i=1 ; i<=np(); i++) starts[i]+=starts[i-1];
57 }
58 
59 /**
60  * create a Distribution from local sizes (dim = np )
61  * (local context)
62  */
63 Distribution::Distribution(const unsigned int * const sizes, MPI_Comm comm)
64 :communicator(comm),
65  lsizes(NULL)
66 {
67  int ierr;
69  ASSERT( ! ierr , "Can not get MPI rank.\n" );
71  ASSERT( ! ierr , "Can not get MPI size.\n" );
72 // TODO: zavest odchytavani vyjimek a pouzivat new a delete
73  starts=(unsigned int *)xmalloc((np()+1)*sizeof(unsigned int));
74  starts[0]=0;
75  for(unsigned int i=0 ; i<np(); i++) starts[i+1]=starts[i]+sizes[i];
76 }
77 
78 /**
79  * constructor from existing PETSC vector
80  * (collective context)
81  */
82 Distribution::Distribution(const Vec &petsc_vector)
83 :communicator(PETSC_COMM_WORLD),
84  lsizes(NULL)
85 {
86  int ierr;
88  ASSERT( ! ierr , "Can not get MPI rank.\n" );
90  ASSERT( ! ierr , "Can not get MPI size.\n" );
91 
92  const PetscInt *petsc_starts;
93  VecGetOwnershipRanges(petsc_vector,&petsc_starts);
94  ASSERT( ! ierr , "Can not get vector ownership range.\n" );
95 
96  starts=(unsigned int *)xmalloc((np()+1)*sizeof(int));
97  for(unsigned int i=0 ; i<=np(); i++) starts[i]=petsc_starts[i];
98 }
99 
100 /**
101  * construct from given global size
102  * (collective context)
103  */
104 Distribution::Distribution(const DistributionType &type, unsigned int global_size, MPI_Comm comm)
105 :communicator(comm),
106  lsizes(NULL)
107 {
108  int ierr;
110  ASSERT( ! ierr , "Can not get MPI rank.\n" );
112  ASSERT( ! ierr , "Can not get MPI size.\n" );
113  ASSERT( num_of_procs > 0, "MPI size is not positive, possibly broken MPI communicator.\n");
114 
115  if (type.type_ == Block) {
116  unsigned int reminder, per_proc;
117 
118  reminder=global_size % np(); per_proc=global_size / np();
119  // set perproc rows to each proc, but for first "reminder" procs set one row more
120  starts=(unsigned int *)xmalloc((np()+1)*sizeof(unsigned int));
121  starts[0]=0;
122  for(unsigned int i=0; i<np(); i++)
123  starts[i+1]=starts[i]+per_proc+(i<reminder?1:0);
124 
125  } else if (type.type_ == Localized) {
126 
127  starts=(unsigned int *)xmalloc((np()+1)*sizeof(unsigned int));
128  starts[0]=0;
129  for(unsigned int i=1; i<=np(); i++) starts[i]=global_size;
130  }
131  else {
132  ASSERT( 0 , "Cyclic distribution is not yet implemented.\n");
133  }
134  }
135 /**
136  * copy constructor
137  *
138  */
140 : communicator(distr.communicator)
141 {
143  my_proc=distr.my_proc;
144  starts=(unsigned int *)xmalloc((np()+1)*sizeof(unsigned int));
145  memcpy(starts,distr.starts,(np()+1) * sizeof(unsigned int));
146  lsizes=NULL;
147 }
148 
149 
150 /**
151  * find the proc to which belongs index "idx" in the distribution
152  * use simple linear search, better binary search could be implemented
153  * (local context)
154  */
155 unsigned int Distribution::get_proc(unsigned int idx) const
156 {
157  ASSERT( starts,"Distribution is not initialized.\n");
158  ASSERT(idx < size(), "Index %d greater then distribution size %d.\n", idx, size());
159 
160  for(unsigned int i=0; i<np(); i++) {
161  if (is_on_proc(idx,i)) return (i);
162  }
163  ASSERT( 0 , "Can not find owner of index %d. \n", idx);
164  return (-1);
165 }
166 
167 const unsigned int * Distribution::get_lsizes_array()
168 {
169  if ( lsizes == NULL ) {
170  lsizes=(unsigned int *) xmalloc(np()*sizeof(int));
171  for(unsigned int i=0;i<np();i++) lsizes[i]=lsize(i);
172  }
173 
174  return lsizes;
175 }
176 
177 
178 
179 const unsigned int * Distribution::get_starts_array() const {
180  return starts;
181 }
182 
183 
184 
185 void Distribution::view(std::ostream &stream) const
186 {
187  stream << "[" <<myp() << "]" << "Distribution size: " << size() << " lsize: " << lsize() << " offset: " << begin() << " mpi_size: " << np() << endl;
188  for(unsigned int i=0; i<np();++i)
189  stream << "[" <<myp() << "]" << "proc: " << i << " offset: " << begin(i) << " lsize: " << lsize(i) << endl;
190 }
191 
192 /**
193  * Destructor.
194  */
196 {
197 // TODO: zavest odchytavani vyjimek a pouzivat new a delete
198  xfree(starts);
199  if (lsizes) xfree(lsizes);
200 }
201 
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:39
int num_of_procs
number of procs
unsigned int lsize() const
const unsigned int * get_lsizes_array()
get local sizes array
void view(std::ostream &stream) const
distribution view
Global macros to enhance readability and debugging, general constants.
#define ASSERT(...)
Definition: global_defs.h:121
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
void * xmalloc(size_t size)
Memory allocation with checking.
Definition: system.cc:209
#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
#define xfree(p)
Definition: system.hh:108
bool is_on_proc(unsigned int idx, unsigned int proc) const