This file is indexed.

/usr/include/deal.II/fe/mapping_q1.h is in libdeal.ii-dev 6.3.1-1.1.

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
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
//---------------------------------------------------------------------------
//    $Id: mapping_q1.h 20158 2009-11-24 00:03:32Z kronbichler $
//    Version: $Name$
//
//    Copyright (C) 2000, 2001, 2002, 2003, 2004, 2005, 2006, 2007, 2008, 2009 by the deal.II authors
//
//    This file is subject to QPL and may not be  distributed
//    without copyright and license information. Please refer
//    to the file deal.II/doc/license.html for the  text  and
//    further information on this license.
//
//---------------------------------------------------------------------------
#ifndef __deal2__mapping_q1_h
#define __deal2__mapping_q1_h


#include <base/config.h>
#include <base/table.h>
#include <base/qprojector.h>
#include <grid/tria_iterator.h>
#include <dofs/dof_accessor.h>
#include <fe/mapping.h>

#include <cmath>

DEAL_II_NAMESPACE_OPEN

/*!@addtogroup mapping */
/*@{*/


/**
 * Mapping of general quadrilateral/hexahedra by d-linear shape
 * functions.
 *
 * This function maps the unit cell to a general grid cell with
 * straight lines in $d$ dimensions (remark that in 3D the surfaces
 * may be curved, even if the edges are not). This is the well-known
 * mapping for polyhedral domains.
 *
 * Shape function for this mapping are the same as for the finite
 * element FE_Q of order 1. Therefore, coupling these two yields
 * an isoparametric element.
 *
 * For more information about the <tt>spacedim</tt> template parameter
 * check the documentation of FiniteElement or the one of
 * Triangulation.
 *
 * @author Guido Kanschat, 2000, 2001; Ralf Hartmann, 2000, 2001, 2005
 */
template <int dim, int spacedim=dim>
class MappingQ1 : public Mapping<dim,spacedim>
{
  public:
				     /**
				      * Default constructor.
				      */
    MappingQ1 ();

    virtual Point<spacedim>
    transform_unit_to_real_cell (
      const typename Triangulation<dim,spacedim>::cell_iterator &cell,
      const Point<dim>                                 &p) const;

				     /**
				      * Transforms the point @p p on
				      * the real cell to the point
				      * @p p_unit on the unit cell
				      * @p cell and returns @p p_unit.
				      *
				      * Uses Newton iteration and the
				      * @p transform_unit_to_real_cell
				      * function.
				      */
    virtual Point<dim>
    transform_real_to_unit_cell (
      const typename Triangulation<dim,spacedim>::cell_iterator &cell,
      const Point<spacedim>                            &p) const;

    virtual void
    transform (const VectorSlice<const std::vector<Tensor<1,dim> > > input,
               VectorSlice<std::vector<Tensor<1,spacedim> > > output,
               const typename Mapping<dim,spacedim>::InternalDataBase &internal,
	       const MappingType type) const;

    virtual void
    transform (const VectorSlice<const std::vector<Tensor<2,dim> > > input,
               VectorSlice<std::vector<Tensor<2,spacedim> > > output,
               const typename Mapping<dim,spacedim>::InternalDataBase &internal,
	       const MappingType type) const;

                                     /**
                                      * Return a pointer to a copy of the
                                      * present object. The caller of this
                                      * copy then assumes ownership of it.
                                      */
    virtual
    Mapping<dim,spacedim> * clone () const;

				     /**
				      * Storage for internal data of
				      * d-linear transformation.
				      */
    class InternalData : public Mapping<dim,spacedim>::InternalDataBase
    {
      public:
					 /**
					  * Constructor. Pass the
					  * number of shape functions.
					  */
	InternalData(const unsigned int n_shape_functions);

					 /**
					  * Shape function at quadrature
					  * point. Shape functions are
					  * in tensor product order, so
					  * vertices must be reordered
					  * to obtain transformation.
					  */
	double shape (const unsigned int qpoint,
		      const unsigned int shape_nr) const;

					 /**
					  * Shape function at quadrature
					  * point. See above.
					  */
	double &shape (const unsigned int qpoint,
		       const unsigned int shape_nr);

					 /**
					  * Gradient of shape function
					  * in quadrature point. See
					  * above.
					  */
	Tensor<1,dim> derivative (const unsigned int qpoint,
				  const unsigned int shape_nr) const;

					 /**
					  * Gradient of shape function
					  * in quadrature point. See
					  * above.
					  */
	Tensor<1,dim> &derivative (const unsigned int qpoint,
				   const unsigned int shape_nr);

					 /**
					  * Second derivative of shape
					  * function in quadrature
					  * point. See above.
					  */
	Tensor<2,dim> second_derivative (const unsigned int qpoint,
					 const unsigned int shape_nr) const;

					 /**
					  * Second derivative of shape
					  * function in quadrature
					  * point. See above.
					  */
	Tensor<2,dim> &second_derivative (const unsigned int qpoint,
					  const unsigned int shape_nr);

					 /**
					  * Return an estimate (in
					  * bytes) or the memory
					  * consumption of this
					  * object.
					  */
	virtual unsigned int memory_consumption () const;

					 /**
					  * Values of shape
					  * functions. Access by
					  * function @p shape.
					  *
					  * Computed once.
					  */
	std::vector<double> shape_values;

					 /**
					  * Values of shape function
					  * derivatives. Access by
					  * function @p derivative.
					  *
					  * Computed once.
					  */
	std::vector<Tensor<1,dim> > shape_derivatives;

					 /**
					  * Values of shape function
					  * second derivatives. Access
					  * by function
					  * @p second_derivative.
					  *
					  * Computed once.
					  */
	std::vector<Tensor<2,dim> > shape_second_derivatives;

					 /**
					  * Tensors of covariant
					  * transformation at each of
					  * the quadrature points. The
					  * matrix stored is the
					  * inverse of the Jacobian
					  * matrix, which itself is
					  * stored in the
					  * @p contravariant field of
					  * this structure.
					  *
					  * Computed on each cell.
					  */
	std::vector<Tensor<2,spacedim> > covariant;

					 /**
					  * Tensors of contravariant
					  * transformation at each of
					  * the quadrature points. The
					  * contravariant matrix is
					  * the Jacobian of the
					  * transformation,
					  * i.e. $J_{ij}=dx_i/d\hat x_j$.
					  *
					  * Computed on each cell.
					  */
	std::vector<Tensor<2,spacedim> > contravariant;

					 /**
					  * Unit tangential vectors. Used
					  * for the computation of
					  * boundary forms and normal
					  * vectors.
					  *
					  * Filled once.
					  */
        std::vector<std::vector<Tensor<1,spacedim> > > unit_tangentials;

					 /**
					  * Auxiliary vectors for internal use.
					  */
        std::vector<std::vector<Tensor<1,spacedim> > > aux;

					 /**
					  * Stores the support points of
					  * the mapping shape functions on
					  * the @p cell_of_current_support_points.
					  */
	std::vector<Point<spacedim> > mapping_support_points;

					 /**
					  * Stores the cell of which the
					  * @p mapping_support_points are
					  * stored.
					  */
	typename Triangulation<dim,spacedim>::cell_iterator cell_of_current_support_points;

					 /**
					  * Default value of this flag
					  * is @p true. If <tt>*this</tt>
					  * is an object of a derived
					  * class, this flag is set to
					  * @p false.
					  */
	bool is_mapping_q1_data;

					 /**
					  * Number of shape
					  * functions. If this is a Q1
					  * mapping, then it is simply
					  * the number of vertices per
					  * cell. However, since also
					  * derived classes use this
					  * class (e.g. the
					  * Mapping_Q() class),
					  * the number of shape
					  * functions may also be
					  * different.
					  */
	unsigned int n_shape_functions;
    };

                                     /**
                                      * Declare a convenience typedef
                                      * for the class that describes
                                      * offsets into quadrature
                                      * formulas projected onto faces
                                      * and subfaces.
                                      */
    typedef
    typename QProjector<dim>::DataSetDescriptor
    DataSetDescriptor;

				     /**
				      * Implementation of the interface in
				      * Mapping.
				      */
    virtual void
    fill_fe_values (const typename Triangulation<dim,spacedim>::cell_iterator &cell,
		    const Quadrature<dim>                                     &quadrature,
		    typename Mapping<dim,spacedim>::InternalDataBase          &mapping_data,
		    typename std::vector<Point<spacedim> >                    &quadrature_points,
		    std::vector<double>                                       &JxW_values,
		    std::vector<Tensor<2,spacedim> >                          &jacobians,
		    std::vector<Tensor<3,spacedim> >                          &jacobian_grads,
		    std::vector<Tensor<2,spacedim> >                          &inverse_jacobians,
		    std::vector<Point<spacedim> >                             &cell_normal_vectors,
		    CellSimilarity::Similarity                           &cell_similarity) const;

				     /**
				      * Implementation of the interface in
				      * Mapping.
				      */
    virtual void
    fill_fe_face_values (const typename Triangulation<dim,spacedim>::cell_iterator &cell,
			 const unsigned int                               face_no,
			 const Quadrature<dim-1>                          &quadrature,
			 typename Mapping<dim,spacedim>::InternalDataBase &mapping_data,
			 typename std::vector<Point<dim> >                &quadrature_points,
			 std::vector<double>                              &JxW_values,
			 typename std::vector<Tensor<1,dim> >             &boundary_form,
			 typename std::vector<Point<spacedim> >           &normal_vectors) const ;

				     /**
				      * Implementation of the interface in
				      * Mapping.
				      */
    virtual void
    fill_fe_subface_values (const typename Triangulation<dim,spacedim>::cell_iterator &cell,
			    const unsigned int face_no,
			    const unsigned int sub_no,
			    const Quadrature<dim-1>& quadrature,
			    typename Mapping<dim,spacedim>::InternalDataBase &mapping_data,
			    typename std::vector<Point<dim> >        &quadrature_points,
			    std::vector<double>             &JxW_values,
			    typename std::vector<Tensor<1,dim> >        &boundary_form,
			    typename std::vector<Point<spacedim> >        &normal_vectors) const ;

				     /**
				      * Compute shape values and/or
				      * derivatives.
				      *
				      * Calls either the
				      * @p compute_shapes_virtual of
				      * this class or that of the
				      * derived class, depending on
				      * whether
				      * <tt>data.is_mapping_q1_data</tt>
				      * equals @p true or @p false.
				      */
    void compute_shapes (const std::vector<Point<dim> > &unit_points,
			 InternalData &data) const;

				     /**
				      * Do the computations for the
				      * @p get_data functions. Here,
				      * the data vectors of
				      * @p InternalData are
				      * reinitialized to proper size
				      * and shape values are computed.
				      */
    void compute_data (const UpdateFlags flags,
		       const Quadrature<dim> &quadrature,
		       const unsigned int n_orig_q_points,
		       InternalData &data) const;

				     /**
				      * Do the computations for the
				      * @p get_face_data
				      * functions. Here, the data
				      * vectors of @p InternalData
				      * are reinitialized to proper
				      * size and shape values and
				      * derivatives are
				      * computed. Furthermore
				      * @p unit_tangential vectors of
				      * the face are computed.
				      */
    void compute_face_data (const UpdateFlags flags,
			    const Quadrature<dim> &quadrature,
			    const unsigned int n_orig_q_points,
			    InternalData &data) const;

				     /**
				      * Do the computation for the
				      * <tt>fill_*</tt> functions.
				      */
    void compute_fill (const typename Triangulation<dim,spacedim>::cell_iterator &cell,
		       const unsigned int      npts,
		       const DataSetDescriptor data_set,
		       const CellSimilarity::Similarity cell_similarity,
		       InternalData           &data,
		       std::vector<Point<spacedim> > &quadrature_points) const;

				     /**
				      * Do the computation for the
				      * <tt>fill_*</tt> functions.
				      */
    void compute_fill_face (const typename Triangulation<dim,spacedim>::cell_iterator &cell,
			    const unsigned int      face_no,
			    const unsigned int      subface_no,
			    const unsigned int      npts,
			    const DataSetDescriptor data_set,
			    const std::vector<double>   &weights,
			    InternalData           &mapping_data,
			    std::vector<Point<dim> >    &quadrature_points,
			    std::vector<double>         &JxW_values,
			    std::vector<Tensor<1,dim> > &boundary_form,
			    std::vector<Point<spacedim> > &normal_vectors) const;

				     /**
				      * Compute shape values and/or
				      * derivatives.
				      */
    virtual void compute_shapes_virtual (const std::vector<Point<dim> > &unit_points,
					 InternalData &data) const;

				     /**
				      * Transforms a point @p p on
				      * the unit cell to the point
				      * @p p_real on the real cell
				      * @p cell and returns @p p_real.
				      *
				      * This function is called by
				      * @p transform_unit_to_real_cell
				      * and multiply (through the
				      * Newton iteration) by
				      * @p transform_real_to_unit_cell_internal.
				      *
				      * Takes a reference to an
				      * @p InternalData that must
				      * already include the shape
				      * values at point @p p and the
				      * mapping support points of the
				      * cell.
				      *
				      * This @p InternalData argument
				      * avoids multiple computations
				      * of the shape values at point
				      * @p p and especially multiple
				      * computations of the mapping
				      * support points.
				      */
    Point<spacedim> transform_unit_to_real_cell_internal (const InternalData &mdata) const;

				     /**
				      * Transforms the point @p p on
				      * the real cell to the point
				      * @p p_unit on the unit cell
				      * @p cell by a Newton
				      * iteration.
				      *
				      * Takes a reference to an
				      * @p InternalData that is
				      * assumed to be previously
				      * created by the @p get_data
				      * function with @p UpdateFlags
				      * including
				      * @p update_transformation_values
				      * and
				      * @p update_transformation_gradients
				      * and a one point Quadrature
				      * including the given point
				      * @p p_unit.  Hence this
				      * function assumes that
				      * @p mdata already includes the
				      * transformation shape values
				      * and gradients computed at
				      * @p p_unit.
				      *
				      * These assumptions should be
				      * fulfilled by the calling
				      * function. That is up to now
				      * only the function
				      * @p transform_real_to_unit_cell
				      * and its overloaded versions.
				      * @p mdata will be changed by
				      * this function.
				      */
    void transform_real_to_unit_cell_internal (const typename Triangulation<dim,spacedim>::cell_iterator &cell,
					       const Point<spacedim> &p,
					       InternalData &mdata,
					       Point<dim> &p_unit) const;

				     /**
				      * Always returns @p true because
				      * MappingQ1 preserves vertex locations.
				      */
    virtual
    bool preserves_vertex_locations () const;

  private:
				     /**
				      * Implementation of the interface in
				      * Mapping.
				      *
				      * Description of effects:
				      * <ul>
				      * <li> if @p update_quadrature_points
				      * is required, the output will
				      * contain
				      * @p update_transformation_values. This
				      * computes the values of the
				      * transformation basis
				      * polynomials at the unit cell
				      * quadrature points.
				      * <li> if any of
				      * @p update_covariant_transformation,
				      * @p update_contravariant_transformation,
				      * @p update_JxW_values,
				      * @p update_boundary_forms,
				      * @p update_normal_vectors is
				      * required, the output will
				      * contain
				      * @p update_transformation_gradients
				      * to compute derivatives of the
				      * transformation basis
				      * polynomials.
				      * </ul>
				      */
    virtual UpdateFlags update_once (const UpdateFlags flags) const;

				     /**
				      * Implementation of the interface in
				      * Mapping.
				      *
				      * Description of effects if
				      * @p flags contains:
				      * <ul>
				      * <li> <code>update_quadrature_points</code> is
				      * copied to the output to
				      * compute the quadrature points
				      * on the real cell.
				      * <li> <code>update_JxW_values</code> is
				      * copied and requires
				      * @p update_boundary_forms on
				      * faces. The latter, because the
				      * surface element is just the
				      * norm of the boundary form.
				      * <li> <code>update_normal_vectors</code>
				      * is copied and requires
				      * @p update_boundary_forms. The
				      * latter, because the normal
				      * vector is the normalized
				      * boundary form.
				      * <li>
				      * <code>update_covariant_transformation</code>
				      * is copied and requires
				      * @p update_contravariant_transformation,
				      * since it is computed as the
				      * inverse of the latter.
				      * <li> <code>update_JxW_values</code> is
				      * copied and requires
				      * <code>update_contravariant_transformation</code>,
				      * since it is computed as one
				      * over determinant of the
				      * latter.
				      * <li> <code>update_boundary_forms</code>
				      * is copied and requires
				      * <code>update_contravariant_transformation</code>,
				      * since the boundary form is
				      * computed as the contravariant
				      * image of the normal vector to
				      * the unit cell.
				      * </ul>
				      */
    virtual UpdateFlags update_each (const UpdateFlags flags) const;

    virtual
    typename Mapping<dim,spacedim>::InternalDataBase *
    get_data (const UpdateFlags,
	      const Quadrature<dim>& quadrature) const;

    virtual
    typename Mapping<dim,spacedim>::InternalDataBase *
    get_face_data (const UpdateFlags flags,
		   const Quadrature<dim-1>& quadrature) const;

    virtual
    typename Mapping<dim,spacedim>::InternalDataBase *
    get_subface_data (const UpdateFlags flags,
		      const Quadrature<dim-1>& quadrature) const;

				     /**
				      * Computes the support points of
				      * the mapping. For @p MappingQ1
				      * these are the
				      * vertices. However, other
				      * classes may override this
				      * function. In particular, the
				      * MappingQ1Eulerian class does
				      * exactly this by not computing
				      * the support points from the
				      * geometry of the current cell
				      * but instead evaluating an
				      * externally given displacement
				      * field in addition to the
				      * geometry of the cell.
				      */
    virtual void compute_mapping_support_points(
      const typename Triangulation<dim,spacedim>::cell_iterator &cell,
      std::vector<Point<spacedim> > &a) const;

				     /**
				      * Number of shape functions. Is
				      * simply the number of vertices
				      * per cell for the Q1 mapping.
				      */
    static const unsigned int n_shape_functions = GeometryInfo<dim>::vertices_per_cell;
};


/**
 * In order to avoid creation of static MappingQ1 objects at several
 * places in the library (in particular in backward compatibility
 * functions), we define a static MappingQ1 objects once and for all
 * places where it is needed.
 */
template <int dim, int spacedim=dim>
struct StaticMappingQ1
{
    static MappingQ1<dim, spacedim> mapping;
};


/*@}*/

/*----------------------------------------------------------------------*/

#ifndef DOXYGEN

template<int dim, int spacedim>
inline
double
MappingQ1<dim,spacedim>::InternalData::shape (const unsigned int qpoint,
				     const unsigned int shape_nr) const
{
  Assert(qpoint*n_shape_functions + shape_nr < shape_values.size(),
	 ExcIndexRange(qpoint*n_shape_functions + shape_nr, 0,
		       shape_values.size()));
  return shape_values [qpoint*n_shape_functions + shape_nr];
}



template<int dim, int spacedim>
inline
double &
MappingQ1<dim,spacedim>::InternalData::shape (const unsigned int qpoint,
				     const unsigned int shape_nr)
{
  Assert(qpoint*n_shape_functions + shape_nr < shape_values.size(),
	 ExcIndexRange(qpoint*n_shape_functions + shape_nr, 0,
		       shape_values.size()));
  return shape_values [qpoint*n_shape_functions + shape_nr];
}


template<int dim, int spacedim>
inline
Tensor<1,dim>
MappingQ1<dim,spacedim>::InternalData::derivative (const unsigned int qpoint,
					  const unsigned int shape_nr) const
{
  Assert(qpoint*n_shape_functions + shape_nr < shape_derivatives.size(),
	 ExcIndexRange(qpoint*n_shape_functions + shape_nr, 0,
		       shape_derivatives.size()));
  return shape_derivatives [qpoint*n_shape_functions + shape_nr];
}



template<int dim, int spacedim>
inline
Tensor<1,dim> &
MappingQ1<dim,spacedim>::InternalData::derivative (const unsigned int qpoint,
					  const unsigned int shape_nr)
{
  Assert(qpoint*n_shape_functions + shape_nr < shape_derivatives.size(),
	 ExcIndexRange(qpoint*n_shape_functions + shape_nr, 0,
		       shape_derivatives.size()));
  return shape_derivatives [qpoint*n_shape_functions + shape_nr];
}


template <int dim, int spacedim>
inline
Tensor<2,dim>
MappingQ1<dim,spacedim>::InternalData::second_derivative (const unsigned int qpoint,
			      			          const unsigned int shape_nr) const
{
  Assert(qpoint*n_shape_functions + shape_nr < shape_second_derivatives.size(),
	 ExcIndexRange(qpoint*n_shape_functions + shape_nr, 0,
		       shape_second_derivatives.size()));
  return shape_second_derivatives [qpoint*n_shape_functions + shape_nr];
}



template <int dim, int spacedim>
inline
Tensor<2,dim> &
MappingQ1<dim,spacedim>::InternalData::second_derivative (const unsigned int qpoint,
						 const unsigned int shape_nr)
{
  Assert(qpoint*n_shape_functions + shape_nr < shape_second_derivatives.size(),
	 ExcIndexRange(qpoint*n_shape_functions + shape_nr, 0,
		       shape_second_derivatives.size()));
  return shape_second_derivatives [qpoint*n_shape_functions + shape_nr];
}



template <int dim, int spacedim>
inline
bool
MappingQ1<dim,spacedim>::preserves_vertex_locations () const
{
  return true;
}

#endif // DOXYGEN

/* -------------- declaration of explicit specializations ------------- */


DEAL_II_NAMESPACE_CLOSE

#endif