/usr/include/deal.II/fe/fe_poly.templates.h is in libdeal.ii-dev 8.1.0-4.
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 483 484 485 486 487 488 489 490 491 492 493 494 495 496 497 498 499 500 501 502 503 504 505 506 | // ---------------------------------------------------------------------
// $Id: fe_poly.templates.h 30036 2013-07-18 16:55:32Z maier $
//
// Copyright (C) 2006 - 2013 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.
//
// ---------------------------------------------------------------------
#include <deal.II/base/qprojector.h>
#include <deal.II/base/polynomial_space.h>
#include <deal.II/base/tensor_product_polynomials.h>
#include <deal.II/fe/fe_values.h>
#include <deal.II/fe/fe_poly.h>
DEAL_II_NAMESPACE_OPEN
template <class POLY, int dim, int spacedim>
FE_Poly<POLY,dim,spacedim>::FE_Poly (const POLY &poly_space,
const FiniteElementData<dim> &fe_data,
const std::vector<bool> &restriction_is_additive_flags,
const std::vector<ComponentMask> &nonzero_components):
FiniteElement<dim,spacedim> (fe_data,
restriction_is_additive_flags,
nonzero_components),
poly_space(poly_space)
{
AssertDimension(dim, POLY::dimension);
}
template <class POLY, int dim, int spacedim>
unsigned int
FE_Poly<POLY,dim,spacedim>::get_degree () const
{
return this->degree;
}
template <class POLY, int dim, int spacedim>
double
FE_Poly<POLY,dim,spacedim>::shape_value (const unsigned int i,
const Point<dim> &p) const
{
Assert (i<this->dofs_per_cell, ExcIndexRange(i,0,this->dofs_per_cell));
return poly_space.compute_value(i, p);
}
template <class POLY, int dim, int spacedim>
double
FE_Poly<POLY,dim,spacedim>::shape_value_component (const unsigned int i,
const Point<dim> &p,
const unsigned int component) const
{
Assert (i<this->dofs_per_cell, ExcIndexRange(i,0,this->dofs_per_cell));
Assert (component == 0, ExcIndexRange (component, 0, 1));
return poly_space.compute_value(i, p);
}
template <class POLY, int dim, int spacedim>
Tensor<1,dim>
FE_Poly<POLY,dim,spacedim>::shape_grad (const unsigned int i,
const Point<dim> &p) const
{
Assert (i<this->dofs_per_cell, ExcIndexRange(i,0,this->dofs_per_cell));
return poly_space.compute_grad(i, p);
}
template <class POLY, int dim, int spacedim>
Tensor<1,dim>
FE_Poly<POLY,dim,spacedim>::shape_grad_component (const unsigned int i,
const Point<dim> &p,
const unsigned int component) const
{
Assert (i<this->dofs_per_cell, ExcIndexRange(i,0,this->dofs_per_cell));
Assert (component == 0, ExcIndexRange (component, 0, 1));
return poly_space.compute_grad(i, p);
}
template <class POLY, int dim, int spacedim>
Tensor<2,dim>
FE_Poly<POLY,dim,spacedim>::shape_grad_grad (const unsigned int i,
const Point<dim> &p) const
{
Assert (i<this->dofs_per_cell, ExcIndexRange(i,0,this->dofs_per_cell));
return poly_space.compute_grad_grad(i, p);
}
template <class POLY, int dim, int spacedim>
Tensor<2,dim>
FE_Poly<POLY,dim,spacedim>::shape_grad_grad_component (const unsigned int i,
const Point<dim> &p,
const unsigned int component) const
{
Assert (i<this->dofs_per_cell, ExcIndexRange(i,0,this->dofs_per_cell));
Assert (component == 0, ExcIndexRange (component, 0, 1));
return poly_space.compute_grad_grad(i, p);
}
//---------------------------------------------------------------------------
// Auxiliary functions
//---------------------------------------------------------------------------
template <class POLY, int dim, int spacedim>
UpdateFlags
FE_Poly<POLY,dim,spacedim>::update_once (const UpdateFlags flags) const
{
// for this kind of elements, only
// the values can be precomputed
// once and for all. set this flag
// if the values are requested at
// all
return (update_default | (flags & update_values));
}
template <class POLY, int dim, int spacedim>
UpdateFlags
FE_Poly<POLY,dim,spacedim>::update_each (const UpdateFlags flags) const
{
UpdateFlags out = update_default;
if (flags & update_gradients)
out |= update_gradients | update_covariant_transformation;
if (flags & update_hessians)
out |= update_hessians | update_covariant_transformation;
if (flags & update_cell_normal_vectors)
out |= update_cell_normal_vectors | update_JxW_values;
return out;
}
//---------------------------------------------------------------------------
// Data field initialization
//---------------------------------------------------------------------------
template <class POLY, int dim, int spacedim>
typename Mapping<dim,spacedim>::InternalDataBase *
FE_Poly<POLY,dim,spacedim>::get_data (const UpdateFlags update_flags,
const Mapping<dim,spacedim> &mapping,
const Quadrature<dim> &quadrature) const
{
// generate a new data object and
// initialize some fields
InternalData *data = new InternalData;
// check what needs to be
// initialized only once and what
// on every cell/face/subface we
// visit
data->update_once = update_once(update_flags);
data->update_each = update_each(update_flags);
data->update_flags = data->update_once | data->update_each;
const UpdateFlags flags(data->update_flags);
const unsigned int n_q_points = quadrature.size();
// some scratch arrays
std::vector<double> values(0);
std::vector<Tensor<1,dim> > grads(0);
std::vector<Tensor<2,dim> > grad_grads(0);
// initialize fields only if really
// necessary. otherwise, don't
// allocate memory
if (flags & update_values)
{
values.resize (this->dofs_per_cell);
data->shape_values.resize (this->dofs_per_cell,
std::vector<double> (n_q_points));
}
if (flags & update_gradients)
{
grads.resize (this->dofs_per_cell);
data->shape_gradients.resize (this->dofs_per_cell,
std::vector<Tensor<1,dim> > (n_q_points));
}
// if second derivatives through
// finite differencing is required,
// then initialize some objects for
// that
if (flags & update_hessians)
data->initialize_2nd (this, mapping, quadrature);
// next already fill those fields
// of which we have information by
// now. note that the shape
// gradients are only those on the
// unit cell, and need to be
// transformed when visiting an
// actual cell
if (flags & (update_values | update_gradients))
for (unsigned int i=0; i<n_q_points; ++i)
{
poly_space.compute(quadrature.point(i),
values, grads, grad_grads);
if (flags & update_values)
for (unsigned int k=0; k<this->dofs_per_cell; ++k)
data->shape_values[k][i] = values[k];
if (flags & update_gradients)
for (unsigned int k=0; k<this->dofs_per_cell; ++k)
data->shape_gradients[k][i] = grads[k];
}
return data;
}
//---------------------------------------------------------------------------
// Fill data of FEValues
//---------------------------------------------------------------------------
template <class POLY, int dim, int spacedim>
void
FE_Poly<POLY,dim,spacedim>::fill_fe_values
(const Mapping<dim,spacedim> &mapping,
const typename Triangulation<dim,spacedim>::cell_iterator &cell,
const Quadrature<dim> &quadrature,
typename Mapping<dim,spacedim>::InternalDataBase &mapping_data,
typename Mapping<dim,spacedim>::InternalDataBase &fedata,
FEValuesData<dim,spacedim> &data,
CellSimilarity::Similarity &cell_similarity) const
{
// convert data object to internal
// data for this class. fails with
// an exception if that is not
// possible
Assert (dynamic_cast<InternalData *> (&fedata) != 0, ExcInternalError());
InternalData &fe_data = static_cast<InternalData &> (fedata);
const UpdateFlags flags(fe_data.current_update_flags());
for (unsigned int k=0; k<this->dofs_per_cell; ++k)
{
if (flags & update_values)
for (unsigned int i=0; i<quadrature.size(); ++i)
data.shape_values(k,i) = fe_data.shape_values[k][i];
if (flags & update_gradients && cell_similarity != CellSimilarity::translation)
mapping.transform(fe_data.shape_gradients[k], data.shape_gradients[k],
mapping_data, mapping_covariant);
}
if (flags & update_hessians && cell_similarity != CellSimilarity::translation)
this->compute_2nd (mapping, cell, QProjector<dim>::DataSetDescriptor::cell(),
mapping_data, fe_data, data);
}
template <class POLY, int dim, int spacedim>
void
FE_Poly<POLY,dim,spacedim>::
fill_fe_face_values (const Mapping<dim,spacedim> &mapping,
const typename Triangulation<dim,spacedim>::cell_iterator &cell,
const unsigned int face,
const Quadrature<dim-1> &quadrature,
typename Mapping<dim,spacedim>::InternalDataBase &mapping_data,
typename Mapping<dim,spacedim>::InternalDataBase &fedata,
FEValuesData<dim,spacedim> &data) const
{
// convert data object to internal
// data for this class. fails with
// an exception if that is not
// possible
Assert (dynamic_cast<InternalData *> (&fedata) != 0, ExcInternalError());
InternalData &fe_data = static_cast<InternalData &> (fedata);
// offset determines which data set
// to take (all data sets for all
// faces are stored contiguously)
const typename QProjector<dim>::DataSetDescriptor offset
= QProjector<dim>::DataSetDescriptor::face (face,
cell->face_orientation(face),
cell->face_flip(face),
cell->face_rotation(face),
quadrature.size());
const UpdateFlags flags(fe_data.update_once | fe_data.update_each);
for (unsigned int k=0; k<this->dofs_per_cell; ++k)
{
if (flags & update_values)
for (unsigned int i=0; i<quadrature.size(); ++i)
data.shape_values(k,i) = fe_data.shape_values[k][i+offset];
if (flags & update_gradients)
mapping.transform(make_slice(fe_data.shape_gradients[k], offset, quadrature.size()),
data.shape_gradients[k],
mapping_data, mapping_covariant);
}
if (flags & update_hessians)
this->compute_2nd (mapping, cell, offset, mapping_data, fe_data, data);
}
//codimension 1
// template <>
// inline
// void
// FE_Poly<TensorProductPolynomials<1>,1,2>::fill_fe_subface_values (const Mapping<1,2> &,
// const Triangulation<1,2>::cell_iterator &,
// const unsigned int,
// const unsigned int,
// const Quadrature<0> &,
// Mapping<1,2>::InternalDataBase &,
// Mapping<1,2>::InternalDataBase &,
// FEValuesData<1,2> &) const
// {
// AssertThrow(false, ExcNotImplemented());
// }
// template <>
// inline
// void
// FE_Poly<TensorProductPolynomials<2>,2,3>::fill_fe_subface_values (const Mapping<2,3> &,
// const Triangulation<2,3>::cell_iterator &,
// const unsigned int,
// const unsigned int,
// const Quadrature<1> &,
// Mapping<2,3>::InternalDataBase &,
// Mapping<2,3>::InternalDataBase &,
// FEValuesData<2,3> &) const
// {
// AssertThrow(false, ExcNotImplemented());
// }
// template <>
// inline
// void
// FE_Poly<PolynomialSpace<1>,1,2>::fill_fe_subface_values (const Mapping<1,2> &,
// const Triangulation<1,2>::cell_iterator &,
// const unsigned int,
// const unsigned int,
// const Quadrature<0> &,
// Mapping<1,2>::InternalDataBase &,
// Mapping<1,2>::InternalDataBase &,
// FEValuesData<1,2> &) const
// {
// AssertThrow(false, ExcNotImplemented());
// }
// template <>
// inline
// void
// FE_Poly<PolynomialSpace<2>,2,3>::fill_fe_subface_values (const Mapping<2,3> &,
// const Triangulation<2,3>::cell_iterator &,
// const unsigned int,
// const unsigned int,
// const Quadrature<1> &,
// Mapping<2,3>::InternalDataBase &,
// Mapping<2,3>::InternalDataBase &,
// FEValuesData<2,3> &) const
// {
// AssertThrow(false, ExcNotImplemented());
// }
template <class POLY, int dim, int spacedim>
void
FE_Poly<POLY,dim,spacedim>::fill_fe_subface_values (const Mapping<dim,spacedim> &mapping,
const typename Triangulation<dim,spacedim>::cell_iterator &cell,
const unsigned int face,
const unsigned int subface,
const Quadrature<dim-1> &quadrature,
typename Mapping<dim,spacedim>::InternalDataBase &mapping_data,
typename Mapping<dim,spacedim>::InternalDataBase &fedata,
FEValuesData<dim,spacedim> &data) const
{
// convert data object to internal
// data for this class. fails with
// an exception if that is not
// possible
Assert (dynamic_cast<InternalData *> (&fedata) != 0, ExcInternalError());
InternalData &fe_data = static_cast<InternalData &> (fedata);
// offset determines which data set
// to take (all data sets for all
// sub-faces are stored contiguously)
const typename QProjector<dim>::DataSetDescriptor offset
= QProjector<dim>::DataSetDescriptor::subface (face, subface,
cell->face_orientation(face),
cell->face_flip(face),
cell->face_rotation(face),
quadrature.size(),
cell->subface_case(face));
const UpdateFlags flags(fe_data.update_once | fe_data.update_each);
for (unsigned int k=0; k<this->dofs_per_cell; ++k)
{
if (flags & update_values)
for (unsigned int i=0; i<quadrature.size(); ++i)
data.shape_values(k,i) = fe_data.shape_values[k][i+offset];
if (flags & update_gradients)
mapping.transform(make_slice(fe_data.shape_gradients[k], offset, quadrature.size()),
data.shape_gradients[k],
mapping_data, mapping_covariant);
}
if (flags & update_hessians)
this->compute_2nd (mapping, cell, offset, mapping_data, fe_data, data);
}
namespace internal
{
template <class POLY>
inline
std::vector<unsigned int>
get_poly_space_numbering (const POLY &)
{
Assert (false, ExcNotImplemented());
return std::vector<unsigned int>();
}
template <class POLY>
inline
std::vector<unsigned int>
get_poly_space_numbering_inverse (const POLY &)
{
Assert (false, ExcNotImplemented());
return std::vector<unsigned int>();
}
template <int dim, typename POLY>
inline
std::vector<unsigned int>
get_poly_space_numbering (const TensorProductPolynomials<dim,POLY> &poly)
{
return poly.get_numbering();
}
template <int dim, typename POLY>
inline
std::vector<unsigned int>
get_poly_space_numbering_inverse (const TensorProductPolynomials<dim,POLY> &poly)
{
return poly.get_numbering_inverse();
}
}
template <class POLY, int dim, int spacedim>
std::vector<unsigned int>
FE_Poly<POLY,dim,spacedim>::get_poly_space_numbering () const
{
return internal::get_poly_space_numbering (poly_space);
}
template <class POLY, int dim, int spacedim>
std::vector<unsigned int>
FE_Poly<POLY,dim,spacedim>::get_poly_space_numbering_inverse () const
{
return internal::get_poly_space_numbering_inverse (poly_space);
}
DEAL_II_NAMESPACE_CLOSE
|