DataTransferKit - Multiphysics Solution Transfer Services  2.0
DTK_SplineCoefficientMatrix_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_SPLINECOEFFICIENTMATRIX_IMPL_HPP
42 #define DTK_SPLINECOEFFICIENTMATRIX_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>> &operator_map,
59  const Teuchos::ArrayView<const double> &source_centers,
60  const Teuchos::ArrayView<const SupportId> &source_center_gids,
61  const Teuchos::ArrayView<const double> &dist_source_centers,
62  const Teuchos::ArrayView<const SupportId> &dist_source_center_gids,
63  const SplineInterpolationPairing<DIM> &source_pairings, const Basis &basis )
64 {
65  DTK_CHECK( 0 == source_centers.size() % DIM );
66  DTK_CHECK( source_centers.size() / DIM == source_center_gids.size() );
67  DTK_CHECK( 0 == dist_source_centers.size() % DIM );
68  DTK_CHECK( dist_source_centers.size() / DIM ==
69  dist_source_center_gids.size() );
70 
71  // Get the number of source centers.
72  unsigned num_source_centers = source_center_gids.size();
73 
74  // Create the P matrix.
75  int offset = DIM + 1;
76  Teuchos::RCP<Tpetra::MultiVector<double, int, SupportId>> P_vec =
77  Tpetra::createMultiVector<double, int, SupportId>( operator_map,
78  offset );
79  int di = 0;
80  for ( unsigned i = 0; i < num_source_centers; ++i )
81  {
82  P_vec->replaceGlobalValue( source_center_gids[i], 0, 1.0 );
83  di = DIM * i;
84  for ( int d = 0; d < DIM; ++d )
85  {
86  P_vec->replaceGlobalValue( source_center_gids[i], d + 1,
87  source_centers[di + d] );
88  }
89  }
90  d_P = Teuchos::rcp(
91  new PolynomialMatrix( P_vec, operator_map, operator_map ) );
92 
93  // Create the M matrix.
94  Teuchos::ArrayRCP<SupportId> children_per_parent =
95  source_pairings.childrenPerParent();
96  SupportId max_entries_per_row = *std::max_element(
97  children_per_parent.begin(), children_per_parent.end() );
98  d_M = Teuchos::rcp( new Tpetra::CrsMatrix<double, int, SupportId>(
99  operator_map, max_entries_per_row ) );
100  Teuchos::Array<SupportId> M_indices( max_entries_per_row );
101  Teuchos::Array<double> values( max_entries_per_row );
102  int dj = 0;
103  Teuchos::ArrayView<const unsigned> source_neighbors;
104  double dist = 0.0;
105  int nsn = 0;
106  double radius = 0.0;
107  for ( unsigned i = 0; i < num_source_centers; ++i )
108  {
109  // Get the source points neighboring this source point.
110  di = DIM * i;
111  source_neighbors = source_pairings.childCenterIds( i );
112  nsn = source_neighbors.size();
113  radius = source_pairings.parentSupportRadius( i );
114 
115  // Add the local basis contributions.
116  for ( int j = 0; j < nsn; ++j )
117  {
118  dj = DIM * source_neighbors[j];
119  M_indices[j] = dist_source_center_gids[source_neighbors[j]];
120 
121  dist = EuclideanDistance<DIM>::distance( &source_centers[di],
122  &dist_source_centers[dj] );
123 
124  values[j] = BP::evaluateValue( basis, radius, dist );
125  }
126  d_M->insertGlobalValues( source_center_gids[i], M_indices( 0, nsn ),
127  values( 0, nsn ) );
128  }
129  d_M->fillComplete();
130 
131  DTK_ENSURE( d_M->isFillComplete() );
132 }
133 
134 //---------------------------------------------------------------------------//
135 
136 } // end namespace DataTransferKit
137 
138 //---------------------------------------------------------------------------//
139 
140 #endif // end DTK_SPLINECOEFFICIENTMATRIX_IMPL_HPP
141 
142 //---------------------------------------------------------------------------//
143 // end DTK_SplineCoefficientMatrix_impl.hpp
144 //---------------------------------------------------------------------------//
Spline interpolation coefficient matrix.
SplineCoefficientMatrix(const Teuchos::RCP< const Tpetra::Map< int, SupportId >> &operator_map, const Teuchos::ArrayView< const double > &source_centers, const Teuchos::ArrayView< const SupportId > &source_center_gids, const Teuchos::ArrayView< const double > &dist_source_centers, const Teuchos::ArrayView< const SupportId > &dist_source_center_gids, const SplineInterpolationPairing< DIM > &source_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.
DTK_BasicEntitySet.cpp.