This file is indexed.

/usr/include/deal.II/base/tensor_product_polynomials_bubbles.h is in libdeal.ii-dev 8.4.2-2+b1.

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
// ---------------------------------------------------------------------
// $Id$
//
// Copyright (C) 2012 - 2016 by the deal.II authors
//
// This file is part of the deal.II library.
//
// The deal.II library is free software; you can use it, redistribute
// it, and/or modify it under the terms of the GNU Lesser General
// Public License as published by the Free Software Foundation; either
// version 2.1 of the License, or (at your option) any later version.
// The full text of the license can be found in the file LICENSE at
// the top level of the deal.II distribution.
//
// ---------------------------------------------------------------------

#ifndef dealii__tensor_product_polynomials_bubbles_h
#define dealii__tensor_product_polynomials_bubbles_h


#include <deal.II/base/config.h>
#include <deal.II/base/exceptions.h>
#include <deal.II/base/tensor.h>
#include <deal.II/base/point.h>
#include <deal.II/base/polynomial.h>
#include <deal.II/base/tensor_product_polynomials.h>
#include <deal.II/base/utilities.h>

#include <vector>

DEAL_II_NAMESPACE_OPEN


/**
 * @addtogroup Polynomials
 * @{
 */

/**
 * Tensor product of given polynomials and bubble functions of form
 * $(2*x_j-1)^{degree-1}\prod_{i=0}^{dim-1}(x_i(1-x_i))$. This class inherits
 * most of its functionality from TensorProductPolynomials. The bubble
 * enrichments are added for the last indices. index.
 *
 * @author Daniel Arndt, 2015
 */
template <int dim>
class TensorProductPolynomialsBubbles : public TensorProductPolynomials<dim>
{
public:
  /**
   * Access to the dimension of this object, for checking and automatic
   * setting of dimension in other classes.
   */
  static const unsigned int dimension = dim;

  /**
   * Constructor. <tt>pols</tt> is a vector of objects that should be derived
   * or otherwise convertible to one-dimensional polynomial objects. It will
   * be copied element by element into a private variable.
   */
  template <class Pol>
  TensorProductPolynomialsBubbles (const std::vector<Pol> &pols);

  /**
   * Computes the value and the first and second derivatives of each tensor
   * product polynomial at <tt>unit_point</tt>.
   *
   * The size of the vectors must either be equal 0 or equal n(). In the first
   * case, the function will not compute these values.
   *
   * If you need values or derivatives of all tensor product polynomials then
   * use this function, rather than using any of the compute_value(),
   * compute_grad() or compute_grad_grad() functions, see below, in a loop
   * over all tensor product polynomials.
   */
  void compute (const Point<dim>            &unit_point,
                std::vector<double>         &values,
                std::vector<Tensor<1,dim> > &grads,
                std::vector<Tensor<2,dim> > &grad_grads,
                std::vector<Tensor<3,dim> > &third_derivatives,
                std::vector<Tensor<4,dim> > &fourth_derivatives) const;

  /**
   * Computes the value of the <tt>i</tt>th tensor product polynomial at
   * <tt>unit_point</tt>. Here <tt>i</tt> is given in tensor product
   * numbering.
   *
   * Note, that using this function within a loop over all tensor product
   * polynomials is not efficient, because then each point value of the
   * underlying (one-dimensional) polynomials is (unnecessarily) computed
   * several times.  Instead use the compute() function with
   * <tt>values.size()==</tt>n() to get the point values of all tensor
   * polynomials all at once and in a much more efficient way.
   */
  double compute_value (const unsigned int i,
                        const Point<dim> &p) const;

  /**
   * Computes the order @p order derivative of the <tt>i</tt>th tensor product
   * polynomial at <tt>unit_point</tt>. Here <tt>i</tt> is given in tensor
   * product numbering.
   *
   * Note, that using this function within a loop over all tensor product
   * polynomials is not efficient, because then each derivative value of the
   * underlying (one-dimensional) polynomials is (unnecessarily) computed
   * several times.  Instead use the compute() function, see above, with the
   * size of the appropriate parameter set to n() to get the point value of
   * all tensor polynomials all at once and in a much more efficient way.
   */
  template <int order>
  Tensor<order,dim> compute_derivative (const unsigned int i,
                                        const Point<dim> &p) const;

  /**
   * Computes the grad of the <tt>i</tt>th tensor product polynomial at
   * <tt>unit_point</tt>. Here <tt>i</tt> is given in tensor product
   * numbering.
   *
   * Note, that using this function within a loop over all tensor product
   * polynomials is not efficient, because then each derivative value of the
   * underlying (one-dimensional) polynomials is (unnecessarily) computed
   * several times.  Instead use the compute() function, see above, with
   * <tt>grads.size()==</tt>n() to get the point value of all tensor
   * polynomials all at once and in a much more efficient way.
   */
  Tensor<1,dim> compute_grad (const unsigned int i,
                              const Point<dim> &p) const;

  /**
   * Computes the second derivative (grad_grad) of the <tt>i</tt>th tensor
   * product polynomial at <tt>unit_point</tt>. Here <tt>i</tt> is given in
   * tensor product numbering.
   *
   * Note, that using this function within a loop over all tensor product
   * polynomials is not efficient, because then each derivative value of the
   * underlying (one-dimensional) polynomials is (unnecessarily) computed
   * several times.  Instead use the compute() function, see above, with
   * <tt>grad_grads.size()==</tt>n() to get the point value of all tensor
   * polynomials all at once and in a much more efficient way.
   */
  Tensor<2,dim> compute_grad_grad (const unsigned int i,
                                   const Point<dim> &p) const;

  /**
   * Returns the number of tensor product polynomials plus the bubble
   * enrichments. For <i>n</i> 1d polynomials this is <i>n<sup>dim</sup>+1</i>
   * if the maximum degree of the polynomials is one and
   * <i>n<sup>dim</sup>+dim</i> otherwise.
   */
  unsigned int n () const;
};

/** @} */


/* ---------------- template and inline functions ---------- */

#ifndef DOXYGEN

template <int dim>
template <class Pol>
inline
TensorProductPolynomialsBubbles<dim>::
TensorProductPolynomialsBubbles(const std::vector<Pol> &pols)
  :
  TensorProductPolynomials<dim>(pols)
{
  const unsigned int q_degree = this->polynomials.size()-1;
  const unsigned int n_bubbles = ((q_degree<=1)?1:dim);
  // append index for renumbering
  for (unsigned int i=0; i<n_bubbles; ++i)
    {
      this->index_map.push_back(i+this->n_tensor_pols);
      this->index_map_inverse.push_back(i+this->n_tensor_pols);
    }
}



template <int dim>
inline
unsigned int
TensorProductPolynomialsBubbles<dim>::n() const
{
  return this->n_tensor_pols+dim;
}



template <>
inline
unsigned int
TensorProductPolynomialsBubbles<0>::n() const
{
  return numbers::invalid_unsigned_int;
}

template <int dim>
template <int order>
Tensor<order,dim>
TensorProductPolynomialsBubbles<dim>::compute_derivative (const unsigned int i,
                                                          const Point<dim> &p) const
{
  const unsigned int q_degree = this->polynomials.size()-1;
  const unsigned int max_q_indices = this->n_tensor_pols;
  const unsigned int n_bubbles = ((q_degree<=1)?1:dim);
  (void)n_bubbles;
  Assert (i<max_q_indices+n_bubbles, ExcInternalError());

  // treat the regular basis functions
  if (i<max_q_indices)
    return this->TensorProductPolynomials<dim>::template compute_derivative<order>(i,p);

  const unsigned int comp = i - this->n_tensor_pols;

  Tensor<order,dim> derivative;
  switch (order)
    {
    case 1:
    {
      Tensor<1,dim> &derivative_1 = *reinterpret_cast<Tensor<1,dim>*>(&derivative);

      for (unsigned int d=0; d<dim ; ++d)
        {
          derivative_1[d] = 1.;
          //compute grad(4*\prod_{i=1}^d (x_i(1-x_i)))(p)
          for (unsigned j=0; j<dim; ++j)
            derivative_1[d] *= (d==j ? 4*(1-2*p(j)) : 4*p(j)*(1-p(j)));
          // and multiply with (2*x_i-1)^{r-1}
          for (unsigned int i=0; i<q_degree-1; ++i)
            derivative_1[d]*=2*p(comp)-1;
        }

      if (q_degree>=2)
        {
          //add \prod_{i=1}^d 4*(x_i(1-x_i))(p)
          double value=1.;
          for (unsigned int j=0; j < dim; ++j)
            value*=4*p(j)*(1-p(j));
          //and multiply with grad(2*x_i-1)^{r-1}
          double tmp=value*2*(q_degree-1);
          for (unsigned int i=0; i<q_degree-2; ++i)
            tmp*=2*p(comp)-1;
          derivative_1[comp]+=tmp;
        }

      return derivative;
    }
    case 2:
    {
      Tensor<2,dim> &derivative_2 = *reinterpret_cast<Tensor<2,dim>*>(&derivative);

      double v [dim+1][3];
      {
        for (unsigned int c=0; c<dim; ++c)
          {
            v[c][0] = 4*p(c)*(1-p(c));
            v[c][1] = 4*(1-2*p(c));
            v[c][2] = -8;
          }

        double tmp=1.;
        for (unsigned int i=0; i<q_degree-1; ++i)
          tmp *= 2*p(comp)-1;
        v[dim][0] = tmp;

        if (q_degree>=2)
          {
            double tmp = 2*(q_degree-1);
            for (unsigned int i=0; i<q_degree-2; ++i)
              tmp *= 2*p(comp)-1;
            v[dim][1] = tmp;
          }
        else
          v[dim][1] = 0.;

        if (q_degree>=3)
          {
            double tmp=4*(q_degree-2)*(q_degree-1);
            for (unsigned int i=0; i<q_degree-3; ++i)
              tmp *= 2*p(comp)-1;
            v[dim][2] = tmp;
          }
        else
          v[dim][2] = 0.;
      }

      //calculate (\partial_j \partial_k \psi) * monomial
      Tensor<2,dim> grad_grad_1;
      for (unsigned int d1=0; d1<dim; ++d1)
        for (unsigned int d2=0; d2<dim; ++d2)
          {
            grad_grad_1[d1][d2] = v[dim][0];
            for (unsigned int x=0; x<dim; ++x)
              {
                unsigned int derivative=0;
                if (d1==x || d2==x)
                  {
                    if (d1==d2)
                      derivative=2;
                    else
                      derivative=1;
                  }
                grad_grad_1[d1][d2] *= v[x][derivative];
              }
          }

      //calculate (\partial_j  \psi) *(\partial_k monomial)
      // and (\partial_k  \psi) *(\partial_j monomial)
      Tensor<2,dim> grad_grad_2;
      Tensor<2,dim> grad_grad_3;
      for (unsigned int d=0; d<dim; ++d)
        {
          grad_grad_2[d][comp] = v[dim][1];
          grad_grad_3[comp][d] = v[dim][1];
          for (unsigned int x=0; x<dim; ++x)
            {
              grad_grad_2[d][comp] *= v[x][d==x];
              grad_grad_3[comp][d] *= v[x][d==x];
            }
        }

      //calculate \psi *(\partial j \partial_k monomial) and sum
      double psi_value = 1.;
      for (unsigned int x=0; x<dim; ++x)
        psi_value *= v[x][0];

      for (unsigned int d1=0; d1<dim; ++d1)
        for (unsigned int d2=0; d2<dim; ++d2)
          derivative_2[d1][d2] = grad_grad_1[d1][d2]
                                 +grad_grad_2[d1][d2]
                                 +grad_grad_3[d1][d2];
      derivative_2[comp][comp]+=psi_value*v[dim][2];

      return derivative;
    }
    default:
    {
      Assert (false, ExcNotImplemented());
      return derivative;
    }
    }
}


#endif // DOXYGEN
DEAL_II_NAMESPACE_CLOSE

#endif