DataTransferKit - Multiphysics Solution Transfer Services  2.0
DTK_MovingLeastSquareReconstructionOperator_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_MOVINGLEASTSQUARERECONSTRUCTIONOPERATOR_IMPL_HPP
42 #define DTK_MOVINGLEASTSQUARERECONSTRUCTIONOPERATOR_IMPL_HPP
43 
44 #include "DTK_BasicEntityPredicates.hpp"
46 #include "DTK_DBC.hpp"
47 #include "DTK_LocalMLSProblem.hpp"
49 #include "DTK_PredicateComposition.hpp"
51 
52 #include <Teuchos_ArrayRCP.hpp>
53 #include <Teuchos_CommHelpers.hpp>
54 #include <Teuchos_ParameterList.hpp>
55 #include <Teuchos_Ptr.hpp>
56 
57 #include <Tpetra_MultiVector.hpp>
58 
59 namespace DataTransferKit
60 {
61 //---------------------------------------------------------------------------//
62 // Constructor.
63 template <class Basis, int DIM>
66  const Teuchos::RCP<const TpetraMap> &domain_map,
67  const Teuchos::RCP<const TpetraMap> &range_map,
68  const Teuchos::ParameterList &parameters )
69  : Base( domain_map, range_map )
70  , d_use_knn( false )
71  , d_knn( 0 )
72  , d_radius( 0.0 )
73  , d_domain_entity_dim( 0 )
74  , d_range_entity_dim( 0 )
75 {
76  // Determine if we are doing kNN search or radius search.
77  if ( parameters.isParameter( "Type of Search" ) )
78  {
79  if ( "Radius" == parameters.get<std::string>( "Type of Search" ) )
80  {
81  d_use_knn = false;
82  }
83  else if ( "Nearest Neighbor" ==
84  parameters.get<std::string>( "Type of Search" ) )
85  {
86  d_use_knn = true;
87  }
88  else
89  {
90  // Otherwise we got an invalid search type.
91  DTK_INSIST( false );
92  }
93  }
94 
95  // If we are doing kNN support get the number of neighbors.
96  if ( d_use_knn )
97  {
98  DTK_REQUIRE( parameters.isParameter( "Num Neighbors" ) );
99  d_knn = parameters.get<int>( "Num Neighbors" );
100  }
101 
102  // Otherwise we are doing the radius search so get the basis radius.
103  else
104  {
105  DTK_REQUIRE( parameters.isParameter( "RBF Radius" ) );
106  d_radius = parameters.get<double>( "RBF Radius" );
107  }
108 
109  // Get the topological dimension of the domain and range entities. This
110  // map will use their centroids for the point cloud.
111  if ( parameters.isParameter( "Domain Entity Dimension" ) )
112  {
113  d_domain_entity_dim = parameters.get<int>( "Domain Entity Dimension" );
114  }
115  if ( parameters.isParameter( "Range Entity Dimension" ) )
116  {
117  d_range_entity_dim = parameters.get<int>( "Range Entity Dimension" );
118  }
119 }
120 
121 //---------------------------------------------------------------------------//
122 // Setup the map operator.
123 template <class Basis, int DIM>
125  const Teuchos::RCP<FunctionSpace> &domain_space,
126  const Teuchos::RCP<FunctionSpace> &range_space )
127 {
128  DTK_REQUIRE( Teuchos::nonnull( domain_space ) );
129  DTK_REQUIRE( Teuchos::nonnull( range_space ) );
130 
131  // Extract the Support maps.
132  const Teuchos::RCP<const typename Base::TpetraMap> domain_map =
133  this->getDomainMap();
134  const Teuchos::RCP<const typename Base::TpetraMap> range_map =
135  this->getRangeMap();
136 
137  // Get the parallel communicator.
138  Teuchos::RCP<const Teuchos::Comm<int>> comm = domain_map->getComm();
139 
140  // Determine if we have range and domain data on this process.
141  bool nonnull_domain = Teuchos::nonnull( domain_space->entitySet() );
142  bool nonnull_range = Teuchos::nonnull( range_space->entitySet() );
143 
144  // Extract the source nodes and their ids.
145  Teuchos::ArrayRCP<double> source_centers;
146  Teuchos::ArrayRCP<GO> source_support_ids;
147  getNodeCoordsAndIds( domain_space, d_domain_entity_dim, source_centers,
148  source_support_ids );
149 
150  // Extract the target nodes and their ids.
151  Teuchos::ArrayRCP<double> target_centers;
152  Teuchos::ArrayRCP<GO> target_support_ids;
153  getNodeCoordsAndIds( range_space, d_range_entity_dim, target_centers,
154  target_support_ids );
155 
156  // Calculate an approximate neighborhood distance for the local target
157  // centers. If using kNN, compute an approximation. If doing a radial
158  // search, use the radius. We will use these distances to expand the local
159  // bounding box to ensure we find all of our neighbors in parallel.
160  double target_proximity = 0.0;
161  if ( d_use_knn )
162  {
163  // Get the local bounding box.
164  Teuchos::Tuple<double, 6> local_box;
165  range_space->entitySet()->localBoundingBox( local_box );
166 
167  // Calculate the largest span of the cardinal directions.
168  target_proximity = local_box[3] - local_box[0];
169  for ( int d = 1; d < DIM; ++d )
170  {
171  target_proximity =
172  std::max( target_proximity, local_box[d + 3] - local_box[d] );
173  }
174 
175  // Take the proximity to be 10% of the largest distance.
176  target_proximity *= 0.1;
177  }
178  else
179  {
180  target_proximity = d_radius;
181  }
182 
183  // Gather the source centers that are in the proximity of the target
184  // centers on this proc.
185  Teuchos::Array<double> dist_sources;
186  CenterDistributor<DIM> distributor( comm, source_centers(),
187  target_centers(), target_proximity,
188  dist_sources );
189 
190  // Gather the global ids of the source centers that are within the proximity
191  // of
192  // the target centers on this proc.
193  Teuchos::Array<GO> dist_source_support_ids( distributor.getNumImports() );
194  Teuchos::ArrayView<const GO> source_support_ids_view = source_support_ids();
195  distributor.distribute( source_support_ids_view,
196  dist_source_support_ids() );
197 
198  // Build the source/target pairings.
199  SplineInterpolationPairing<DIM> pairings( dist_sources, target_centers(),
200  d_use_knn, d_knn, d_radius );
201 
202  // Build the basis.
203  Teuchos::RCP<Basis> basis = BP::create();
204 
205  // Build the interpolation matrix.
206  Teuchos::ArrayRCP<SupportId> children_per_parent =
207  pairings.childrenPerParent();
208  SupportId max_entries_per_row = *std::max_element(
209  children_per_parent.begin(), children_per_parent.end() );
210  d_coupling_matrix = Teuchos::rcp( new Tpetra::CrsMatrix<Scalar, LO, GO>(
211  range_map, max_entries_per_row ) );
212  Teuchos::ArrayView<const double> target_view;
213  Teuchos::Array<GO> indices( max_entries_per_row );
214  Teuchos::ArrayView<const double> values;
215  Teuchos::ArrayView<const unsigned> pair_gids;
216  int nn = 0;
217  int local_num_tgt = target_support_ids.size();
218  for ( int i = 0; i < local_num_tgt; ++i )
219  {
220  // If there is no support for this target center then do not build a
221  // local basis.
222  if ( 0 < pairings.childCenterIds( i ).size() )
223  {
224  // Get a view of this target center.
225  target_view = target_centers( i * DIM, DIM );
226 
227  // Build the local interpolation problem.
228  LocalMLSProblem<Basis, DIM> local_problem(
229  target_view, pairings.childCenterIds( i ), dist_sources, *basis,
230  pairings.parentSupportRadius( i ) );
231 
232  // Get MLS shape function values for this target point.
233  values = local_problem.shapeFunction();
234  nn = values.size();
235 
236  // Populate the interpolation matrix row.
237  pair_gids = pairings.childCenterIds( i );
238  for ( int j = 0; j < nn; ++j )
239  {
240  indices[j] = dist_source_support_ids[pair_gids[j]];
241  }
242  d_coupling_matrix->insertGlobalValues( target_support_ids[i],
243  indices( 0, nn ), values );
244  }
245  }
246  d_coupling_matrix->fillComplete( domain_map, range_map );
247  DTK_ENSURE( d_coupling_matrix->isFillComplete() );
248 }
249 
250 //---------------------------------------------------------------------------//
251 // Apply the operator.
252 template <class Basis, int DIM>
254  const TpetraMultiVector &X, TpetraMultiVector &Y, Teuchos::ETransp mode,
255  double alpha, double beta ) const
256 {
257  d_coupling_matrix->apply( X, Y, mode, alpha, beta );
258 }
259 
260 //---------------------------------------------------------------------------//
261 // Transpose apply option.
262 template <class Basis, int DIM>
264  DIM>::hasTransposeApplyImpl() const
265 {
266  return true;
267 }
268 
269 //---------------------------------------------------------------------------//
270 // Extract node coordinates and ids from an iterator.
271 template <class Basis, int DIM>
273  const Teuchos::RCP<FunctionSpace> &space, const int entity_dim,
274  Teuchos::ArrayRCP<double> &centers,
275  Teuchos::ArrayRCP<GO> &support_ids ) const
276 {
277  // Get an iterator over the local nodes.
278  EntityIterator iterator;
279  if ( Teuchos::nonnull( space->entitySet() ) )
280  {
281  LocalEntityPredicate local_predicate(
282  space->entitySet()->communicator()->getRank() );
283  PredicateFunction predicate = PredicateComposition::And(
284  space->selectFunction(), local_predicate.getFunction() );
285  iterator = space->entitySet()->entityIterator( entity_dim, predicate );
286  }
287 
288  // Extract the coordinates and support ids of the nodes.
289  int local_num_node = iterator.size();
290  centers = Teuchos::ArrayRCP<double>( DIM * local_num_node );
291  support_ids = Teuchos::ArrayRCP<GO>( local_num_node );
292  Teuchos::Array<SupportId> node_supports;
293  EntityIterator begin = iterator.begin();
294  EntityIterator end = iterator.end();
295  int entity_counter = 0;
296  for ( EntityIterator entity = begin; entity != end;
297  ++entity, ++entity_counter )
298  {
299  space->shapeFunction()->entitySupportIds( *entity, node_supports );
300  DTK_CHECK( 1 == node_supports.size() );
301  support_ids[entity_counter] = node_supports[0];
302  space->localMap()->centroid( *entity,
303  centers( DIM * entity_counter, DIM ) );
304  }
305 }
306 
307 //---------------------------------------------------------------------------//
308 
309 } // end namespace DataTransferKit
310 
311 //---------------------------------------------------------------------------//
312 
313 #endif // end DTK_MOVINGLEASTSQUARERECONSTRUCTIONOPERATOR_IMPL_HPP
314 
315 //---------------------------------------------------------------------------//
316 // end DTK_MovingLeastSquareReconstructionOperator_impl.hpp
317 //---------------------------------------------------------------------------//
MovingLeastSquareReconstructionOperator(const Teuchos::RCP< const TpetraMap > &domain_map, const Teuchos::RCP< const TpetraMap > &range_map, const Teuchos::ParameterList &parameters)
Constructor.
Local moving least square problem.
Local source/parent center pairings.
Parallel moving least square interpolator.
Parallel moving least square interpolator MapOperator implementation.
void applyImpl(const TpetraMultiVector &X, TpetraMultiVector &Y, Teuchos::ETransp mode=Teuchos::NO_TRANS, double alpha=Teuchos::ScalarTraits< double >::one(), double beta=Teuchos::ScalarTraits< double >::zero()) const override
Apply the operator.
Entity iterator interface.
unsigned long int SupportId
Support id type.
Definition: DTK_Types.hpp:57
static Teuchos::RCP< RadialBasis > create()
Creation method.
Assertions and Design-by-Contract for error handling.
std::function< bool(Entity)> PredicateFunction
Predicate function typedef.
Definition: DTK_Types.hpp:64
Teuchos::ArrayView< const unsigned > childCenterIds(const unsigned parent_id) const
Given a parent center local id get the ids of the child centers within the given radius.
Global center distributor.
Local moving least square problem about a single target center using quadratic polynomials.
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.