44 template <
typename ordinal_type,
typename value_type>
52 pce_basis(pce[0].basis()),
53 pce_sz(pce_basis->size()),
56 verbose(params.get(
"Verbose", false)),
57 rank_threshold(params.get(
"Rank Threshold", 1.0e-12)),
58 orthogonalization_method(params.get(
"Orthogonalization Method",
63 template <
typename ordinal_type,
typename value_type>
75 if (pce[i].standard_deviation() > 1.0e-15)
83 quad->getQuadPoints();
85 quad->getBasisAtQuadPoints();
93 A(i,
j) = basis_values[i][
j];
99 pce_norms[
j] += (*pce2[
j])[i]*(*pce2[
j])[i]*pce_basis->norm_squared(i);
113 F(i,
j) = pce2[
j]->evaluate(points[i], basis_values[i]);
116 sz = buildReducedBasis(max_p, rank_threshold, A, F, weights, terms, num_terms,
124 if (quad_params.
isParameter(
"Reduced Quadrature Method") &&
125 quad_params.
get<std::string>(
"Reduced Quadrature Method") ==
"Q2") {
128 value_type rank_threshold2 = quad_params.
get(
"Q2 Rank Threshold",
132 buildReducedBasis(2*max_p, rank_threshold2, A, F, weights, terms2,
133 num_terms2, Qp2, Q2);
138 norms.resize(sz, 1.0);
141 template <
typename ordinal_type,
typename value_type>
147 template <
typename ordinal_type,
typename value_type>
155 template <
typename ordinal_type,
typename value_type>
163 template <
typename ordinal_type,
typename value_type>
171 template <
typename ordinal_type,
typename value_type>
179 template <
typename ordinal_type,
typename value_type>
187 template <
typename ordinal_type,
typename value_type>
196 template <
typename ordinal_type,
typename value_type>
205 template <
typename ordinal_type,
typename value_type>
213 template <
typename ordinal_type,
typename value_type>
222 template <
typename ordinal_type,
typename value_type>
227 os <<
"Gram-Schmidt basis of order " << p <<
", dimension " << d
228 <<
", and size " << sz <<
". Matrix coefficients:\n";
229 os << Qp << std::endl;
230 os <<
"Basis vector norms (squared):\n\t";
232 os << norms[i] <<
" ";
236 template <
typename ordinal_type,
typename value_type>
258 template <
typename ordinal_type,
typename value_type>
280 template <
typename ordinal_type,
typename value_type>
KOKKOS_INLINE_FUNCTION PCE< Storage > sqrt(const PCE< Storage > &a)
virtual Teuchos::RCP< const Stokhos::UserDefinedQuadrature< ordinal_type, value_type > > createReducedQuadrature(const Teuchos::SerialDenseMatrix< ordinal_type, value_type > &Q, const Teuchos::SerialDenseMatrix< ordinal_type, value_type > &Q2, const Teuchos::SerialDenseMatrix< ordinal_type, value_type > &F, const Teuchos::Array< value_type > &weights) const
Get reduced quadrature object.
virtual Teuchos::RCP< Stokhos::Sparse3Tensor< ordinal_type, value_type > > computeTripleProductTensor() const
Compute triple product tensor.
T & get(ParameterList &l, const std::string &name)
int multiply(ETransp transa, ETransp transb, ScalarType alpha, const SerialDenseMatrix< OrdinalType, ScalarType > &A, const SerialDenseMatrix< OrdinalType, ScalarType > &B, ScalarType beta)
virtual ~GSReducedPCEBasisBase()
Destructor.
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
RCP< ParameterList > sublist(const RCP< ParameterList > ¶mList, const std::string &name, bool mustAlreadyExist=false, const std::string &docString="")
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 Teuchos::Array< value_type > & norm_squared() const
Return array storing norm-squared of each basis polynomial.
bool isParameter(const std::string &name) const
void setup(ordinal_type p, const Teuchos::Array< Stokhos::OrthogPolyApprox< ordinal_type, value_type > > &pce, const Teuchos::RCP< const Stokhos::Quadrature< ordinal_type, value_type > > &quad)
virtual void print(std::ostream &os) const
Print basis to stream os.
Abstract base class for quadrature methods.
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.
GSReducedPCEBasisBase(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::ParameterList ¶ms=Teuchos::ParameterList())
Constructor.
virtual ordinal_type size() const
Return total size of basis.
virtual value_type evaluateZero(ordinal_type i) const
Evaluate basis polynomial i at zero.
Generate a basis from a given set of PCE expansions that is orthogonal with respect to the product me...
ordinal_type order() const
Return order of basis.
void push_back(const value_type &x)
ordinal_type dimension() const
Return dimension of basis.
#define TEUCHOS_ASSERT(assertion_test)
virtual Teuchos::RCP< Stokhos::Sparse3Tensor< ordinal_type, value_type > > computeLinearTripleProductTensor() const
Compute linear triple product tensor where k = 0,1,..,d.
virtual Teuchos::RCP< const Stokhos::Quadrature< ordinal_type, value_type > > getReducedQuadrature() const
Get reduced quadrature object.
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.