DataTransferKit - Multiphysics Solution Transfer Services  2.0
DTK_SplineInterpolationOperator_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_SPLINEINTERPOLATIONOPERATOR_IMPL_HPP
42 #define DTK_SPLINEINTERPOLATIONOPERATOR_IMPL_HPP
43 
44 #include "DTK_BasicEntityPredicates.hpp"
46 #include "DTK_DBC.hpp"
47 #include "DTK_PredicateComposition.hpp"
53 
54 #include <Teuchos_ArrayRCP.hpp>
55 #include <Teuchos_CommHelpers.hpp>
56 #include <Teuchos_ParameterList.hpp>
57 #include <Teuchos_Ptr.hpp>
58 #include <Teuchos_XMLParameterListCoreHelpers.hpp>
59 
60 #include <BelosPseudoBlockGmresSolMgr.hpp>
61 
62 #include <Thyra_DefaultAddedLinearOp.hpp>
63 #include <Thyra_DefaultMultipliedLinearOp.hpp>
64 #include <Thyra_DefaultScaledAdjointLinearOp.hpp>
65 #include <Thyra_LinearOpWithSolveFactoryHelpers.hpp>
66 #include <Thyra_TpetraThyraWrappers.hpp>
67 
68 #include <Stratimikos_DefaultLinearSolverBuilder.hpp>
69 
70 namespace DataTransferKit
71 {
72 //---------------------------------------------------------------------------//
73 // Constructor.
74 template <class Basis, int DIM>
75 SplineInterpolationOperator<Basis, DIM>::SplineInterpolationOperator(
76  const Teuchos::RCP<const TpetraMap> &domain_map,
77  const Teuchos::RCP<const TpetraMap> &range_map,
78  const Teuchos::ParameterList &parameters )
79  : Base( domain_map, range_map )
80  , d_use_knn( false )
81  , d_knn( 0 )
82  , d_radius( 0.0 )
83  , d_domain_entity_dim( 0 )
84  , d_range_entity_dim( 0 )
85 {
86  // Determine if we are doing kNN search or radius search.
87  if ( parameters.isParameter( "Type of Search" ) )
88  {
89  if ( "Radius" == parameters.get<std::string>( "Type of Search" ) )
90  {
91  d_use_knn = false;
92  }
93  else if ( "Nearest Neighbor" ==
94  parameters.get<std::string>( "Type of Search" ) )
95  {
96  d_use_knn = true;
97  }
98  else
99  {
100  // Otherwise we got an invalid search type.
101  DTK_INSIST( false );
102  }
103  }
104 
105  // If we are doing kNN support get the number of neighbors.
106  if ( d_use_knn )
107  {
108  DTK_REQUIRE( parameters.isParameter( "Num Neighbors" ) );
109  d_knn = parameters.get<int>( "Num Neighbors" );
110  }
111 
112  // Otherwise we are doing the radius search so get the basis radius.
113  else
114  {
115  DTK_REQUIRE( parameters.isParameter( "RBF Radius" ) );
116  d_radius = parameters.get<double>( "RBF Radius" );
117  }
118 
119  // Get the topological dimension of the domain and range entities. This
120  // map will use their centroids for the point cloud.
121  if ( parameters.isParameter( "Domain Entity Dimension" ) )
122  {
123  d_domain_entity_dim = parameters.get<int>( "Domain Entity Dimension" );
124  }
125  if ( parameters.isParameter( "Range Entity Dimension" ) )
126  {
127  d_range_entity_dim = parameters.get<int>( "Range Entity Dimension" );
128  }
129 
130  // Get the stratimikos parameters if they exist.
131  if ( parameters.isSublist( "Stratimikos" ) )
132  {
133  d_stratimikos_list = Teuchos::rcp(
134  new Teuchos::ParameterList( parameters.sublist( "Stratimikos" ) ) );
135  }
136 }
137 
138 //---------------------------------------------------------------------------//
139 // Setup the map operator.
140 template <class Basis, int DIM>
141 void SplineInterpolationOperator<Basis, DIM>::setupImpl(
142  const Teuchos::RCP<FunctionSpace> &domain_space,
143  const Teuchos::RCP<FunctionSpace> &range_space )
144 {
145  DTK_REQUIRE( Teuchos::nonnull( domain_space ) );
146  DTK_REQUIRE( Teuchos::nonnull( range_space ) );
147 
148  // Extract the Support maps.
149  const Teuchos::RCP<const typename Base::TpetraMap> domain_map =
150  this->getDomainMap();
151  const Teuchos::RCP<const typename Base::TpetraMap> range_map =
152  this->getRangeMap();
153 
154  // Prolongation operator.
155  Teuchos::RCP<const Root> S;
156 
157  // Coefficient matrix polynomial component.
158  Teuchos::RCP<const Root> P;
159 
160  // Coefficient matrix basis component.
161  Teuchos::RCP<const Root> M;
162 
163  // Evaluation matrix polynomial component.
164  Teuchos::RCP<const Root> Q;
165 
166  // Evaluation matrix basis component.
167  Teuchos::RCP<const Root> N;
168 
169  // Build the concrete operators.
170  buildConcreteOperators( domain_space, range_space, S, P, M, Q, N );
171 
172  // Create an abstract wrapper for S.
173  Teuchos::RCP<const Thyra::VectorSpaceBase<Scalar>>
174  thyra_range_vector_space_S =
175  Thyra::createVectorSpace<Scalar>( S->getRangeMap() );
176  Teuchos::RCP<const Thyra::VectorSpaceBase<Scalar>>
177  thyra_domain_vector_space_S =
178  Thyra::createVectorSpace<Scalar>( S->getDomainMap() );
179  Teuchos::RCP<const Thyra::TpetraLinearOp<Scalar, LO, GO>> thyra_S =
180  Teuchos::rcp( new Thyra::TpetraLinearOp<Scalar, LO, GO>() );
181  Teuchos::rcp_const_cast<Thyra::TpetraLinearOp<Scalar, LO, GO>>( thyra_S )
182  ->constInitialize( thyra_range_vector_space_S,
183  thyra_domain_vector_space_S, S );
184 
185  // Create an abstract wrapper for P.
186  Teuchos::RCP<const Thyra::VectorSpaceBase<Scalar>>
187  thyra_range_vector_space_P =
188  Thyra::createVectorSpace<Scalar>( P->getRangeMap() );
189  Teuchos::RCP<const Thyra::VectorSpaceBase<Scalar>>
190  thyra_domain_vector_space_P =
191  Thyra::createVectorSpace<Scalar>( P->getDomainMap() );
192  Teuchos::RCP<const Thyra::TpetraLinearOp<Scalar, LO, GO>> thyra_P =
193  Teuchos::rcp( new Thyra::TpetraLinearOp<Scalar, LO, GO>() );
194  Teuchos::rcp_const_cast<Thyra::TpetraLinearOp<Scalar, LO, GO>>( thyra_P )
195  ->constInitialize( thyra_range_vector_space_P,
196  thyra_domain_vector_space_P, P );
197 
198  // Create an abstract wrapper for M.
199  Teuchos::RCP<const Thyra::VectorSpaceBase<Scalar>>
200  thyra_range_vector_space_M =
201  Thyra::createVectorSpace<Scalar>( M->getRangeMap() );
202  Teuchos::RCP<const Thyra::VectorSpaceBase<Scalar>>
203  thyra_domain_vector_space_M =
204  Thyra::createVectorSpace<Scalar>( M->getDomainMap() );
205  Teuchos::RCP<const Thyra::TpetraLinearOp<Scalar, LO, GO>> thyra_M =
206  Teuchos::rcp( new Thyra::TpetraLinearOp<Scalar, LO, GO>() );
207  Teuchos::rcp_const_cast<Thyra::TpetraLinearOp<Scalar, LO, GO>>( thyra_M )
208  ->constInitialize( thyra_range_vector_space_M,
209  thyra_domain_vector_space_M, M );
210 
211  // Create an abstract wrapper for Q.
212  Teuchos::RCP<const Thyra::VectorSpaceBase<Scalar>>
213  thyra_range_vector_space_Q =
214  Thyra::createVectorSpace<Scalar>( Q->getRangeMap() );
215  Teuchos::RCP<const Thyra::VectorSpaceBase<Scalar>>
216  thyra_domain_vector_space_Q =
217  Thyra::createVectorSpace<Scalar>( Q->getDomainMap() );
218  Teuchos::RCP<const Thyra::TpetraLinearOp<Scalar, LO, GO>> thyra_Q =
219  Teuchos::rcp( new Thyra::TpetraLinearOp<Scalar, LO, GO>() );
220  Teuchos::rcp_const_cast<Thyra::TpetraLinearOp<Scalar, LO, GO>>( thyra_Q )
221  ->constInitialize( thyra_range_vector_space_Q,
222  thyra_domain_vector_space_Q, Q );
223 
224  // Create an abstract wrapper for N.
225  Teuchos::RCP<const Thyra::VectorSpaceBase<Scalar>>
226  thyra_range_vector_space_N =
227  Thyra::createVectorSpace<Scalar>( N->getRangeMap() );
228  Teuchos::RCP<const Thyra::VectorSpaceBase<Scalar>>
229  thyra_domain_vector_space_N =
230  Thyra::createVectorSpace<Scalar>( N->getDomainMap() );
231  Teuchos::RCP<const Thyra::TpetraLinearOp<Scalar, LO, GO>> thyra_N =
232  Teuchos::rcp( new Thyra::TpetraLinearOp<Scalar, LO, GO>() );
233  Teuchos::rcp_const_cast<Thyra::TpetraLinearOp<Scalar, LO, GO>>( thyra_N )
234  ->constInitialize( thyra_range_vector_space_N,
235  thyra_domain_vector_space_N, N );
236 
237  // COUPLING MATRIX ASSEMBLY: A = (Q + N)*[(P + M + P^T)^-1]*S
238  // Create a transpose of P.
239  Teuchos::RCP<const Thyra::LinearOpBase<Scalar>> thyra_P_T =
240  Thyra::transpose<Scalar>( thyra_P );
241 
242  // Create a composite operator C = (P + M + P^T)
243  Teuchos::RCP<const Thyra::LinearOpBase<Scalar>> thyra_PpM =
244  Thyra::add<Scalar>( thyra_P, thyra_M );
245  Teuchos::RCP<const Thyra::LinearOpBase<Scalar>> thyra_C =
246  Thyra::add<Scalar>( thyra_PpM, thyra_P_T );
247 
248  // If we didnt get stratimikos parameters from the input list, create some
249  // here.
250  if ( Teuchos::is_null( d_stratimikos_list ) )
251  {
252  d_stratimikos_list = Teuchos::parameterList( "Stratimikos" );
253 
254  d_stratimikos_list->set( "Linear Solver Type", "Belos" );
255  d_stratimikos_list->set( "Preconditioner Type", "None" );
256 
257  auto &linear_solver_types_list =
258  d_stratimikos_list->sublist( "Linear Solver Types" );
259  auto &belos_list = linear_solver_types_list.sublist( "Belos" );
260  belos_list.set( "Solver Type", "Pseudo Block GMRES" );
261  auto &solver_types_list = belos_list.sublist( "Solver Types" );
262  auto &gmres_list = solver_types_list.sublist( "Pseudo Block GMRES" );
263  gmres_list.set( "Convergence Tolerance", 1.0e-10 );
264  gmres_list.set( "Verbosity",
265  Belos::Errors + Belos::Warnings + Belos::TimingDetails +
266  Belos::FinalSummary + Belos::StatusTestDetails );
267  gmres_list.set( "Output Frequency", 1 );
268  }
269 
270  // Create the inverse of the composite operator C.
271  Stratimikos::DefaultLinearSolverBuilder builder;
272  builder.setParameterList( d_stratimikos_list );
273  Teuchos::RCP<Thyra::LinearOpWithSolveFactoryBase<Scalar>> factory =
274  Thyra::createLinearSolveStrategy( builder );
275  Teuchos::RCP<const Thyra::LinearOpBase<Scalar>> thyra_C_inv =
276  Thyra::inverse<Scalar>( *factory, thyra_C );
277 
278  // Create the composite operator B = (Q + N);
279  Teuchos::RCP<const Thyra::LinearOpBase<Scalar>> thyra_B =
280  Thyra::add<Scalar>( thyra_Q, thyra_N );
281 
282  // Create the coupling matrix A = (B * C^-1 * S).
283  d_coupling_matrix =
284  Thyra::multiply<Scalar>( thyra_B, thyra_C_inv, thyra_S );
285  DTK_ENSURE( Teuchos::nonnull( d_coupling_matrix ) );
286 }
287 
288 //---------------------------------------------------------------------------//
289 // Apply the operator.
290 template <class Basis, int DIM>
292  const TpetraMultiVector &X, TpetraMultiVector &Y, Teuchos::ETransp mode,
293  Scalar alpha, Scalar beta ) const
294 {
295  DTK_REQUIRE( Teuchos::NO_TRANS == mode );
296  Teuchos::RCP<const Thyra::MultiVectorBase<Scalar>> thyra_X =
297  Thyra::createConstMultiVector<Scalar>( Teuchos::rcpFromRef( X ) );
298  Teuchos::RCP<Thyra::MultiVectorBase<Scalar>> thyra_Y =
299  Thyra::createMultiVector<Scalar>( Teuchos::rcpFromRef( Y ) );
300  d_coupling_matrix->apply( Thyra::NOTRANS, *thyra_X, thyra_Y.ptr(), alpha,
301  beta );
302 }
303 
304 //---------------------------------------------------------------------------//
305 // Transpose apply option.
306 template <class Basis, int DIM>
308 {
309  return false;
310 }
311 
312 //---------------------------------------------------------------------------//
316 template <class Basis, int DIM>
318  const Teuchos::RCP<FunctionSpace> &domain_space,
319  const Teuchos::RCP<FunctionSpace> &range_space, Teuchos::RCP<const Root> &S,
320  Teuchos::RCP<const Root> &P, Teuchos::RCP<const Root> &M,
321  Teuchos::RCP<const Root> &Q, Teuchos::RCP<const Root> &N ) const
322 {
323  // Extract the Support maps.
324  const Teuchos::RCP<const typename Base::TpetraMap> domain_map =
325  this->getDomainMap();
326  const Teuchos::RCP<const typename Base::TpetraMap> range_map =
327  this->getRangeMap();
328 
329  // Determine if we have range and domain data on this process.
330  bool nonnull_domain = Teuchos::nonnull( domain_space->entitySet() );
331  bool nonnull_range = Teuchos::nonnull( range_space->entitySet() );
332 
333  // Get the parallel communicator.
334  Teuchos::RCP<const Teuchos::Comm<int>> comm = domain_map->getComm();
335 
336  // Extract the source nodes and their ids.
337  Teuchos::ArrayRCP<double> source_centers;
338  Teuchos::ArrayRCP<GO> source_support_ids;
339  getNodeCoordsAndIds( domain_space, d_domain_entity_dim, source_centers,
340  source_support_ids );
341 
342  // Extract the target nodes and their ids.
343  Teuchos::ArrayRCP<double> target_centers;
344  Teuchos::ArrayRCP<GO> target_support_ids;
345  getNodeCoordsAndIds( range_space, d_range_entity_dim, target_centers,
346  target_support_ids );
347 
348  // Calculate an approximate neighborhood distance for the local source
349  // centers. If using kNN, compute an approximation. If doing a radial
350  // search, use the radius.
351  double source_proximity = 0.0;
352  if ( d_use_knn )
353  {
354  // Get the local bounding box.
355  Teuchos::Tuple<double, 6> local_box;
356  domain_space->entitySet()->localBoundingBox( local_box );
357 
358  // Calculate the largest span of the cardinal directions.
359  source_proximity = local_box[3] - local_box[0];
360  for ( int d = 1; d < DIM; ++d )
361  {
362  source_proximity =
363  std::max( source_proximity, local_box[d + 3] - local_box[d] );
364  }
365 
366  // Take the proximity to be 10% of the largest distance.
367  source_proximity *= 0.1;
368  }
369  else
370  {
371  source_proximity = d_radius;
372  }
373 
374  // Gather the source centers that are in the proximity of the source
375  // centers on this proc.
376  Teuchos::Array<double> dist_sources;
377  CenterDistributor<DIM> source_distributor( comm, source_centers(),
378  source_centers(),
379  source_proximity, dist_sources );
380 
381  // Distribute the global source ids.
382  Teuchos::Array<GO> dist_source_support_ids(
383  source_distributor.getNumImports() );
384  Teuchos::ArrayView<const GO> source_support_ids_view = source_support_ids();
385  Teuchos::ArrayView<GO> dist_gids_view = dist_source_support_ids();
386  source_distributor.distribute( source_support_ids_view, dist_gids_view );
387 
388  // Build the source/source pairings.
389  SplineInterpolationPairing<DIM> source_pairings(
390  dist_sources(), source_centers(), d_use_knn, d_knn, d_radius );
391 
392  // Build the basis.
393  Teuchos::RCP<Basis> basis = BP::create();
394 
395  // PROLONGATION OPERATOR.
396  GO offset = comm->getRank() ? 0 : DIM + 1;
397  S = Teuchos::rcp( new SplineProlongationOperator( offset, domain_map ) );
398 
399  // Get the operator map.
400  Teuchos::RCP<const Tpetra::Map<int, GO>> prolongated_map = S->getRangeMap();
401 
402  // COEFFICIENT OPERATORS.
403  // Build the coefficient operators.
405  prolongated_map, source_centers(), source_support_ids(), dist_sources(),
406  dist_source_support_ids(), source_pairings, *basis );
407  P = C.getP();
408  M = C.getM();
409 
410  // Cleanup.
411  dist_source_support_ids.clear();
412 
413  // EVALUATION OPERATORS.
414  // Calculate an approximate neighborhood distance for the local target
415  // centers. If using kNN, compute an approximation. If doing a radial
416  // search, use the radius.
417  double target_proximity = 0.0;
418  if ( d_use_knn )
419  {
420  // Get the local bounding box.
421  Teuchos::Tuple<double, 6> local_box;
422  range_space->entitySet()->localBoundingBox( local_box );
423 
424  // Calculate the largest span of the cardinal directions.
425  target_proximity = local_box[3] - local_box[0];
426  for ( int d = 1; d < DIM; ++d )
427  {
428  target_proximity =
429  std::max( target_proximity, local_box[d + 3] - local_box[d] );
430  }
431 
432  // Take the proximity to be 10% of the largest distance.
433  target_proximity *= 0.1;
434  }
435  else
436  {
437  target_proximity = d_radius;
438  }
439 
440  // Gather the source centers that are in the proximity of the target
441  // centers on this proc.
442  CenterDistributor<DIM> target_distributor( comm, source_centers(),
443  target_centers(),
444  target_proximity, dist_sources );
445 
446  // Distribute the global source ids.
447  dist_source_support_ids.resize( target_distributor.getNumImports() );
448  source_support_ids_view = source_support_ids();
449  dist_gids_view = dist_source_support_ids();
450  target_distributor.distribute( source_support_ids_view, dist_gids_view );
451 
452  // Build the source/target pairings.
453  SplineInterpolationPairing<DIM> target_pairings(
454  dist_sources(), target_centers(), d_use_knn, d_knn, d_radius );
455 
456  // Build the transformation operators.
458  prolongated_map, range_map, target_centers(), target_support_ids(),
459  dist_sources(), dist_source_support_ids(), target_pairings, *basis );
460  N = B.getN();
461  Q = B.getQ();
462 
463  DTK_ENSURE( Teuchos::nonnull( S ) );
464  DTK_ENSURE( Teuchos::nonnull( P ) );
465  DTK_ENSURE( Teuchos::nonnull( M ) );
466  DTK_ENSURE( Teuchos::nonnull( Q ) );
467  DTK_ENSURE( Teuchos::nonnull( N ) );
468 }
469 
470 //---------------------------------------------------------------------------//
471 // Extract node coordinates and ids from an iterator.
472 template <class Basis, int DIM>
474  const Teuchos::RCP<FunctionSpace> &space, const int entity_dim,
475  Teuchos::ArrayRCP<double> &centers,
476  Teuchos::ArrayRCP<GO> &support_ids ) const
477 {
478  // Get an iterator over the local nodes.
479  EntityIterator iterator;
480  if ( Teuchos::nonnull( space->entitySet() ) )
481  {
482  LocalEntityPredicate local_predicate(
483  space->entitySet()->communicator()->getRank() );
484  PredicateFunction predicate = PredicateComposition::And(
485  space->selectFunction(), local_predicate.getFunction() );
486  iterator = space->entitySet()->entityIterator( entity_dim, predicate );
487  }
488 
489  // Extract the coordinates and support ids of the nodes.
490  int local_num_node = iterator.size();
491  centers = Teuchos::ArrayRCP<double>( DIM * local_num_node );
492  support_ids = Teuchos::ArrayRCP<GO>( local_num_node );
493  Teuchos::Array<SupportId> node_supports;
494  EntityIterator begin = iterator.begin();
495  EntityIterator end = iterator.end();
496  int entity_counter = 0;
497  for ( EntityIterator entity = begin; entity != end;
498  ++entity, ++entity_counter )
499  {
500  space->shapeFunction()->entitySupportIds( *entity, node_supports );
501  DTK_CHECK( 1 == node_supports.size() );
502  support_ids[entity_counter] = node_supports[0];
503  space->localMap()->centroid( *entity,
504  centers( DIM * entity_counter, DIM ) );
505  }
506 }
507 
508 //---------------------------------------------------------------------------//
509 
510 } // end namespace DataTransferKit
511 
512 //---------------------------------------------------------------------------//
513 
514 #endif // end DTK_SPLINEINTERPOLATIONOPERATOR_IMPL_HPP
515 
516 //---------------------------------------------------------------------------//
517 // end DTK_SplineInterpolationOperator_impl.hpp
518 //---------------------------------------------------------------------------//
Local source/parent center pairings.
Sparse spline transformation operator (the A matrix). A = N + Q.
Spline interpolation coefficient matrix.
Entity iterator interface.
static Teuchos::RCP< RadialBasis > create()
Creation method.
Prolongation operator for projecting a vector into the extended spline space.
Assertions and Design-by-Contract for error handling.
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.
Parallel spline interpolator.
std::function< bool(Entity)> PredicateFunction
Predicate function typedef.
Definition: DTK_Types.hpp:64
Spline transformation operator.
Global center distributor.
Spline transformation operator.
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.