Stokhos Package Browser (Single Doxygen Collection) Version of the Day
Loading...
Searching...
No Matches
Stokhos_ProductLanczosPCEBasisImp.hpp
Go to the documentation of this file.
1// @HEADER
2// ***********************************************************************
3//
4// Stokhos Package
5// Copyright (2009) Sandia Corporation
6//
7// Under terms of Contract DE-AC04-94AL85000, there is a non-exclusive
8// license for use of this work by or on behalf of the U.S. Government.
9//
10// Redistribution and use in source and binary forms, with or without
11// modification, are permitted provided that the following conditions are
12// met:
13//
14// 1. Redistributions of source code must retain the above copyright
15// notice, this list of conditions and the following disclaimer.
16//
17// 2. Redistributions in binary form must reproduce the above copyright
18// notice, this list of conditions and the following disclaimer in the
19// documentation and/or other materials provided with the distribution.
20//
21// 3. Neither the name of the Corporation nor the names of the
22// contributors may be used to endorse or promote products derived from
23// this software without specific prior written permission.
24//
25// THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
26// EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
27// IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
28// PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
29// CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
30// EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
31// PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
32// PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
33// LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
34// NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
35// SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
36//
37// Questions? Contact Eric T. Phipps (etphipp@sandia.gov).
38//
39// ***********************************************************************
40// @HEADER
41
42#include "Stokhos_SDMUtils.hpp"
47
48template <typename ordinal_type, typename value_type>
51 ordinal_type p,
52 const Teuchos::Array< Stokhos::OrthogPolyApprox<ordinal_type, value_type> >& pce,
53 const Teuchos::RCP<const Stokhos::Quadrature<ordinal_type, value_type> >& quad,
54 const Teuchos::RCP< const Stokhos::Sparse3Tensor<ordinal_type, value_type> >& Cijk,
55 const Teuchos::ParameterList& params_) :
56 name("Product Lanczos PCE Basis"),
57 params(params_)
58{
59 Teuchos::RCP<const Stokhos::OrthogPolyBasis<ordinal_type,value_type> > pce_basis = pce[0].basis();
60 ordinal_type pce_sz = pce_basis->size();
61
62 // Check if basis is a product basis
63 Teuchos::RCP<const Stokhos::ProductBasis<ordinal_type,value_type> > prod_basis = Teuchos::rcp_dynamic_cast<const Stokhos::ProductBasis<ordinal_type,value_type> >(pce_basis);
64 Teuchos::Array< Teuchos::RCP<const OneDOrthogPolyBasis<ordinal_type,value_type> > > coord_bases;
65 if (prod_basis != Teuchos::null)
66 coord_bases = prod_basis->getCoordinateBases();
67
68 // Build Lanczos basis for each pce
69 bool project = params.get("Project", true);
70 bool normalize = params.get("Normalize", true);
71 bool limit_integration_order = params.get("Limit Integration Order", false);
72 bool use_stieltjes = params.get("Use Old Stieltjes Method", false);
73 Teuchos::Array< Teuchos::RCP<const Stokhos::OneDOrthogPolyBasis<int,double > > > coordinate_bases;
74 Teuchos::Array<int> is_invariant(pce.size(),-2);
75 for (ordinal_type i=0; i<pce.size(); i++) {
76
77 // Check for pce's lying in invariant subspaces, which are pce's that
78 // depend on only a single dimension. In this case use the corresponding
79 // original coordinate basis. Convention is: -2 -- not invariant, -1 --
80 // constant, i >= 0 pce depends only on dimension i
81 if (prod_basis != Teuchos::null)
82 is_invariant[i] = isInvariant(pce[i]);
83 if (is_invariant[i] >= 0) {
84 coordinate_bases.push_back(coord_bases[is_invariant[i]]);
85 }
86
87 // Exclude constant pce's from the basis since they don't represent
88 // stochastic dimensions
89 else if (is_invariant[i] != -1) {
90 if (use_stieltjes) {
91 coordinate_bases.push_back(
92 Teuchos::rcp(
94 p, Teuchos::rcp(&(pce[i]),false), quad, false,
95 normalize, project, Cijk)));
96 }
97 else {
98 if (project)
99 coordinate_bases.push_back(
100 Teuchos::rcp(
102 p, Teuchos::rcp(&(pce[i]),false), Cijk,
103 normalize, limit_integration_order)));
104 else
105 coordinate_bases.push_back(
106 Teuchos::rcp(
108 p, Teuchos::rcp(&(pce[i]),false), quad,
109 normalize, limit_integration_order)));
110 }
111 }
112 }
113 ordinal_type d = coordinate_bases.size();
114
115 // Build tensor product basis
117 Teuchos::rcp(
119 coordinate_bases,
120 params.get("Cijk Drop Tolerance", 1.0e-15),
121 params.get("Use Old Cijk Algorithm", false)));
122
123 // Build reduced quadrature
124 Teuchos::ParameterList sg_params;
125 sg_params.sublist("Basis").set< Teuchos::RCP< const Stokhos::OrthogPolyBasis<ordinal_type,value_type> > >("Stochastic Galerkin Basis", tensor_lanczos_basis);
126 sg_params.sublist("Quadrature") = params.sublist("Reduced Quadrature");
127 reduced_quad =
129
130 // Build Psi matrix -- Psi_ij = Psi_i(x^j)*w_j/<Psi_i^2>
131 const Teuchos::Array<value_type>& weights = quad->getQuadWeights();
132 const Teuchos::Array< Teuchos::Array<value_type> >& points =
133 quad->getQuadPoints();
134 const Teuchos::Array< Teuchos::Array<value_type> >& basis_vals =
135 quad->getBasisAtQuadPoints();
136 ordinal_type nqp = weights.size();
137 SDM Psi(pce_sz, nqp);
138 for (ordinal_type i=0; i<pce_sz; i++)
139 for (ordinal_type k=0; k<nqp; k++)
140 Psi(i,k) = basis_vals[k][i]*weights[k]/pce_basis->norm_squared(i);
141
142 // Build Phi matrix -- Phi_ij = Phi_i(y(x^j))
143 ordinal_type sz = tensor_lanczos_basis->size();
144 Teuchos::Array<value_type> red_basis_vals(sz);
145 Teuchos::Array<value_type> pce_vals(d);
146 Phi.shape(sz, nqp);
147 for (int k=0; k<nqp; k++) {
148 ordinal_type jdx = 0;
149 for (int j=0; j<pce.size(); j++) {
150
151 // Exclude constant pce's
152 if (is_invariant[j] != -1) {
153
154 // Use the identity mapping for invariant subspaces
155 if (is_invariant[j] >= 0)
156 pce_vals[jdx] = points[k][is_invariant[j]];
157 else
158 pce_vals[jdx] = pce[j].evaluate(points[k], basis_vals[k]);
159 jdx++;
160
161 }
162
163 }
164 tensor_lanczos_basis->evaluateBases(pce_vals, red_basis_vals);
165 for (int i=0; i<sz; i++)
166 Phi(i,k) = red_basis_vals[i];
167 }
168
169 bool verbose = params.get("Verbose", false);
170
171 // Compute matrix A mapping reduced space to original
172 A.shape(pce_sz, sz);
173 ordinal_type ret =
174 A.multiply(Teuchos::NO_TRANS, Teuchos::TRANS, 1.0, Psi, Phi, 0.0);
175 TEUCHOS_ASSERT(ret == 0);
176 //print_matlab(std::cout << "A = " << std::endl, A);
177
178 // Compute pseudo-inverse of A mapping original space to reduced
179 // A = U*S*V^T -> A^+ = V*S^+*U^T = (S^+*V^T)^T*U^T where
180 // S^+ is a diagonal matrix comprised of the inverse of the diagonal of S
181 // for each nonzero, and zero otherwise
182 Teuchos::Array<value_type> sigma;
183 SDM U, Vt;
184 value_type rank_threshold = params.get("Rank Threshold", 1.0e-12);
185 ordinal_type rank = svd_threshold(rank_threshold, A, sigma, U, Vt);
186 Ainv.shape(sz, pce_sz);
187 TEUCHOS_ASSERT(rank == Vt.numRows());
188 for (ordinal_type i=0; i<Vt.numRows(); i++)
189 for (ordinal_type j=0; j<Vt.numCols(); j++)
190 Vt(i,j) = Vt(i,j) / sigma[i];
191 ret = Ainv.multiply(Teuchos::TRANS, Teuchos::TRANS, 1.0, Vt, U, 0.0);
192 TEUCHOS_ASSERT(ret == 0);
193 //print_matlab(std::cout << "Ainv = " << std::endl, Ainv);
194
195 if (verbose) {
196 std::cout << "rank = " << rank << std::endl;
197
198 std::cout << "diag(S) = [";
199 for (ordinal_type i=0; i<rank; i++)
200 std::cout << sigma[i] << " ";
201 std::cout << "]" << std::endl;
202
203 // Check A = U*S*V^T
204 SDM SVt(rank, Vt.numCols());
205 for (ordinal_type i=0; i<Vt.numRows(); i++)
206 for (ordinal_type j=0; j<Vt.numCols(); j++)
207 SVt(i,j) = Vt(i,j) * sigma[i] * sigma[i]; // since we divide by sigma
208 // above
209 SDM err_A(pce_sz,sz);
210 err_A.assign(A);
211 ret = err_A.multiply(Teuchos::NO_TRANS, Teuchos::NO_TRANS, -1.0, U, SVt,
212 1.0);
213 TEUCHOS_ASSERT(ret == 0);
214 std::cout << "||A - U*S*V^T||_infty = " << err_A.normInf() << std::endl;
215 //print_matlab(std::cout << "A - U*S*V^T = " << std::endl, err_A);
216
217 // Check Ainv*A == I
218 SDM err(sz,sz);
219 err.putScalar(0.0);
220 for (ordinal_type i=0; i<sz; i++)
221 err(i,i) = 1.0;
222 ret = err.multiply(Teuchos::NO_TRANS, Teuchos::NO_TRANS, 1.0, Ainv, A,
223 -1.0);
224 TEUCHOS_ASSERT(ret == 0);
225 std::cout << "||Ainv*A - I||_infty = " << err.normInf() << std::endl;
226 //print_matlab(std::cout << "Ainv*A-I = " << std::endl, err);
227
228 // Check A*Ainv == I
229 SDM err2(pce_sz,pce_sz);
230 err2.putScalar(0.0);
231 for (ordinal_type i=0; i<pce_sz; i++)
232 err2(i,i) = 1.0;
233 ret = err2.multiply(Teuchos::NO_TRANS, Teuchos::NO_TRANS, 1.0, A, Ainv, -1.0);
234 TEUCHOS_ASSERT(ret == 0);
235 std::cout << "||A*Ainv - I||_infty = " << err2.normInf() << std::endl;
236 //print_matlab(std::cout << "A*Ainv-I = " << std::endl, err2);
237 }
238}
239
240template <typename ordinal_type, typename value_type>
245
246template <typename ordinal_type, typename value_type>
247ordinal_type
249order() const
250{
251 return tensor_lanczos_basis->order();
252}
253
254template <typename ordinal_type, typename value_type>
255ordinal_type
257dimension() const
258{
259 return tensor_lanczos_basis->dimension();
260}
261
262template <typename ordinal_type, typename value_type>
263ordinal_type
265size() const
266{
267 return tensor_lanczos_basis->size();
268}
269
270template <typename ordinal_type, typename value_type>
271const Teuchos::Array<value_type>&
273norm_squared() const
274{
275 return tensor_lanczos_basis->norm_squared();
276}
277
278template <typename ordinal_type, typename value_type>
279const value_type&
281norm_squared(ordinal_type i) const
282{
283 return tensor_lanczos_basis->norm_squared(i);
284}
285
286template <typename ordinal_type, typename value_type>
287Teuchos::RCP< Stokhos::Sparse3Tensor<ordinal_type, value_type> >
290
291{
292 Teuchos::RCP< Stokhos::Sparse3Tensor<ordinal_type, value_type> > Cijk =
293 tensor_lanczos_basis->computeTripleProductTensor();
294 //std::cout << *Cijk << std::endl;
295 return Cijk;
296}
297
298template <typename ordinal_type, typename value_type>
299Teuchos::RCP< Stokhos::Sparse3Tensor<ordinal_type, value_type> >
302
303{
304 Teuchos::RCP< Stokhos::Sparse3Tensor<ordinal_type, value_type> > Cijk =
305 tensor_lanczos_basis->computeLinearTripleProductTensor();
306 //std::cout << *Cijk << std::endl;
307 return Cijk;
308}
309
310template <typename ordinal_type, typename value_type>
311value_type
313evaluateZero(ordinal_type i) const
314{
315 return tensor_lanczos_basis->evaluateZero(i);
316}
317
318template <typename ordinal_type, typename value_type>
319void
321evaluateBases(const Teuchos::ArrayView<const value_type>& point,
322 Teuchos::Array<value_type>& basis_vals) const
323{
324 return tensor_lanczos_basis->evaluateBases(point, basis_vals);
325}
326
327template <typename ordinal_type, typename value_type>
328void
330print(std::ostream& os) const
331{
332 tensor_lanczos_basis->print(os);
333}
334
335template <typename ordinal_type, typename value_type>
336const std::string&
342
343template <typename ordinal_type, typename value_type>
346term(ordinal_type i) const
347{
348 return tensor_lanczos_basis->term(i);
349}
350
351template <typename ordinal_type, typename value_type>
352ordinal_type
355{
356 return tensor_lanczos_basis->index(term);
357}
358
359template <typename ordinal_type, typename value_type>
360Teuchos::Array< Teuchos::RCP<const Stokhos::OneDOrthogPolyBasis<ordinal_type, value_type> > >
362getCoordinateBases() const
363{
364 return tensor_lanczos_basis->getCoordinateBases();
365}
366
367template <typename ordinal_type, typename value_type>
370getMaxOrders() const
371{
372 return tensor_lanczos_basis->getMaxOrders();
373}
374
375template <typename ordinal_type, typename value_type>
376void
378transformToOriginalBasis(const value_type *in, value_type *out,
379 ordinal_type ncol, bool transpose) const
380{
381 ordinal_type pce_sz = A.numRows();
382 ordinal_type sz = A.numCols();
383 if (transpose) {
384 SDM zbar(Teuchos::View, const_cast<value_type*>(in), ncol, ncol, sz);
385 SDM z(Teuchos::View, out, ncol, ncol, pce_sz);
386 ordinal_type ret =
387 z.multiply(Teuchos::NO_TRANS, Teuchos::TRANS, 1.0, zbar, A, 0.0);
388 TEUCHOS_ASSERT(ret == 0);
389 }
390 else {
391 SDM zbar(Teuchos::View, const_cast<value_type*>(in), sz, sz, ncol);
392 SDM z(Teuchos::View, out, pce_sz, pce_sz, ncol);
393 ordinal_type ret =
394 z.multiply(Teuchos::NO_TRANS, Teuchos::NO_TRANS, 1.0, A, zbar, 0.0);
395 TEUCHOS_ASSERT(ret == 0);
396 }
397}
398
399template <typename ordinal_type, typename value_type>
400void
402transformFromOriginalBasis(const value_type *in, value_type *out,
403 ordinal_type ncol, bool transpose) const
404{
405 ordinal_type pce_sz = A.numRows();
406 ordinal_type sz = A.numCols();
407 if (transpose) {
408 SDM z(Teuchos::View, const_cast<value_type*>(in), ncol, ncol, pce_sz);
409 SDM zbar(Teuchos::View, out, ncol, ncol, sz);
410 ordinal_type ret =
411 zbar.multiply(Teuchos::NO_TRANS, Teuchos::TRANS, 1.0, z, Ainv, 0.0);
412 TEUCHOS_ASSERT(ret == 0);
413 }
414 else {
415 SDM z(Teuchos::View, const_cast<value_type*>(in), pce_sz, pce_sz, ncol);
416 SDM zbar(Teuchos::View, out, sz, sz, ncol);
417 ordinal_type ret =
418 zbar.multiply(Teuchos::NO_TRANS, Teuchos::NO_TRANS, 1.0, Ainv, z, 0.0);
419 TEUCHOS_ASSERT(ret == 0);
420 }
421}
422
423template <typename ordinal_type, typename value_type>
424Teuchos::RCP<const Stokhos::Quadrature<ordinal_type, value_type> >
430
431template <typename ordinal_type, typename value_type>
432ordinal_type
435{
436 Teuchos::RCP<const Stokhos::OrthogPolyBasis<ordinal_type,value_type> > basis =
437 pce.basis();
438 ordinal_type dim = basis->dimension();
439 value_type tol = 1.0e-15;
440
441 // Check if basis is a product basis
442 Teuchos::RCP<const Stokhos::ProductBasis<ordinal_type,value_type> > prod_basis = Teuchos::rcp_dynamic_cast<const Stokhos::ProductBasis<ordinal_type,value_type> >(basis);
443 if (prod_basis == Teuchos::null)
444 return -2;
445
446 // Build list of dimensions pce depends on by looping over each dimension,
447 // computing norm of pce with just that dimension -- note we don't include
448 // the constant term
449 Teuchos::Array<ordinal_type> dependent_dims;
450 tmp_pce.reset(basis);
451 for (ordinal_type i=0; i<dim; i++) {
452 ordinal_type p = prod_basis->getCoordinateBases()[i]->order();
453 tmp_pce.init(0.0);
454 for (ordinal_type j=1; j<=p; j++)
455 tmp_pce.term(i,j) = pce.term(i,j);
456 value_type nrm = tmp_pce.two_norm();
457 if (nrm > tol) dependent_dims.push_back(i);
458 }
459
460 // If dependent_dims has length 1, pce a function of a single variable,
461 // which is an invariant subspace
462 if (dependent_dims.size() == 1)
463 return dependent_dims[0];
464
465 // If dependent_dims has length 0, pce is constant
466 else if (dependent_dims.size() == 0)
467 return -1;
468
469 // Otherwise pce depends on more than one variable
470 return -2;
471}
Multivariate orthogonal polynomial basis generated from a total-order complete-polynomial tensor prod...
Generates three-term recurrence using the Lanczos procedure applied to a polynomial chaos expansion i...
Generates three-term recurrence using the Lanczos procedure applied to a polynomial chaos expansion i...
A multidimensional index.
Class to store coefficients of a projection onto an orthogonal polynomial basis.
reference term(ordinal_type dimension, ordinal_type order)
Get coefficient term for given dimension and order.
Teuchos::RCP< const Stokhos::OrthogPolyBasis< ordinal_type, value_type > > basis() const
Return basis.
virtual value_type evaluateZero(ordinal_type i) const
Evaluate basis polynomial i at zero.
virtual ordinal_type index(const MultiIndex< ordinal_type > &term) const
Get index of the multivariate polynomial given orders of each coordinate.
virtual Teuchos::RCP< Stokhos::Sparse3Tensor< ordinal_type, value_type > > computeLinearTripleProductTensor() const
Compute linear triple product tensor where k = 0,1,..,d.
Teuchos::Array< Teuchos::RCP< const OneDOrthogPolyBasis< ordinal_type, value_type > > > getCoordinateBases() const
Return coordinate bases.
virtual void print(std::ostream &os) const
Print basis to stream os.
virtual const Teuchos::Array< value_type > & norm_squared() const
Return array storing norm-squared of each basis polynomial.
Teuchos::ParameterList params
Algorithm parameters.
SDM Phi
Values of transformed basis at quadrature points.
virtual const std::string & getName() const
Return string name of basis.
SDM A
Transition matrix from reduced basis to original.
virtual Teuchos::RCP< Stokhos::Sparse3Tensor< ordinal_type, value_type > > computeTripleProductTensor() const
Compute triple product tensor.
virtual Teuchos::RCP< const Stokhos::Quadrature< ordinal_type, value_type > > getReducedQuadrature() const
Get reduced quadrature object.
ordinal_type isInvariant(const Stokhos::OrthogPolyApprox< ordinal_type, value_type > &pce) const
Teuchos::RCP< Stokhos::CompletePolynomialBasis< ordinal_type, value_type > > tensor_lanczos_basis
Product Lanczos basis.
Teuchos::RCP< const Stokhos::Quadrature< ordinal_type, value_type > > reduced_quad
Reduced quadrature object.
Teuchos::SerialDenseMatrix< ordinal_type, value_type > SDM
ordinal_type order() const
Return order of basis.
virtual void transformToOriginalBasis(const value_type *in, value_type *out, ordinal_type ncol=1, bool transpose=false) const
Transform coefficients to original basis from this basis.
virtual void evaluateBases(const Teuchos::ArrayView< const value_type > &point, Teuchos::Array< value_type > &basis_vals) const
Evaluate basis polynomials at given point point.
virtual const MultiIndex< ordinal_type > & term(ordinal_type i) const
Get orders of each coordinate polynomial given an index i.
virtual void transformFromOriginalBasis(const value_type *in, value_type *out, ordinal_type ncol=1, bool transpose=false) const
Transform coefficients from original basis to this basis.
ProductLanczosPCEBasis(ordinal_type p, const Teuchos::Array< Stokhos::OrthogPolyApprox< ordinal_type, value_type > > &pce, const Teuchos::RCP< const Stokhos::Quadrature< ordinal_type, value_type > > &quad, const Teuchos::RCP< const Stokhos::Sparse3Tensor< ordinal_type, value_type > > &Cijk, const Teuchos::ParameterList &params=Teuchos::ParameterList())
Constructor.
ordinal_type dimension() const
Return dimension of basis.
SDM Ainv
Projection matrix from original matrix to reduced.
virtual ordinal_type size() const
Return total size of basis.
virtual MultiIndex< ordinal_type > getMaxOrders() const
Return maximum order allowable for each coordinate basis.
static Teuchos::RCP< const Stokhos::Quadrature< ordinal_type, value_type > > create(Teuchos::ParameterList &sgParams)
Generate quadrature object.
Abstract base class for quadrature methods.
Data structure storing a sparse 3-tensor C(i,j,k) in a a compressed format.
Generates three-term recurrence using the Discretized Stieltjes procedure applied to a polynomial cha...
ordinal_type svd_threshold(const scalar_type &rank_threshold, const Teuchos::SerialDenseMatrix< ordinal_type, scalar_type > &A, Teuchos::Array< scalar_type > &s, Teuchos::SerialDenseMatrix< ordinal_type, scalar_type > &U, Teuchos::SerialDenseMatrix< ordinal_type, scalar_type > &Vt)