DataTransferKit - Multiphysics Solution Transfer Services  2.0
DTK_CenterDistributor_impl.hpp
Go to the documentation of this file.
1 //---------------------------------------------------------------------------//
2 /*
3  Copyright (c) 2012, Stuart R. Slattery
4  All rights reserved.
5 
6  Redistribution and use in source and binary forms, with or without
7  modification, are permitted provided that the following conditions are
8  met:
9 
10  *: Redistributions of source code must retain the above copyright
11  notice, this list of conditions and the following disclaimer.
12 
13  *: Redistributions in binary form must reproduce the above copyright
14  notice, this list of conditions and the following disclaimer in the
15  documentation and/or other materials provided with the distribution.
16 
17  *: Neither the name of the University of Wisconsin - Madison nor the
18  names of its contributors may be used to endorse or promote products
19  derived from this software without specific prior written permission.
20 
21  THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS
22  "AS IS" AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT
23  LIMITED TO, THE IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR
24  A PARTICULAR PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL THE COPYRIGHT
25  HOLDER OR CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL,
26  SPECIAL, EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT
27  LIMITED TO, PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE,
28  DATA, OR PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY
29  THEORY OF LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT
30  (INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE
31  OF THIS SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
32 */
33 //---------------------------------------------------------------------------//
39 //---------------------------------------------------------------------------//
40 
41 #ifndef DTK_CENTERDISTRIBUTOR_IMPL_HPP
42 #define DTK_CENTERDISTRIBUTOR_IMPL_HPP
43 
44 #include <algorithm>
45 #include <limits>
46 
47 #include <Teuchos_Array.hpp>
48 #include <Teuchos_CommHelpers.hpp>
49 
50 #include "DTK_DBC.hpp"
51 
52 namespace DataTransferKit
53 {
54 //---------------------------------------------------------------------------//
58 template <int DIM>
60  const Teuchos::RCP<const Teuchos::Comm<int>> &comm,
61  const Teuchos::ArrayView<const double> &source_centers,
62  const Teuchos::ArrayView<const double> &target_centers, const double radius,
63  Teuchos::Array<double> &target_decomp_source_centers )
64  : d_distributor( new Tpetra::Distributor( comm ) )
65 {
66  DTK_REQUIRE( 0 == source_centers.size() % DIM );
67  DTK_REQUIRE( 0 == target_centers.size() % DIM );
68 
69  // Build the import/export data.
70  Teuchos::Array<int> export_procs;
71  {
72  // Compute the radius to expand the local domain with.
73  double radius_tol = 1.0e-2;
74  double radius_expand = radius * ( 1.0 + radius_tol );
75 
76  // Gather the bounding domains for each target proc.
77  CloudDomain<DIM> local_target_domain =
78  localCloudDomain( target_centers );
79  local_target_domain.expand( radius_expand );
80  Teuchos::Array<CloudDomain<DIM>> global_target_domains(
81  comm->getSize() );
82  Teuchos::gatherAll<int, CloudDomain<DIM>>(
83  *comm, 1, &local_target_domain, global_target_domains.size(),
84  global_target_domains.getRawPtr() );
85 
86  // Get those that are neighbors to this source proc.
87  CloudDomain<DIM> local_source_domain =
88  localCloudDomain( source_centers );
89  Teuchos::Array<CloudDomain<DIM>> neighbor_target_domains;
90  Teuchos::Array<int> neighbor_ranks;
91  for ( unsigned i = 0; i < global_target_domains.size(); ++i )
92  {
93  if ( local_source_domain.checkForIntersection(
94  global_target_domains[i] ) )
95  {
96  neighbor_target_domains.push_back( global_target_domains[i] );
97  neighbor_ranks.push_back( i );
98  }
99  }
100  global_target_domains.clear();
101 
102  // Find the procs to which the sources will be sent.
103  Teuchos::ArrayView<const double> source_point;
104  for ( unsigned source_id = 0; source_id < source_centers.size() / DIM;
105  ++source_id )
106  {
107  source_point = source_centers.view( DIM * source_id, DIM );
108  for ( unsigned b = 0; b < neighbor_target_domains.size(); ++b )
109  {
110  if ( neighbor_target_domains[b].pointInDomain( source_point ) )
111  {
112  export_procs.push_back( neighbor_ranks[b] );
113  d_export_ids.push_back( source_id );
114  }
115  }
116  }
117  }
118  DTK_CHECK( d_export_ids.size() == export_procs.size() );
119  d_num_exports = d_export_ids.size();
120 
121  // Create the communication plan.
122  Teuchos::ArrayView<int> export_procs_view = export_procs();
123  d_num_imports = d_distributor->createFromSends( export_procs_view );
124  export_procs.clear();
125 
126  // Unroll the coordinates to handle cases where single source centers
127  // may have multiple destinations.
128  Teuchos::Array<unsigned>::const_iterator export_id_it;
129  Teuchos::Array<double> src_coords( d_num_exports * DIM );
130  for ( int n = 0; n < d_num_exports; ++n )
131  {
132  src_coords( DIM * n, DIM )
133  .assign( source_centers( DIM * d_export_ids[n], DIM ) );
134  }
135 
136  // Move the source center coordinates to the target decomposition.
137  Teuchos::ArrayView<const double> src_coords_view = src_coords();
138  target_decomp_source_centers.resize( d_num_imports * DIM );
139  d_distributor->doPostsAndWaits( src_coords_view, DIM,
140  target_decomp_source_centers() );
141 }
142 
143 //---------------------------------------------------------------------------//
148 template <int DIM>
149 template <class T>
151  const Teuchos::ArrayView<const T> &source_decomp_data,
152  const Teuchos::ArrayView<T> &target_decomp_data ) const
153 {
154  DTK_REQUIRE( d_num_imports == target_decomp_data.size() );
155 
156  // Unroll the source data to handle cases where single data points may
157  // have multiple destinations.
158  Teuchos::Array<unsigned>::const_iterator export_id_it;
159  Teuchos::Array<T> src_data( d_num_exports );
160  typename Teuchos::Array<T>::iterator src_it;
161  for ( export_id_it = d_export_ids.begin(), src_it = src_data.begin();
162  export_id_it != d_export_ids.end(); ++export_id_it, ++src_it )
163  {
164  DTK_CHECK( *export_id_it < source_decomp_data.size() );
165  *src_it = source_decomp_data[*export_id_it];
166  }
167 
168  // Distribute.
169  Teuchos::ArrayView<const T> src_view = src_data();
170  d_distributor->doPostsAndWaits( src_view, 1, target_decomp_data );
171 }
172 
173 //---------------------------------------------------------------------------//
177 template <int DIM>
179  const Teuchos::ArrayView<const double> &centers ) const
180 {
181  Teuchos::Array<double> bounds( 2 * DIM, 0.0 );
182 
183  if ( centers.size() > 0 )
184  {
185  for ( int d = 0; d < DIM; ++d )
186  {
187  bounds[2 * d] = std::numeric_limits<double>::max();
188  bounds[2 * d + 1] = std::numeric_limits<double>::min();
189  }
190  }
191 
192  Teuchos::ArrayView<const double>::const_iterator center_it;
193  for ( center_it = centers.begin(); center_it != centers.end(); )
194  {
195  for ( int d = 0; d < DIM; ++d )
196  {
197  bounds[2 * d] = std::min( bounds[2 * d], *( center_it + d ) );
198  bounds[2 * d + 1] =
199  std::max( bounds[2 * d + 1], *( center_it + d ) );
200  }
201  center_it += DIM;
202  }
203 
204  return CloudDomain<DIM>( bounds.getRawPtr() );
205 }
206 
207 //---------------------------------------------------------------------------//
208 
209 } // end namespace DataTransferKit
210 
211 //---------------------------------------------------------------------------//
212 
213 #endif // end DTK_CENTERDISTRIBUTOR_IMPL_HPP
214 
215 //---------------------------------------------------------------------------//
216 // end DTK_CenterDistributor_impl.hpp
217 //---------------------------------------------------------------------------//
CenterDistributor(const Teuchos::RCP< const Teuchos::Comm< int >> &comm, const Teuchos::ArrayView< const double > &source_centers, const Teuchos::ArrayView< const double > &target_centers, const double radius, Teuchos::Array< double > &target_decomp_source_centers)
Constructor.
Assertions and Design-by-Contract for error handling.
Axis-aligned Cartesian cloud domain container.
void distribute(const Teuchos::ArrayView< const T > &source_decomp_data, const Teuchos::ArrayView< T > &target_decomp_data) const
Given a set of scalar values at the given source centers in the source decomposition, distribute them to the target decomposition.
DTK_BasicEntitySet.cpp.