DataTransferKit - Multiphysics Solution Transfer Services  2.0
DTK_SplineEvaluationMatrix_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_SPLINEEVALUATIONMATRIX_IMPL_HPP
42 #define DTK_SPLINEEVALUATIONMATRIX_IMPL_HPP
43 
44 #include "DTK_DBC.hpp"
47 
48 #include <Teuchos_Array.hpp>
49 
50 namespace DataTransferKit
51 {
52 //---------------------------------------------------------------------------//
56 template <class Basis, int DIM>
58  const Teuchos::RCP<const Tpetra::Map<int, SupportId>> &domain_map,
59  const Teuchos::RCP<const Tpetra::Map<int, SupportId>> &range_map,
60  const Teuchos::ArrayView<const double> &target_centers,
61  const Teuchos::ArrayView<const SupportId> &target_center_gids,
62  const Teuchos::ArrayView<const double> &dist_source_centers,
63  const Teuchos::ArrayView<const SupportId> &dist_source_center_gids,
64  const SplineInterpolationPairing<DIM> &target_pairings, const Basis &basis )
65 {
66  DTK_CHECK( 0 == target_centers.size() % DIM );
67  DTK_CHECK( target_centers.size() / DIM == target_center_gids.size() );
68  DTK_CHECK( 0 == dist_source_centers.size() % DIM );
69  DTK_CHECK( dist_source_centers.size() / DIM ==
70  dist_source_center_gids.size() );
71 
72  // Get the number of target centers.
73  unsigned num_target_centers = target_center_gids.size();
74 
75  // Create the Q matrix.
76  int offset = DIM + 1;
77  Teuchos::RCP<Tpetra::MultiVector<double, int, SupportId>> Q_vec =
78  Tpetra::createMultiVector<double, int, SupportId>( range_map, offset );
79  int di = 0;
80  for ( unsigned i = 0; i < num_target_centers; ++i )
81  {
82  Q_vec->replaceGlobalValue( target_center_gids[i], 0, 1.0 );
83  di = DIM * i;
84  for ( int d = 0; d < DIM; ++d )
85  {
86  Q_vec->replaceGlobalValue( target_center_gids[i], d + 1,
87  target_centers[di + d] );
88  }
89  }
90  d_Q = Teuchos::rcp( new PolynomialMatrix( Q_vec, domain_map, range_map ) );
91 
92  // Create the N matrix.
93  Teuchos::ArrayRCP<SupportId> children_per_parent =
94  target_pairings.childrenPerParent();
95  SupportId max_entries_per_row = *std::max_element(
96  children_per_parent.begin(), children_per_parent.end() );
97  d_N = Teuchos::rcp( new Tpetra::CrsMatrix<double, int, SupportId>(
98  range_map, max_entries_per_row ) );
99  Teuchos::Array<SupportId> N_indices( max_entries_per_row );
100  Teuchos::Array<double> values( max_entries_per_row );
101  int dj = 0;
102  Teuchos::ArrayView<const unsigned> target_neighbors;
103  double dist = 0.0;
104  int ntn = 0;
105  double radius = 0.0;
106  for ( unsigned i = 0; i < num_target_centers; ++i )
107  {
108  di = DIM * i;
109 
110  // Get the source points neighboring this target point.
111  target_neighbors = target_pairings.childCenterIds( i );
112  ntn = target_neighbors.size();
113  radius = target_pairings.parentSupportRadius( i );
114 
115  // Add the local basis contributions.
116  for ( int j = 0; j < ntn; ++j )
117  {
118  dj = DIM * target_neighbors[j];
119 
120  N_indices[j] = dist_source_center_gids[target_neighbors[j]];
121 
122  dist = EuclideanDistance<DIM>::distance( &target_centers[di],
123  &dist_source_centers[dj] );
124 
125  values[j] = BP::evaluateValue( basis, radius, dist );
126  }
127 
128  d_N->insertGlobalValues( target_center_gids[i], N_indices( 0, ntn ),
129  values( 0, ntn ) );
130  }
131  d_N->fillComplete( domain_map, range_map );
132 
133  DTK_ENSURE( d_N->isFillComplete() );
134 }
135 
136 //---------------------------------------------------------------------------//
137 
138 } // end namespace DataTransferKit
139 
140 //---------------------------------------------------------------------------//
141 
142 #endif // end DTK_SPLINEEVALUATIONMATRIX_IMPL_HPP
143 
144 //---------------------------------------------------------------------------//
145 // end DTK_SplineEvaluationMatrix_impl.hpp
146 //---------------------------------------------------------------------------//
SplineEvaluationMatrix(const Teuchos::RCP< const Tpetra::Map< int, SupportId >> &domain_map, const Teuchos::RCP< const Tpetra::Map< int, SupportId >> &range_map, const Teuchos::ArrayView< const double > &target_centers, const Teuchos::ArrayView< const SupportId > &target_center_gids, const Teuchos::ArrayView< const double > &dist_source_centers, const Teuchos::ArrayView< const SupportId > &dist_source_center_gids, const SplineInterpolationPairing< DIM > &target_pairings, const Basis &basis)
Constructor.
unsigned long int SupportId
Support id type.
Definition: DTK_Types.hpp:57
Assertions and Design-by-Contract for error handling.
Vector apply implementation for polynomial matrices.
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.
Euclidean distance function.
Spline transformation operator.
DTK_BasicEntitySet.cpp.