Rocstar  1.0
Rocstar multiphysics simulation application
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
RFC_Window_transfer.h
Go to the documentation of this file.
1 /* *******************************************************************
2  * Rocstar Simulation Suite *
3  * Copyright@2015, Illinois Rocstar LLC. All rights reserved. *
4  * *
5  * Illinois Rocstar LLC *
6  * Champaign, IL *
7  * www.illinoisrocstar.com *
8  * sales@illinoisrocstar.com *
9  * *
10  * License: See LICENSE file in top level of distribution package or *
11  * http://opensource.org/licenses/NCSA *
12  *********************************************************************/
13 /* *******************************************************************
14  * THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, *
15  * EXPRESS OR IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES *
16  * OF MERCHANTABILITY, FITNESS FOR A PARTICULAR PURPOSE AND *
17  * NONINFRINGEMENT. IN NO EVENT SHALL THE CONTRIBUTORS OR *
18  * COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER *
19  * LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, *
20  * ARISING FROM, OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE *
21  * USE OR OTHER DEALINGS WITH THE SOFTWARE. *
22  *********************************************************************/
23 // $Id: RFC_Window_transfer.h,v 1.14 2008/12/06 08:43:27 mtcampbe Exp $
24 
25 //===============================================================
26 // This file defines RFC_Window_transfer and RFC_Pane_transfer, the
27 // basic data structures for data transfer.
28 // Author: Xiangmin Jiao
29 // Date: June. 15, 2001
30 //===============================================================
31 
32 #ifndef RFC_WINDOW_TRANSFER_H
33 #define RFC_WINDOW_TRANSFER_H
34 
35 #include "RFC_Window_base.h"
36 #include "Vector_n.h"
37 #include "commpi.h"
38 
40 
42 template < class _Tag> class RFC_Data;
43 template < class _Tag> class RFC_Data_const;
44 struct Tag_facial {};
45 struct Tag_nodal {};
50 
51 // RFC_Pane_transfer is built based on RFC_Pane_base with extension of
52 // extra attribute for storing buffers for data transfer.
54 public:
57  friend class RFC_Window_transfer;
58 
59  RFC_Pane_transfer( COM::Pane *b, int c);
60  virtual ~RFC_Pane_transfer();
61 
63  const RFC_Window_transfer *window() const { return _window; }
64 
65  Real *pointer( int i) {
66  if ( i>=0) {
67  if ( is_master())
68  return Base::pointer(i);
69  else {
71  return &_data_buf.front();
72  }
73  }
74  else
75  return &_buffer[-i-1][0];
76  }
77  const Real *pointer( int i) const
78  { return const_cast<Self*>(this)->pointer(i); }
79  const Real *coordinates() const {
80  if ( is_master())
81  return Base::coordinates();
82  else {
83  return &_coor_buf.front();
84  }
85  }
86  bool need_recv( int i) const {
87  // Check whether an element is on the receiver list
88  if ( _to_recv == NULL) return true;
89  RFC_assertion( i>=1 && i<=size_of_faces());
90  return _to_recv[ i-1];
91  }
92 
93  Real *get_emm( int f) {
94  RFC_assertion( f>=1 && f<=size_of_faces());
95  return &_emm_buffer[ _emm_offset[f-1]];
96  }
97  const Real *get_emm( int f) const {
98  RFC_assertion( f>=1 && f<=size_of_faces());
99  return &_emm_buffer[ _emm_offset[f-1]];
100  }
101 
102  bool is_master() const { return _base->window()!=NULL; }
103 
104 private:
105  // Data member
106  RFC_Window_transfer *_window; // Point to its parent window.
107 
108  std::vector< std::vector<Real> > _buffer; // Buffer for PCG
109  std::vector< int> _emm_offset; // Element mass matrix
110  std::vector< Real> _emm_buffer;
111 
113  std::vector< Real> _coor_buf;
114  std::vector< Real> _data_buf;
115 
116  std::vector< MPI_Request> _msg_requests;
117  int *_to_recv;
118 
119  //======= Data members to reduce communication volume for data transfer
120  // List of faces and nodes to be received from owner process.
121  // These arrays are allocated only if this pane is cached copy.
122  std::vector<int> _recv_faces;
123  std::vector<int> _recv_nodes;
124 
125  // Maps from process ids to the lists of faces and nodes to be sent.
126  // These fields are used only if this pane is the master copy.
127  std::map<int, std::vector<int> > _send_faces;
128  std::map<int, std::vector<int> > _send_nodes;
129  std::map<int, std::vector<Real> > _send_buffers;
130 };
131 
132 // A window is a collection of panes.
133 class RFC_Window_transfer : public RFC_Window_derived< RFC_Pane_transfer> {
134 public:
137 
138  RFC_Window_transfer( COM::Window *b, int color,
139  MPI_Comm com, const char *pre=NULL,
140  const char *format=NULL);
141  virtual ~RFC_Window_transfer();
142 
143  RFC_Pane_transfer &pane( const int pid);
144  const RFC_Pane_transfer &pane( const int pid) const;
145 
146  // ======== Functions for supporting data transfer algorithms
147  void init_facial_buffers( const Facial_data &nd, int);
149  void delete_facial_buffers();
150 
151  void init_nodal_buffers( const Nodal_data &nd, int n, bool init_emm);
152  Nodal_data nodal_buffer( int);
153  void delete_nodal_buffers();
154 
155  // Set _to_recv tags for the next data transfer algorithm.
156  // If tag is NULL, reset the tags to NULL.
157  void set_tags( const COM::Attribute *tag);
158 
159  // ============ Communication subroutines for source panes ==================
161  bool replicated() const { return _replicated; }
162 
164  void replicate_metadata( int *pane_ids, int n);
166  void clear_replicated_data();
168  void incident_panes( std::vector<int> &pane_ids);
169 
171  void incident_faces( std::map< int, std::vector<int> > &) const;
172 
175  void replicate_metadata( const RFC_Window_transfer &opp_win);
176 
179  void replicate_data( const Facial_data_const &data, bool replicate_coor);
180  void replicate_data( const Nodal_data_const &data, bool replicate_coor);
181 
182  //============= communication subroutines for target panes ==================
183  void reduce_to_all( Nodal_data &, MPI_Op);
185 
186  // ===== Lower level communication routines =============
188  void barrier() const;
190  bool is_root() const;
191  int comm_rank() const
192  { return COMMPI_Initialized()?COMMPI_Comm_rank( _comm):0; }
193  int comm_size() const
194  { return COMMPI_Initialized()?COMMPI_Comm_size( _comm):1; }
195 
196  void wait_all( int n, MPI_Request *requests);
197  void wait_any( int n, MPI_Request *requests, int *index, MPI_Status *stat=NULL);
198 
199  void allreduce( Array_n &arr, MPI_Op op) const
200  { allreduce( arr.begin(), arr.end()-arr.begin(), op); }
201  void allreduce( Real *x, MPI_Op op) const
202  { allreduce( x, 1, op); }
203 
204 private:
205  void allreduce( Real *, int n, MPI_Op op) const;
206 
207  // Initialize the data structures _pane_map and _num_panes.
208  void counts_to_displs( const std::vector<int> &counts,
209  std::vector<int> &displs) const {
210  RFC_assertion( !displs.empty() && counts.size()>=displs.size()-1);
211  displs[0]=0;
212  for ( int i=1, size=displs.size(); i<size; ++i) {
213  displs[i] = displs[i-1]+counts[i-1];
214  }
215  }
216 
217  void init_send_buffer( int pane_id, int to_rank);
218  void init_recv_buffer( int pane_id, int from_rank);
219 
220 private:
221  int _buf_dim;
222  MPI_Comm _comm;
223  std::map< int, std::pair<int, int> > _pane_map;
224  std::vector<int> _num_panes;
225  std::map< int, RFC_Pane_transfer*> _replic_panes;
227 
228  std::set< std::pair<int, RFC_Pane_transfer*> > _panes_to_send; //<to_rank, p>
229  const std::string _prefix;
230  const int _IO_format;
231 };
232 
233 //================================================================
234 // Wrapper classes for accessing data fields of a window.
235 //================================================================
236 template <class _Tag>
237 class RFC_Data_const {
238 public:
239  typedef Vector_n Value;
242  typedef unsigned int Size;
243  typedef _Tag Tag;
244 
245  RFC_Data_const() : _id(0), _dim(0) {}
246  RFC_Data_const( const RFC_Window_transfer &win, const char *aname) {
247  const COM::Attribute *a = win.attribute( aname);
248  _id = a->id(); _dim = a->size_of_components();
249  }
250  RFC_Data_const( const COM::Attribute *a)
251  : _id(a->id()), _dim(a->size_of_components()) {}
252  RFC_Data_const( int i, int d)
253  : _id( i), _dim(d) {}
254 
255  int id() const { return _id; }
256  Size dimension() const { return _dim; }
257  bool is_valid( const RFC_Pane_transfer *p, int i) const
258  { return true; }
259 
260  // Get the data associated with node v
261  Value_const get_value( const RFC_Pane_transfer *p, int v) const {
262  RFC_assertion( v>=1);
263  return Value_const( p->pointer( _id) + (v-1)*_dim, _dim);
264  }
265 
266  Value_const get_value( const Real *buf, int v) const {
267  RFC_assertion( v>=1);
268  return Value_const( buf + (v-1)*_dim, _dim);
269  }
270 
271  Tag tag() const { return Tag(); }
272 protected:
273  int _id;
274  unsigned int _dim;
275 };
276 
277 
278 // Wrapper for node centered data.
279 template < class _Tag>
280 class RFC_Data : public RFC_Data_const<_Tag> {
281 public:
283  typedef Vector_n Value;
286  typedef unsigned int Size;
287  typedef _Tag Tag;
288 
289  RFC_Data() {}
290  RFC_Data( RFC_Window_transfer &win, const char *aname) : Base( win, aname) {}
291  RFC_Data( COM::Attribute *a) : Base( a) {}
292  RFC_Data( int i, int d) : Base(i, d) {}
293 
294  // Get the data associated with node v
296  RFC_assertion( v>=1);
297  return Value_nonconst( p->pointer( _id) + (v-1)*_dim, _dim);
298  }
299  // Get the data associated with node v
300  Value_const get_value( const RFC_Pane_transfer *p, int v) const {
301  RFC_assertion( v>=1);
302  return Value_const( p->pointer( _id) + (v-1)*_dim, _dim);
303  }
304  void set_value( RFC_Pane_transfer *p, int v, const Value_const &vec) {
305  RFC_assertion( _dim == vec.dimension() && v>=1);
306  std::copy( vec.begin(), vec.end(), p->pointer( _id) + (v-1)*_dim);
307  }
308 
309  Value_nonconst get_value( Real *buf, int v) const {
310  RFC_assertion( v>=1);
311  return Value_nonconst( buf + (v-1)*_dim, _dim);
312  }
313  Value_const get_value( const Real *buf, int v) const {
314  RFC_assertion( v>=1);
315  return Value_const( buf + (v-1)*_dim, _dim);
316  }
317  void set_value( Real *buf, int v, const Value_const &vec) {
318  RFC_assertion( _dim == vec.dimension() && v>=1);
319  std::copy( vec.begin(), vec.end(), buf + (v-1)*_dim);
320  }
321 
322  Tag tag() const { return Tag(); }
323 protected:
326 };
327 
329 
330 #endif
331 
332 
333 
334 
335 
336 
int COMMPI_Comm_rank(MPI_Comm c)
Definition: commpi.h:162
RFC_Data(int i, int d)
void replicate_data(const Facial_data_const &data, bool replicate_coor)
Replicate the given data from remote processes onto local process.
RFC_Data_const< Tag_facial > Facial_data_const
RFC_Pane_transfer & pane(const int pid)
std::vector< int > _num_panes
std::map< int, RFC_Pane_transfer * > _replic_panes
RFC_Window_transfer(COM::Window *b, int color, MPI_Comm com, const char *pre=NULL, const char *format=NULL)
const NT & d
void wait_any(int n, MPI_Request *requests, int *index, MPI_Status *stat=NULL)
Size dimension() const
const Real * coordinates() const
RFC_Data_const(const RFC_Window_transfer &win, const char *aname)
std::vector< int > _recv_faces
void counts_to_displs(const std::vector< int > &counts, std::vector< int > &displs) const
bool is_valid(const RFC_Pane_transfer *p, int i) const
const RFC_Window_transfer * window() const
Size dimension() const
Definition: Vector_n.h:79
void set_tags(const COM::Attribute *tag)
const Attribute * attribute(const char *f) const
Retrieve an attribute object from the base using the attribute name.
void wait_all(int n, MPI_Request *requests)
std::vector< Real > _coor_buf
Value_const get_value(const RFC_Pane_transfer *p, int v) const
Tag tag() const
double Real
Definition: mapbasic.h:322
Value_const get_value(const Real *buf, int v) const
Base * _base
Reference to its base object.
void barrier() const
Block until all processes of the window have reached here.
RFC_Pane_transfer(COM::Pane *b, int c)
std::map< int, std::vector< Real > > _send_buffers
int color() const
The color of the window for overlay or for data transfer (BLUE or GREEN).
Facial_data facial_buffer(int)
RFC_Pane_transfer Self
int size_of_faces() const
The total number of faces in the pane.
Value_const get_value(const Real *buf, int v) const
const std::string _prefix
RFC_Data(COM::Attribute *a)
RFC_Data_const< Tag_nodal > Nodal_data_const
*********************************************************************Illinois Open Source License ****University of Illinois NCSA **Open Source License University of Illinois All rights reserved ****Developed free of to any person **obtaining a copy of this software and associated documentation to deal with the Software without including without limitation the rights to and or **sell copies of the and to permit persons to whom the **Software is furnished to do subject to the following this list of conditions and the following disclaimers ****Redistributions in binary form must reproduce the above **copyright this list of conditions and the following **disclaimers in the documentation and or other materials **provided with the distribution ****Neither the names of the Center for Simulation of Advanced the University of nor the names of its **contributors may be used to endorse or promote products derived **from this Software without specific prior written permission ****THE SOFTWARE IS PROVIDED AS WITHOUT WARRANTY OF ANY **EXPRESS OR INCLUDING BUT NOT LIMITED TO THE WARRANTIES **OF FITNESS FOR A PARTICULAR PURPOSE AND **NONINFRINGEMENT IN NO EVENT SHALL THE CONTRIBUTORS OR **COPYRIGHT HOLDERS BE LIABLE FOR ANY DAMAGES OR OTHER WHETHER IN AN ACTION OF TORT OR **ARISING OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE **USE OR OTHER DEALINGS WITH THE SOFTWARE v
Definition: roccomf90.h:20
Real * pointer(int i)
Get the address of a given attribute with id i.
RFC_Data_const(const COM::Attribute *a)
int COMMPI_Comm_size(MPI_Comm c)
Definition: commpi.h:165
The base implementation of RFC_Pane.
RFC_Data_const< _Tag > Base
#define RFC_END_NAME_SPACE
Definition: rfc_basic.h:29
Const_pointer begin() const
Definition: Vector_n.h:81
*********************************************************************Illinois Open Source License ****University of Illinois NCSA **Open Source License University of Illinois All rights reserved ****Developed free of to any person **obtaining a copy of this software and associated documentation to deal with the Software without including without limitation the rights to ** copy
Definition: roccomf90.h:20
void clear_replicated_data()
Clear all the replicate data but keep metadata.
void reduce_maxabs_to_all(Nodal_data &)
Const_pointer end() const
Definition: Vector_n.h:82
Array_n Value_nonconst
void incident_faces(std::map< int, std::vector< int > > &) const
Obtain the list of incident faces in each pane of the opposite mesh.
RFC_Data< Tag_facial > Facial_data
std::map< int, std::pair< int, int > > _pane_map
bool replicated() const
Returns whether replication has been performed.
RFC_Window_transfer * window()
RFC_Window_transfer Self
RFC_Data(RFC_Window_transfer &win, const char *aname)
blockLoc i
Definition: read.cpp:79
RFC_Data< Tag_nodal > Nodal_data
#define RFC_BEGIN_NAME_SPACE
Definition: rfc_basic.h:28
RFC_Data_const(int i, int d)
Const_pointer end() const
Definition: Vector_n.h:122
void int int REAL * x
Definition: read.cpp:74
void set_value(RFC_Pane_transfer *p, int v, const Value_const &vec)
Contains declarations of MPI subroutines used in Roccom.
Value_nonconst get_value(RFC_Pane_transfer *p, int v)
const Real * get_emm(int f) const
const NT & n
Nodal_data nodal_buffer(int)
void incident_panes(std::vector< int > &pane_ids)
Determine the remote panes that are incident with local panes.
void init_nodal_buffers(const Nodal_data &nd, int n, bool init_emm)
void allreduce(Real *x, MPI_Op op) const
Value_nonconst get_value(Real *buf, int v) const
void set_value(Real *buf, int v, const Value_const &vec)
const Real * coordinates() const
Reusable implementation for derived class of RFC_Window_base.
void init_facial_buffers(const Facial_data &nd, int)
void replicate_metadata(int *pane_ids, int n)
Replicate the metadata of a remote pane only the local process.
unsigned int Size
std::vector< Real > _emm_buffer
std::vector< std::vector< Real > > _buffer
Const_pointer begin() const
Definition: Vector_n.h:121
Array_n_const Value_const
std::map< int, std::vector< int > > _send_nodes
const Real * pointer(int i) const
Value_const get_value(const RFC_Pane_transfer *p, int v) const
int COMMPI_Initialized()
Definition: commpi.h:168
void init_recv_buffer(int pane_id, int from_rank)
std::vector< int > _emm_offset
Array_n_const Value_const
std::set< std::pair< int, RFC_Pane_transfer * > > _panes_to_send
std::vector< MPI_Request > _msg_requests
#define RFC_assertion
Definition: rfc_basic.h:65
bool need_recv(int i) const
void reduce_to_all(Nodal_data &, MPI_Op)
RFC_Window_derived< RFC_Pane_transfer > Base
bool is_root() const
Check whether the process has rank 0.
std::map< int, std::vector< int > > _send_faces
RFC_Window_transfer * _window
std::vector< Real > _data_buf
void init_send_buffer(int pane_id, int to_rank)
void allreduce(Array_n &arr, MPI_Op op) const
std::vector< int > _recv_nodes