DataTransferKit - Multiphysics Solution Transfer Services  2.0
DTK_LocalMLSProblem_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_LOCALMLSPROBLEM_IMPL_HPP
42 #define DTK_LOCALMLSPROBLEM_IMPL_HPP
43 
44 #include <algorithm>
45 #include <functional>
46 #include <limits>
47 
48 #include "DTK_DBC.hpp"
50 #include "DTK_LocalMLSProblem.hpp"
52 
53 #include <Teuchos_LAPACK.hpp>
54 #include <Teuchos_SerialDenseSolver.hpp>
55 #include <Teuchos_SerialDenseVector.hpp>
56 
57 namespace DataTransferKit
58 {
59 //---------------------------------------------------------------------------//
63 template <class Basis, int DIM>
65  const Teuchos::ArrayView<const double> &target_center,
66  const Teuchos::ArrayView<const unsigned> &source_lids,
67  const Teuchos::ArrayView<const double> &source_centers, const Basis &basis,
68  const double radius )
69  : d_shape_function( source_lids.size() )
70 {
71  DTK_REQUIRE( 0 == source_centers.size() % DIM );
72  DTK_REQUIRE( 0 == target_center.size() % DIM );
73 
74  // Number of source centers supporting this target center.
75  int num_sources = source_lids.size();
76  int poly_size = 0;
77  Teuchos::SerialDenseMatrix<int, double> P;
78  Teuchos::SerialDenseVector<int, double> target_poly;
79 
80  // Make Phi.
81  Teuchos::SerialDenseMatrix<int, double> phi( num_sources, num_sources );
82  Teuchos::ArrayView<const double> source_center_view;
83  double dist = 0.0;
84  for ( int i = 0; i < num_sources; ++i )
85  {
86  source_center_view = source_centers( DIM * source_lids[i], DIM );
88  target_center.getRawPtr(), source_center_view.getRawPtr() );
89  phi( i, i ) = BP::evaluateValue( basis, radius, dist );
90  }
91 
92  // Make P.
93  Teuchos::Array<int> poly_ids( 1, 0 );
94  int poly_id = 0;
95  P.reshape( num_sources, poly_id + 1 );
96  for ( int i = 0; i < num_sources; ++i )
97  {
98  source_center_view = source_centers( DIM * source_lids[i], DIM );
99  P( i, poly_id ) = polynomialCoefficient( 0, source_center_view );
100  }
101  ++poly_id;
102 
103  // Add polynomial columns until we are full rank.
104  bool full_rank = false;
105  int num_poly = 10;
106  int total_poly = std::min( num_poly, num_sources );
107  for ( int j = 1; j < total_poly; ++j )
108  {
109  // Add the next column.
110  P.reshape( num_sources, poly_id + 1 );
111  for ( int i = 0; i < num_sources; ++i )
112  {
113  source_center_view = source_centers( DIM * source_lids[i], DIM );
114  P( i, poly_id ) = polynomialCoefficient( j, source_center_view );
115  }
116 
117  // Check for rank deficiency.
118  full_rank = isFullRank( P );
119 
120  // If we are full rank, add this coefficient.
121  if ( full_rank )
122  {
123  poly_ids.push_back( j );
124  ++poly_id;
125  }
126 
127  // If we are rank deficient, remove the last column.
128  else
129  {
130  P.reshape( num_sources, poly_id );
131  }
132  }
133 
134  // Make p.
135  poly_size = poly_ids.size();
136  target_poly.resize( poly_size );
137  for ( int i = 0; i < poly_size; ++i )
138  {
139  target_poly( i ) = polynomialCoefficient( poly_ids[i], target_center );
140  }
141 
142  // Construct b.
143  Teuchos::SerialDenseMatrix<int, double> b( poly_size, num_sources );
144  b.multiply( Teuchos::TRANS, Teuchos::NO_TRANS, 1.0, P, phi, 0.0 );
145 
146  // Construct the A matrix.
147  Teuchos::SerialDenseMatrix<int, double> A( poly_size, poly_size );
148  {
149  // Build A.
150  Teuchos::SerialDenseMatrix<int, double> work( num_sources, poly_size );
151  work.multiply( Teuchos::NO_TRANS, Teuchos::NO_TRANS, 1.0, phi, P, 0.0 );
152  A.multiply( Teuchos::TRANS, Teuchos::NO_TRANS, 1.0, P, work, 0.0 );
153  }
154 
155  // Apply the inverse of the A matrix to b.
156  Teuchos::LAPACK<int, double> lapack;
157  double A_rcond = std::numeric_limits<double>::epsilon();
158  Teuchos::Array<double> work( 4 * A.numCols() );
159  Teuchos::SerialDenseVector<int, double> s( poly_size );
160  int rank = 0;
161  int info = 0;
162 
163  // Estimate the reciprocal of the condition number.
164  Teuchos::Array<int> ipiv( std::min( A.numRows(), A.numCols() ) );
165  Teuchos::SerialDenseMatrix<int, double> LU_A( A );
166  lapack.GETRF( LU_A.numRows(), LU_A.numCols(), LU_A.values(), LU_A.numRows(),
167  ipiv.getRawPtr(), &info );
168  DTK_CHECK( 0 == info );
169 
170  Teuchos::Array<int> iwork( A.numCols() );
171  lapack.GECON( '1', LU_A.numCols(), LU_A.values(), LU_A.numRows(),
172  A.normOne(), &A_rcond, work.getRawPtr(), iwork.getRawPtr(),
173  &info );
174  DTK_CHECK( 0 == info );
175 
176  // Get the optimal work size.
177  lapack.GELSS( A.numRows(), A.numCols(), b.numCols(), A.values(),
178  A.numRows(), b.values(), b.numRows(), s.values(), A_rcond,
179  &rank, work.getRawPtr(), -1, &info );
180  DTK_CHECK( 0 == info );
181 
182  // Apply the inverse of A to b.
183  work.resize( work[0] );
184  lapack.GELSS( A.numRows(), A.numCols(), b.numCols(), A.values(),
185  A.numRows(), b.values(), b.numRows(), s.values(), A_rcond,
186  &rank, work.getRawPtr(), work.size(), &info );
187  DTK_CHECK( 0 == info );
188 
189  // Construct the basis.
190  Teuchos::SerialDenseMatrix<int, double> shape_matrix(
191  Teuchos::View, d_shape_function.getRawPtr(), 1, 1,
192  d_shape_function.size() );
193  shape_matrix.multiply( Teuchos::TRANS, Teuchos::NO_TRANS, 1.0, target_poly,
194  b, 0.0 );
195 }
196 
197 //---------------------------------------------------------------------------//
198 // Get a polynomial coefficient.
199 template <class Basis, int DIM>
201  const int coeff, const Teuchos::ArrayView<const double> &center ) const
202 {
203  switch ( coeff )
204  {
205  // Linear.
206  case ( 0 ):
207  return 1.0;
208  break;
209  case ( 1 ):
210  return center[0];
211  break;
212  case ( 2 ):
213  return center[1];
214  break;
215  case ( 3 ):
216  return center[2];
217  break;
218 
219  // Quadratic
220  case ( 4 ):
221  return center[0] * center[1];
222  break;
223  case ( 5 ):
224  return center[0] * center[2];
225  break;
226  case ( 6 ):
227  return center[1] * center[2];
228  break;
229  case ( 7 ):
230  return center[0] * center[0];
231  break;
232  case ( 8 ):
233  return center[1] * center[1];
234  break;
235  case ( 9 ):
236  return center[2] * center[2];
237  break;
238  }
239  return 0.0;
240 }
241 
242 //---------------------------------------------------------------------------//
243 // Check if a matrix is full rank.
244 template <class Basis, int DIM>
246  const Teuchos::SerialDenseMatrix<int, double> &matrix ) const
247 {
248  // Copy the matrix.
249  Teuchos::SerialDenseMatrix<int, double> A = matrix;
250 
251  // Determine the full rank.
252  int full_rank = std::min( A.numRows(), A.numCols() );
253 
254  // Compute the singular value decomposition.
255  Teuchos::LAPACK<int, double> lapack;
256  Teuchos::Array<double> S( full_rank );
257  Teuchos::SerialDenseMatrix<int, double> U( A.numRows(), A.numRows() );
258  Teuchos::SerialDenseMatrix<int, double> VT( A.numCols(), A.numCols() );
259  Teuchos::Array<double> work( full_rank );
260  Teuchos::Array<double> rwork( full_rank );
261  int info = 0;
262  lapack.GESVD( 'A', 'A', A.numRows(), A.numCols(), A.values(), A.numRows(),
263  S.getRawPtr(), U.values(), U.numRows(), VT.values(),
264  VT.numRows(), work.getRawPtr(), -1, rwork.getRawPtr(),
265  &info );
266  DTK_CHECK( 0 == info );
267 
268  work.resize( work[0] );
269  rwork.resize( work.size() );
270  lapack.GESVD( 'A', 'A', A.numRows(), A.numCols(), A.values(), A.numRows(),
271  S.getRawPtr(), U.values(), U.numRows(), VT.values(),
272  VT.numRows(), work.getRawPtr(), work.size(),
273  rwork.getRawPtr(), &info );
274  DTK_CHECK( 0 == info );
275 
276  // Check the singular values. If they are greater than delta they count.
277  double epsilon = std::numeric_limits<double>::epsilon();
278  double delta = S[0] * epsilon;
279  int rank = std::count_if( S.begin(), S.end(),
280  [=]( double s ) { return ( s > delta ); } );
281 
282  return ( rank == full_rank );
283 }
284 
285 //---------------------------------------------------------------------------//
286 
287 } // end namespace DataTransferKit
288 
289 //---------------------------------------------------------------------------//
290 
291 #endif // end DTK_LOCALMLSPROBLEM_IMPL_HPP
292 
293 //---------------------------------------------------------------------------//
294 // end DTK_LocalMLSProblem_impl.hpp
295 //---------------------------------------------------------------------------//
Local moving least square problem.
static double evaluateValue(const RadialBasis &basis, const double radius, const double x)
Compute the value of the basis at the given value.
Assertions and Design-by-Contract for error handling.
Policy class for spline interpolation basis functions.
Local moving least square problem about a single target center using quadratic polynomials.
Euclidean distance function.
DTK_BasicEntitySet.cpp.