Stokhos Package Browser (Single Doxygen Collection)  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Stokhos_StieltjesGramSchmidtBuilderImp.hpp
Go to the documentation of this file.
1 // @HEADER
2 // *****************************************************************************
3 // Stokhos Package
4 //
5 // Copyright 2009 NTESS and the Stokhos contributors.
6 // SPDX-License-Identifier: BSD-3-Clause
7 // *****************************************************************************
8 // @HEADER
9 
13 
14 template <typename ordinal_type, typename value_type>
19  ordinal_type new_order, bool use_pce_qp, bool normalize) :
20  quad(quad_)
21 {
22  // Create array to store new coordinate bases
23  ordinal_type new_dim = pces.size();
25 
26  // Create Stieltjes basis for each pce
27  for (ordinal_type k=0; k<new_dim; k++) {
28  new_coordinate_bases[k] = Teuchos::rcp(
30  new_order, Teuchos::rcp(&(pces[k]),false), quad, use_pce_qp,
31  normalize)
32  );
33  }
34 
35  // Create tensor product basis from coordinate bases
37  new CompletePolynomialBasis<ordinal_type,value_type>(new_coordinate_bases)
38  );
39 
40  // Use Gram-Schmidt to orthogonalize tensor product bases
41  const Teuchos::Array<value_type>& weights = quad->getQuadWeights();
43  quad->getQuadPoints();
44  ordinal_type nqp = points.size();
49  for (ordinal_type i=0; i<nqp; i++)
50  (*new_points)[i].resize(new_dim);
51  for (ordinal_type k=0; k<new_dim; k++) {
53  Teuchos::Array<value_type> st_weights;
55  new_coordinate_bases[k]->getQuadPoints(new_order+1, st_points, st_weights,
56  st_values);
57 
58  for (ordinal_type i=0; i<nqp; i++)
59  (*new_points)[i][k] = st_points[i];
60  }
63  *new_points,
64  *new_weights,
65  1e-15)
66  );
67 
68  // Create new quadrature object
70  gs_basis;
73  new_points,
74  new_weights)
75  );
76 
77 }
78 
79 template <typename ordinal_type, typename value_type>
83 {
84  return gs_basis;
85 }
86 
87 template <typename ordinal_type, typename value_type>
91 {
92  return gs_quad;
93 }
94 
95 template <typename ordinal_type, typename value_type>
96 void
101 {
102  // Map pce coefficients to tensor basis to Gram-Schmidt basis
103  ordinal_type dim = pces.size();
104  if (new_pces.size() != pces.size())
105  new_pces.resize(dim);
106  for (ordinal_type k=0; k<dim; k++) {
107  OrthogPolyApprox<ordinal_type,value_type> p_tensor(tensor_basis);
108  p_tensor.term(k, 0) = pces[k].mean();
109  p_tensor.term(k, 1) = 1.0;
110  new_pces[k].reset(gs_basis);
111  gs_basis->transformCoeffs(p_tensor.coeff(), new_pces[k].coeff());
112  }
113 }
Teuchos::RCP< const OrthogPolyBasis< ordinal_type, value_type > > tensor_basis
Reduced tensor basis.
Teuchos::RCP< const Quadrature< ordinal_type, value_type > > quad
Quadrature object for original basis.
Teuchos::RCP< UserDefinedQuadrature< ordinal_type, value_type > > gs_quad
Reduced quadrature.
pointer coeff()
Return coefficient array.
StieltjesGramSchmidtBuilder(const Teuchos::RCP< const Quadrature< ordinal_type, value_type > > &quad, const Teuchos::Array< OrthogPolyApprox< ordinal_type, value_type > > &pces, ordinal_type new_order, bool use_pce_qp, bool normalize)
Constructor.
Generates three-term recurrence using the Discretized Stieltjes procedure applied to a polynomial cha...
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
Abstract base class for quadrature methods.
Transforms a non-orthogonal multivariate basis to an orthogonal one using the Gram-Schmit procedure...
Teuchos::RCP< Quadrature< ordinal_type, value_type > > getReducedQuadrature() const
Get reduced quadrature.
Multivariate orthogonal polynomial basis generated from a total-order complete-polynomial tensor prod...
void computeReducedPCEs(const Teuchos::Array< OrthogPolyApprox< ordinal_type, value_type > > &pces, Teuchos::Array< OrthogPolyApprox< ordinal_type, value_type > > &new_pces)
Get reduced PCEs.
size_type size() const
Teuchos::RCP< const OrthogPolyBasis< ordinal_type, value_type > > getReducedBasis() const
Get reduced basis.
Teuchos::RCP< GramSchmidtBasis< ordinal_type, value_type > > gs_basis
Reduced Gram-Schmidt basis.
reference term(ordinal_type dimension, ordinal_type order)
Get coefficient term for given dimension and order.