DataTransferKit - Multiphysics Solution Transfer Services  2.0
DTK_MovingLeastSquareReconstructionOperator.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_HPP
42 #define DTK_MOVINGLEASTSQUARERECONSTRUCTIONOPERATOR_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_CrsMatrix.hpp>
53 
54 namespace DataTransferKit
55 {
56 //---------------------------------------------------------------------------//
62 //---------------------------------------------------------------------------//
63 template <class Basis, int DIM>
64 class MovingLeastSquareReconstructionOperator : virtual public MapOperator
65 {
66  public:
68  typedef MapOperator Base;
70  typedef typename Base::Root Root;
71  typedef typename Root::scalar_type Scalar;
72  typedef typename Root::local_ordinal_type LO;
73  typedef typename Root::global_ordinal_type GO;
74  typedef typename Base::TpetraMultiVector TpetraMultiVector;
75  typedef typename Base::TpetraMap TpetraMap;
78 
89  const Teuchos::RCP<const TpetraMap> &domain_map,
90  const Teuchos::RCP<const TpetraMap> &range_map,
91  const Teuchos::ParameterList &parameters );
92 
93  protected:
94  /*
95  * \brief Setup the map operator from a domain entity set and a range
96  * entity set.
97  *
98  * \param domain_function The function that contains the data that will be
99  * sent to the range. Must always be nonnull but the pointers it contains
100  * may be null of no entities are on-process.
101  *
102  * \param range_space The function that will receive the data from the
103  * domain. Must always be nonnull but the pointers it contains to entity
104  * data may be null of no entities are on-process.
105  *
106  * \param parameters Parameters for the setup.
107  */
108  void setupImpl( const Teuchos::RCP<FunctionSpace> &domain_space,
109  const Teuchos::RCP<FunctionSpace> &range_space ) override;
110 
114  void applyImpl(
115  const TpetraMultiVector &X, TpetraMultiVector &Y,
116  Teuchos::ETransp mode = Teuchos::NO_TRANS,
117  double alpha = Teuchos::ScalarTraits<double>::one(),
118  double beta = Teuchos::ScalarTraits<double>::zero() ) const override;
119 
120  /*
121  * \brief Transpose apply option.
122  */
123  bool hasTransposeApplyImpl() const override;
124 
125  private:
126  // Extract node coordinates and ids from an iterator.
127  void getNodeCoordsAndIds( const Teuchos::RCP<FunctionSpace> &space,
128  const int entity_dim,
129  Teuchos::ArrayRCP<double> &centers,
130  Teuchos::ArrayRCP<GO> &support_ids ) const;
131 
132  private:
133  // Flag for search type. True if kNN, false if radius.
134  bool d_use_knn;
135 
136  // k-nearest-neighbors for support.
137  int d_knn;
138 
139  // Basis radius for support.
140  double d_radius;
141 
142  // Domain entity topological dimension. Default is 0 (vertex).
143  int d_domain_entity_dim;
144 
145  // Range entity topological dimension. Default is 0 (vertex).
146  int d_range_entity_dim;
147 
148  // Coupling matrix.
149  Teuchos::RCP<Tpetra::CrsMatrix<Scalar, LO, GO>> d_coupling_matrix;
150 };
151 
152 //---------------------------------------------------------------------------//
153 
154 } // end namespace DataTransferKit
155 
156 //---------------------------------------------------------------------------//
157 
158 #endif // end DTK_MOVINGLEASTSQUARERECONSTRUCTIONOPERATOR_HPP
159 
160 //---------------------------------------------------------------------------//
161 // end DTK_MovingLeastSquareReconstructionOperator.hpp
162 //---------------------------------------------------------------------------//
MovingLeastSquareReconstructionOperator(const Teuchos::RCP< const TpetraMap > &domain_map, const Teuchos::RCP< const TpetraMap > &range_map, const Teuchos::ParameterList &parameters)
Constructor.
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.
Policy class for spline interpolation basis functions.
DTK_BasicEntitySet.cpp.