DataTransferKit - Multiphysics Solution Transfer Services  2.0
DTK_SplineInterpolationOperator.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_HPP
42 #define DTK_SPLINEINTERPOLATIONOPERATOR_HPP
43 
44 #include "DTK_MapOperator.hpp"
46 
47 #include <Teuchos_Array.hpp>
48 #include <Teuchos_ArrayView.hpp>
49 #include <Teuchos_Comm.hpp>
50 #include <Teuchos_RCP.hpp>
51 
52 #include <Tpetra_Map.hpp>
53 
54 #include <Thyra_LinearOpBase.hpp>
55 
56 namespace DataTransferKit
57 {
58 //---------------------------------------------------------------------------//
67 //---------------------------------------------------------------------------//
68 template <class Basis, int DIM>
69 class SplineInterpolationOperator : virtual public MapOperator
70 {
71  public:
73  typedef MapOperator Base;
75  typedef typename Base::Root Root;
76  typedef typename Root::scalar_type Scalar;
77  typedef typename Root::local_ordinal_type LO;
78  typedef typename Root::global_ordinal_type GO;
79  typedef typename Base::TpetraMultiVector TpetraMultiVector;
80  typedef typename Base::TpetraMap TpetraMap;
83 
84  /*
85  * \brief Constructor.
86  *
87  * \param domain_map Parallel map for domain vectors this map should be
88  * compatible with.
89  *
90  * \param range_map Parallel map for range vectors this map should be
91  * compatible with.
92  */
94  const Teuchos::RCP<const TpetraMap> &domain_map,
95  const Teuchos::RCP<const TpetraMap> &range_map,
96  const Teuchos::ParameterList &parameters );
97 
98  protected:
99  /*
100  * \brief Setup the map operator from a domain entity set and a range
101  * entity set.
102  *
103  * \param domain_function The function that contains the data that will be
104  * sent to the range. Must always be nonnull but the pointers it contains
105  * may be null of no entities are on-process.
106  *
107  * \param range_space The function that will receive the data from the
108  * domain. Must always be nonnull but the pointers it contains to entity
109  * data may be null of no entities are on-process.
110  *
111  * \param parameters Parameters for the setup.
112  */
113  void setupImpl( const Teuchos::RCP<FunctionSpace> &domain_space,
114  const Teuchos::RCP<FunctionSpace> &range_space ) override;
115 
119  void applyImpl(
120  const TpetraMultiVector &X, TpetraMultiVector &Y,
121  Teuchos::ETransp mode = Teuchos::NO_TRANS,
122  double alpha = Teuchos::ScalarTraits<double>::one(),
123  double beta = Teuchos::ScalarTraits<double>::zero() ) const override;
124 
125  /*
126  * \brief Transpose apply option.
127  */
128  bool hasTransposeApplyImpl() const override;
129 
130  private:
131  // Build the concrete operators.
132  void buildConcreteOperators(
133  const Teuchos::RCP<FunctionSpace> &domain_space,
134  const Teuchos::RCP<FunctionSpace> &range_space,
135  Teuchos::RCP<const Root> &S, Teuchos::RCP<const Root> &P,
136  Teuchos::RCP<const Root> &M, Teuchos::RCP<const Root> &Q,
137  Teuchos::RCP<const Root> &N ) const;
138 
139  private:
140  // Extract node coordinates and ids from an iterator.
141  void getNodeCoordsAndIds( const Teuchos::RCP<FunctionSpace> &space,
142  const int entity_dim,
143  Teuchos::ArrayRCP<double> &centers,
144  Teuchos::ArrayRCP<GO> &support_ids ) const;
145 
146  private:
147  // Flag for search type. True if kNN, false if radius.
148  bool d_use_knn;
149 
150  // k-nearest-neighbors for support.
151  int d_knn;
152 
153  // Basis radius.
154  double d_radius;
155 
156  // Domain entity topological dimension. Default is 0 (vertex).
157  int d_domain_entity_dim;
158 
159  // Range entity topological dimension. Default is 0 (vertex).
160  int d_range_entity_dim;
161 
162  // Stratimikos parameter list.
163  Teuchos::RCP<Teuchos::ParameterList> d_stratimikos_list;
164 
165  // Coupling matrix.
166  Teuchos::RCP<const Thyra::LinearOpBase<double>> d_coupling_matrix;
167 };
168 
169 //---------------------------------------------------------------------------//
170 
171 } // end namespace DataTransferKit
172 
173 //---------------------------------------------------------------------------//
174 
175 #endif // end DTK_SPLINEINTERPOLATIONOPERATOR_HPP
176 
177 //---------------------------------------------------------------------------//
178 // end DTK_SplineInterpolationOperator.hpp
179 //---------------------------------------------------------------------------//
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.
Policy class for spline interpolation basis functions.
DTK_BasicEntitySet.cpp.