This file is indexed.

/usr/include/trilinos/Stokhos_RecurrenceBasisImp.hpp is in libtrilinos-stokhos-dev 12.12.1-5.

This file is owned by root:root, with mode 0o644.

The actual contents of the file can be viewed below.

  1
  2
  3
  4
  5
  6
  7
  8
  9
 10
 11
 12
 13
 14
 15
 16
 17
 18
 19
 20
 21
 22
 23
 24
 25
 26
 27
 28
 29
 30
 31
 32
 33
 34
 35
 36
 37
 38
 39
 40
 41
 42
 43
 44
 45
 46
 47
 48
 49
 50
 51
 52
 53
 54
 55
 56
 57
 58
 59
 60
 61
 62
 63
 64
 65
 66
 67
 68
 69
 70
 71
 72
 73
 74
 75
 76
 77
 78
 79
 80
 81
 82
 83
 84
 85
 86
 87
 88
 89
 90
 91
 92
 93
 94
 95
 96
 97
 98
 99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
// $Id$
// $Source$
// @HEADER
// ***********************************************************************
//
//                           Stokhos Package
//                 Copyright (2009) Sandia Corporation
//
// Under terms of Contract DE-AC04-94AL85000, there is a non-exclusive
// license for use of this work by or on behalf of the U.S. Government.
//
// Redistribution and use in source and binary forms, with or without
// modification, are permitted provided that the following conditions are
// met:
//
// 1. Redistributions of source code must retain the above copyright
// notice, this list of conditions and the following disclaimer.
//
// 2. Redistributions in binary form must reproduce the above copyright
// notice, this list of conditions and the following disclaimer in the
// documentation and/or other materials provided with the distribution.
//
// 3. Neither the name of the Corporation nor the names of the
// contributors may be used to endorse or promote products derived from
// this software without specific prior written permission.
//
// THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
// EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
// IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
// PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
// CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
// EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
// PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
// PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
// LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
// NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
// SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
//
// Questions? Contact Eric T. Phipps (etphipp@sandia.gov).
//
// ***********************************************************************
// @HEADER

#include "Teuchos_LAPACK.hpp"
#include "Teuchos_SerialDenseMatrix.hpp"

template <typename ordinal_type, typename value_type>
Stokhos::RecurrenceBasis<ordinal_type, value_type>::
RecurrenceBasis(const std::string& name_, ordinal_type p_, bool normalize_,
                Stokhos::GrowthPolicy growth_) :
  name(name_),
  p(p_),
  normalize(normalize_),
  growth(growth_),
  quad_zero_tol(1.0e-14),
#ifdef HAVE_STOKHOS_DAKOTA
  sparse_grid_growth_rule(webbur::level_to_order_linear_nn),
#else
  sparse_grid_growth_rule(NULL),
#endif
  alpha(p+1, value_type(0.0)),
  beta(p+1, value_type(0.0)),
  delta(p+1, value_type(0.0)),
  gamma(p+1, value_type(0.0)),
  norms(p+1, value_type(0.0))
{
}

template <typename ordinal_type, typename value_type>
Stokhos::RecurrenceBasis<ordinal_type, value_type>::
RecurrenceBasis(ordinal_type p_, const RecurrenceBasis& basis) :
  name(basis.name),
  p(p_),
  normalize(basis.normalize),
  growth(basis.growth),
  quad_zero_tol(basis.quad_zero_tol),
  sparse_grid_growth_rule(basis.sparse_grid_growth_rule),
  alpha(p+1, value_type(0.0)),
  beta(p+1, value_type(0.0)),
  delta(p+1, value_type(0.0)),
  gamma(p+1, value_type(0.0)),
  norms(p+1, value_type(0.0))
{
}

template <typename ordinal_type, typename value_type>
void
Stokhos::RecurrenceBasis<ordinal_type, value_type>::
setup()
{
  bool is_normalized =
    computeRecurrenceCoefficients(p+1, alpha, beta, delta, gamma);
  if (normalize && !is_normalized) {
    normalizeRecurrenceCoefficients(alpha, beta, delta, gamma);
  }


  // Compute norms -- when polynomials are normalized, this should work
  // out to one (norms[0] == 1, delta[k] == 1, beta[k] == gamma[k]
  norms[0] = beta[0]/(gamma[0]*gamma[0]);
  for (ordinal_type k=1; k<=p; k++) {
    norms[k] = (beta[k]/gamma[k])*(delta[k-1]/delta[k])*norms[k-1];
  }
}

template <typename ordinal_type, typename value_type>
Stokhos::RecurrenceBasis<ordinal_type, value_type>::
~RecurrenceBasis()
{
}

template <typename ordinal_type, typename value_type>
ordinal_type
Stokhos::RecurrenceBasis<ordinal_type, value_type>::
order() const
{
  return p;
}

template <typename ordinal_type, typename value_type>
ordinal_type
Stokhos::RecurrenceBasis<ordinal_type, value_type>::
size() const
{
  return p+1;
}

template <typename ordinal_type, typename value_type>
const Teuchos::Array<value_type>&
Stokhos::RecurrenceBasis<ordinal_type, value_type>::
norm_squared() const
{
  return norms;
}

template <typename ordinal_type, typename value_type>
const value_type&
Stokhos::RecurrenceBasis<ordinal_type, value_type>::
norm_squared(ordinal_type i) const
{
  return norms[i];
}

template <typename ordinal_type, typename value_type>
Teuchos::RCP< Stokhos::Dense3Tensor<ordinal_type, value_type> >
Stokhos::RecurrenceBasis<ordinal_type, value_type>::
computeTripleProductTensor() const
{
  // Compute Cijk = < \Psi_i \Psi_j \Psi_k >
  ordinal_type sz = size();
  Teuchos::RCP< Stokhos::Dense3Tensor<ordinal_type, value_type> > Cijk =
    Teuchos::rcp(new Dense3Tensor<ordinal_type, value_type>(sz));
  Teuchos::Array<value_type> points, weights;
  Teuchos::Array< Teuchos::Array<value_type> > values;
  getQuadPoints(3*p, points, weights, values);

  for (ordinal_type i=0; i<sz; i++) {
    for (ordinal_type j=0; j<sz; j++) {
      for (ordinal_type k=0; k<sz; k++) {
        value_type triple_product = 0;
        for (ordinal_type l=0; l<static_cast<ordinal_type>(points.size());
             l++){
          triple_product +=
            weights[l]*(values[l][i])*(values[l][j])*(values[l][k]);
        }
        (*Cijk)(i,j,k) = triple_product;
      }
    }
  }

  return Cijk;
}

template <typename ordinal_type, typename value_type>
Teuchos::RCP< Stokhos::Sparse3Tensor<ordinal_type, value_type> >
Stokhos::RecurrenceBasis<ordinal_type, value_type>::
computeSparseTripleProductTensor(ordinal_type theOrder) const
{
  // Compute Cijk = < \Psi_i \Psi_j \Psi_k >
  value_type sparse_tol = 1.0e-15;
  ordinal_type sz = size();
  Teuchos::RCP< Stokhos::Sparse3Tensor<ordinal_type, value_type> > Cijk =
    Teuchos::rcp(new Sparse3Tensor<ordinal_type, value_type>());
  Teuchos::Array<value_type> points, weights;
  Teuchos::Array< Teuchos::Array<value_type> > values;
  getQuadPoints(3*p, points, weights, values);

  for (ordinal_type i=0; i<sz; i++) {
    for (ordinal_type j=0; j<sz; j++) {
      for (ordinal_type k=0; k<theOrder; k++) {
        value_type triple_product = 0;
        for (ordinal_type l=0; l<static_cast<ordinal_type>(points.size());
             l++){
          triple_product +=
            weights[l]*(values[l][i])*(values[l][j])*(values[l][k]);
        }
        if (std::abs(triple_product/norms[i]) > sparse_tol)
          Cijk->add_term(i,j,k,triple_product);
      }
    }
  }
  Cijk->fillComplete();

  return Cijk;
}

template <typename ordinal_type, typename value_type>
Teuchos::RCP< Teuchos::SerialDenseMatrix<ordinal_type, value_type> >
Stokhos::RecurrenceBasis<ordinal_type, value_type>::
computeDerivDoubleProductTensor() const
{
  // Compute Bij = < \Psi_i' \Psi_j >
  Teuchos::Array<value_type> points, weights;
  Teuchos::Array< Teuchos::Array<value_type> > values, derivs;
  getQuadPoints(2*p, points, weights, values);
  ordinal_type nqp = weights.size();
  derivs.resize(nqp);
  ordinal_type sz = size();
  for (ordinal_type i=0; i<nqp; i++) {
    derivs[i].resize(sz);
    evaluateBasesAndDerivatives(points[i], values[i], derivs[i]);
  }
  Teuchos::RCP< Teuchos::SerialDenseMatrix<ordinal_type, value_type> > Bij =
    Teuchos::rcp(new Teuchos::SerialDenseMatrix<ordinal_type, value_type>(sz,sz));
  for (ordinal_type i=0; i<sz; i++) {
    for (ordinal_type j=0; j<sz; j++) {
      value_type b = value_type(0.0);
      for (int qp=0; qp<nqp; qp++)
        b += weights[qp]*derivs[qp][i]*values[qp][j];
      (*Bij)(i,j) = b;
    }
  }

  return Bij;
}

template <typename ordinal_type, typename value_type>
void
Stokhos::RecurrenceBasis<ordinal_type,value_type>::
evaluateBases(const value_type& x, Teuchos::Array<value_type>& basis_pts) const
{
  // Evaluate basis polynomials P(x) using 3 term recurrence
  // P_0(x) = 1/gamma[0], P_-1 = 0
  // P_i(x) = [(delta[i-1]*x-alpha[i-1])*P_{i-1}(x) -
  //           beta[i-1]*P_{i-2}(x)]/gamma[i],
  // i=2,3,...
  basis_pts[0] = value_type(1)/gamma[0];
  if (p >= 1)
    basis_pts[1] = (delta[0]*x-alpha[0])*basis_pts[0]/gamma[1];
  for (ordinal_type i=2; i<=p; i++)
    basis_pts[i] = ((delta[i-1]*x-alpha[i-1])*basis_pts[i-1] -
                    beta[i-1]*basis_pts[i-2])/gamma[i];
}

template <typename ordinal_type, typename value_type>
void
Stokhos::RecurrenceBasis<ordinal_type,value_type>::
evaluateBasesAndDerivatives(const value_type& x,
                            Teuchos::Array<value_type>& vals,
                            Teuchos::Array<value_type>& derivs) const
{
  evaluateBases(x, vals);
  derivs[0] = 0.0;
  if (p >= 1)
    derivs[1] = delta[0]/(gamma[0]*gamma[1]);
  for (ordinal_type i=2; i<=p; i++)
    derivs[i] = (delta[i-1]*vals[i-1] + (delta[i-1]*x-alpha[i-1])*derivs[i-1] -
                 beta[i-1]*derivs[i-2])/gamma[i];
}

template <typename ordinal_type, typename value_type>
value_type
Stokhos::RecurrenceBasis<ordinal_type,value_type>::
evaluate(const value_type& x, ordinal_type k) const
{
  // Evaluate basis polynomials P(x) using 3 term recurrence
  // P_0(x) = 1/gamma[0], P_-1 = 0
  // P_i(x) = [(delta[i-1]*x-alpha[i-1])*P_{i-1}(x) -
  //           beta[i-1]*P_{i-2}(x)]/gamma[i],
  // i=2,3,...

  value_type v0 = value_type(1.0)/gamma[0];
  if (k == 0)
    return v0;

  value_type v1 = (x*delta[0]-alpha[0])*v0/gamma[1];
  if (k == 1)
    return v1;

  value_type v2 = value_type(0.0);
  for (ordinal_type i=2; i<=k; i++) {
    v2 = ((delta[i-1]*x-alpha[i-1])*v1 - beta[i-1]*v0)/gamma[i];
    v0 = v1;
    v1 = v2;
  }

  return v2;
}

template <typename ordinal_type, typename value_type>
void
Stokhos::RecurrenceBasis<ordinal_type, value_type>::
print(std::ostream& os) const
{
  os << name << " basis of order " << p << "." << std::endl;

  os << "Alpha recurrence coefficients:\n\t";
  for (ordinal_type i=0; i<=p; i++)
    os << alpha[i] << " ";
  os << std::endl;

  os << "Beta recurrence coefficients:\n\t";
  for (ordinal_type i=0; i<=p; i++)
    os << beta[i] << " ";
  os << std::endl;

  os << "Delta recurrence coefficients:\n\t";
  for (ordinal_type i=0; i<=p; i++)
    os << delta[i] << " ";
  os << std::endl;

  os << "Gamma recurrence coefficients:\n\t";
  for (ordinal_type i=0; i<=p; i++)
    os << gamma[i] << " ";
  os << std::endl;

  os << "Basis polynomial norms (squared):\n\t";
  for (ordinal_type i=0; i<=p; i++)
    os << norms[i] << " ";
  os << std::endl;
}

template <typename ordinal_type, typename value_type>
const std::string&
Stokhos::RecurrenceBasis<ordinal_type, value_type>::
getName() const
{
  return name;
}

template <typename ordinal_type, typename value_type>
void
Stokhos::RecurrenceBasis<ordinal_type,value_type>::
getQuadPoints(ordinal_type quad_order,
              Teuchos::Array<value_type>& quad_points,
              Teuchos::Array<value_type>& quad_weights,
              Teuchos::Array< Teuchos::Array<value_type> >& quad_values) const
{

  //This is a transposition into C++ of Gautschi's code for taking the first
  // N recurrance coefficient and generating a N point quadrature rule.
  // The MATLAB version is available at
  // http://www.cs.purdue.edu/archives/2002/wxg/codes/gauss.m
  ordinal_type num_points =
    static_cast<ordinal_type>(std::ceil((quad_order+1)/2.0));
  Teuchos::Array<value_type> a(num_points,0);
  Teuchos::Array<value_type> b(num_points,0);
  Teuchos::Array<value_type> c(num_points,0);
  Teuchos::Array<value_type> d(num_points,0);

  // If we don't have enough recurrance coefficients, get some more.
  if(num_points > p+1){
    bool is_normalized = computeRecurrenceCoefficients(num_points, a, b, c, d);
    if (!is_normalized)
      normalizeRecurrenceCoefficients(a, b, c, d);
  }
  else {  //else just take the ones we already have.
    for(ordinal_type n = 0; n<num_points; n++){
      a[n] = alpha[n];
      b[n] = beta[n];
      c[n] = delta[n];
      d[n] = gamma[n];
    }
    if (!normalize)
      normalizeRecurrenceCoefficients(a, b, c, d);
  }

  // With normalized coefficients, A is symmetric and tri-diagonal, with
  // diagonal = a, and off-diagonal = b, starting at b[1]
  Teuchos::SerialDenseMatrix<ordinal_type,value_type> eig_vectors(num_points,
                                                                  num_points);
  Teuchos::Array<value_type> workspace(2*num_points);
  ordinal_type info_flag;
  Teuchos::LAPACK<ordinal_type,value_type> my_lapack;

  // compute the eigenvalues (stored in a) and right eigenvectors.
  if (num_points == 1)
    my_lapack.STEQR('I', num_points, &a[0], &b[0], eig_vectors.values(),
                    num_points, &workspace[0], &info_flag);
  else
    my_lapack.STEQR('I', num_points, &a[0], &b[1], eig_vectors.values(),
                    num_points, &workspace[0], &info_flag);

  // eigenvalues are sorted by STEQR
  quad_points.resize(num_points);
  quad_weights.resize(num_points);
  for (ordinal_type i=0; i<num_points; i++) {
    quad_points[i] = a[i];
    if (std::abs(quad_points[i]) < quad_zero_tol)
      quad_points[i] = 0.0;
    quad_weights[i] = beta[0]*eig_vectors[i][0]*eig_vectors[i][0];
  }

  // Evalute basis at gauss points
  quad_values.resize(num_points);
  for (ordinal_type i=0; i<num_points; i++) {
    quad_values[i].resize(p+1);
    evaluateBases(quad_points[i], quad_values[i]);
  }
}

template <typename ordinal_type, typename value_type>
ordinal_type
Stokhos::RecurrenceBasis<ordinal_type,value_type>::
quadDegreeOfExactness(ordinal_type n) const
{
  return ordinal_type(2)*n-ordinal_type(1);
}

template <typename ordinal_type, typename value_type>
ordinal_type
Stokhos::RecurrenceBasis<ordinal_type,value_type>::
coefficientGrowth(ordinal_type n) const
{
  if (growth == SLOW_GROWTH)
    return n;

  // else moderate
  return ordinal_type(2)*n;
}

template <typename ordinal_type, typename value_type>
ordinal_type
Stokhos::RecurrenceBasis<ordinal_type,value_type>::
pointGrowth(ordinal_type n) const
{
  if (growth == SLOW_GROWTH) {
    if (n % ordinal_type(2) == ordinal_type(1))
      return n + ordinal_type(1);
    else
      return n;
  }

  // else moderate
  return n;
}

template <typename ordinal_type, typename value_type>
void
Stokhos::RecurrenceBasis<ordinal_type,value_type>::
getRecurrenceCoefficients(Teuchos::Array<value_type>& a,
                          Teuchos::Array<value_type>& b,
                          Teuchos::Array<value_type>& c,
                          Teuchos::Array<value_type>& g) const
{
  a = alpha;
  b = beta;
  c = delta;
  g = gamma;
}

template <typename ordinal_type, typename value_type>
void
Stokhos::RecurrenceBasis<ordinal_type,value_type>::
normalizeRecurrenceCoefficients(
  Teuchos::Array<value_type>& a,
  Teuchos::Array<value_type>& b,
  Teuchos::Array<value_type>& c,
  Teuchos::Array<value_type>& g) const
{
  ordinal_type n = a.size();
  a[0] = a[0] / c[0];
  b[0] = std::sqrt(b[0]);
  g[0] = b[0];
  for (ordinal_type k=1; k<n; k++) {
    a[k] = a[k] / c[k];
    b[k] = std::sqrt((b[k]*g[k])/(c[k]*c[k-1]));
    g[k] = b[k];
  }
  for (ordinal_type k=0; k<n; k++)
    c[k] = value_type(1);
}