escript  Revision_
ReferenceElementSets.h
Go to the documentation of this file.
1 
2 /*****************************************************************************
3 *
4 * Copyright (c) 2003-2016 by The University of Queensland
5 * http://www.uq.edu.au
6 *
7 * Primary Business: Queensland, Australia
8 * Licensed under the Apache License, version 2.0
9 * http://www.apache.org/licenses/LICENSE-2.0
10 *
11 * Development until 2012 by Earth Systems Science Computational Center (ESSCC)
12 * Development 2012-2013 by School of Earth Sciences
13 * Development from 2014 by Centre for Geoscience Computing (GeoComp)
14 *
15 *****************************************************************************/
16 
17 #ifndef __FINLEY_REFERENCEELEMENTSETS_H__
18 #define __FINLEY_REFERENCEELEMENTSETS_H__
19 
20 #include "ReferenceElements.h"
21 
22 namespace finley {
23 
27  ReferenceElementSet(ElementTypeId id, int order, int reduced_order)
28  {
31  id_info->BasisFunctions);
32  if (order<0)
33  order=std::max(2*bf_info->numOrder, 0);
34 
35  referenceElement.reset(new ReferenceElement(id, order));
36  if (reduced_order<0)
37  reduced_order=std::max(2*(bf_info->numOrder-1), 0);
39  reduced_order));
40 
41  if (referenceElement->getNumNodes() != referenceElementReducedQuadrature->getNumNodes()) {
42  throw escript::ValueError("ReferenceElementSet: numNodes in referenceElement and referenceElementReducedQuadrature don't match.");
43  }
44  }
45 
46  const_ShapeFunction_ptr borrowBasisFunctions(bool reducedShapefunction,
47  bool reducedIntegrationOrder) const
48  {
49  if (reducedShapefunction) {
50  return (reducedIntegrationOrder ?
51  referenceElementReducedQuadrature->LinearBasisFunctions :
52  referenceElement->LinearBasisFunctions);
53  }
54  return (reducedIntegrationOrder ?
55  referenceElementReducedQuadrature->BasisFunctions :
56  referenceElement->BasisFunctions);
57  }
58 
59  const_ShapeFunction_ptr borrowParametrization(bool reducedIntegrationOrder) const
60  {
61  return (reducedIntegrationOrder ?
62  referenceElementReducedQuadrature->Parametrization :
63  referenceElement->Parametrization);
64  }
65 
67  {
68  return (reducedIntOrder ? referenceElementReducedQuadrature :
70  }
71 
72  inline int getNumNodes() const { return referenceElement->getNumNodes(); }
73 
76 };
77 
78 
79 typedef boost::shared_ptr<const ReferenceElementSet> const_ReferenceElementSet_ptr;
80 
81 
82 } // namespace finley
83 
84 #endif // __FINLEY_REFERENCEELEMENTSETS_H__
85 
static const ShapeFunctionInfo * getInfo(ShapeFunctionTypeId id)
Definition: ShapeFunctions.cpp:104
boost::shared_ptr< const ReferenceElementSet > const_ReferenceElementSet_ptr
Definition: ReferenceElementSets.h:79
const_ReferenceElement_ptr borrowReferenceElement(bool reducedIntOrder) const
Definition: ReferenceElementSets.h:66
static const ReferenceElementInfo * getInfo(ElementTypeId id)
returns the element information structure for the given type id
Definition: ReferenceElements.cpp:663
const_ShapeFunction_ptr borrowBasisFunctions(bool reducedShapefunction, bool reducedIntegrationOrder) const
Definition: ReferenceElementSets.h:46
ReferenceElement_ptr referenceElement
Definition: ReferenceElementSets.h:75
Definition: ReferenceElementSets.h:26
ReferenceElementSet(ElementTypeId id, int order, int reduced_order)
Definition: ReferenceElementSets.h:27
A suite of factory methods for creating various finley domains.
Definition: finley/src/Assemble.h:31
ShapeFunctionTypeId BasisFunctions
shape function for the basis functions
Definition: ReferenceElements.h:145
boost::shared_ptr< const ReferenceElement > const_ReferenceElement_ptr
Definition: ReferenceElements.h:213
const_ShapeFunction_ptr borrowParametrization(bool reducedIntegrationOrder) const
Definition: ReferenceElementSets.h:59
boost::shared_ptr< const ShapeFunction > const_ShapeFunction_ptr
Definition: ShapeFunctions.h:99
this struct holds the definition of the reference element
Definition: ReferenceElements.h:120
An exception class that signals an invalid argument value.
Definition: EsysException.h:88
ElementTypeId
Definition: ReferenceElements.h:38
boost::shared_ptr< ReferenceElement > ReferenceElement_ptr
Definition: ReferenceElements.h:212
int getNumNodes() const
Definition: ReferenceElementSets.h:72
this struct holds the definition of the shape functions on an element
Definition: ShapeFunctions.h:57
this struct holds the realization of a reference element
Definition: ReferenceElements.h:176
ReferenceElement_ptr referenceElementReducedQuadrature
Definition: ReferenceElementSets.h:74
int numOrder
order of the shape functions
Definition: ShapeFunctions.h:67