ug4
piecewise_constant.h
Go to the documentation of this file.
1 /*
2  * Copyright (c) 2012-2015: G-CSC, Goethe University Frankfurt
3  * Author: Christian Wehner
4  *
5  * This file is part of UG4.
6  *
7  * UG4 is free software: you can redistribute it and/or modify it under the
8  * terms of the GNU Lesser General Public License version 3 (as published by the
9  * Free Software Foundation) with the following additional attribution
10  * requirements (according to LGPL/GPL v3 §7):
11  *
12  * (1) The following notice must be displayed in the Appropriate Legal Notices
13  * of covered and combined works: "Based on UG4 (www.ug4.org/license)".
14  *
15  * (2) The following notice must be displayed at a prominent place in the
16  * terminal output of covered works: "Based on UG4 (www.ug4.org/license)".
17  *
18  * (3) The following bibliography is recommended for citation and must be
19  * preserved in all covered files:
20  * "Reiter, S., Vogel, A., Heppner, I., Rupp, M., and Wittum, G. A massively
21  * parallel geometric multigrid solver on hierarchically distributed grids.
22  * Computing and visualization in science 16, 4 (2013), 151-164"
23  * "Vogel, A., Reiter, S., Rupp, M., Nägel, A., and Wittum, G. UG4 -- a novel
24  * flexible software system for simulating pde based models on high performance
25  * computers. Computing and visualization in science 16, 4 (2013), 165-179"
26  *
27  * This program is distributed in the hope that it will be useful,
28  * but WITHOUT ANY WARRANTY; without even the implied warranty of
29  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
30  * GNU Lesser General Public License for more details.
31  */
32 
33 #ifndef __H__UG__LIB_DISC__LOCAL_SHAPE_FUNCTION_SET__PIECEWISE_CONSTANT__
34 #define __H__UG__LIB_DISC__LOCAL_SHAPE_FUNCTION_SET__PIECEWISE_CONSTANT__
35 
36 #include "common/util/provider.h"
40 
41 namespace ug{
42 
44 template <typename TRefElem>
46  : public BaseLSFS<PiecewiseConstantLSFS<TRefElem>, TRefElem::dim>
47 {
48  public:
50  static const int dim = TRefElem::dim;
51 
52  public:
55  {
56  const TRefElem& rRef = Provider<TRefElem>::get();
57 
58  bary = rRef.corner(0);
59  for (size_t j=1; j < rRef.num(0); ++j){
60  bary += rRef.corner(j);
61  }
62  bary *= 1./rRef.num(0);
63 
64  m_vLocalDoF = LocalDoF(dim, 0, 0);
65  }
66 
67  public:
69  ReferenceObjectID roid() const {return TRefElem::REFERENCE_OBJECT_ID;}
70 
72  size_t num_dof() const {return 1;};
73 
75  size_t num_dof(ReferenceObjectID type) const
76  {
77  if (type == TRefElem::REFERENCE_OBJECT_ID) return 1;
78  else return 0;
79  }
80 
82  const LocalDoF& local_dof(size_t dof) const {return m_vLocalDoF;}
83 
85  inline bool position(size_t i, MathVector<dim>& pos) const
86  {
87  pos = bary; return true;
88  }
89 
91  bool exact_position_available() const {return true;};
92 
93  public:
95  bool continuous() const {return false;}
96 
98  size_t num_sh() const {return 1;}
99 
101  number shape(const size_t i, const MathVector<dim>& x) const
102  {
103  return 1;
104  }
105 
107  void grad(MathVector<dim>& g, const size_t i, const MathVector<dim>& x) const
108  {
109  TRefElem::check_position(x);
110  VecSet(g, 0.0);
111  }
112 
113  protected:
116 };
117 
118 } //namespace ug
119 
120 #endif /* __H__UG__LIB_DISC__LOCAL_SHAPE_FUNCTION_SET__PIECEWISE_CONSTANT__ */
121 
static interface for trial spaces
Definition: local_shape_function_set.h:176
Definition: local_dof_set.h:58
Elementwise constant shape functions.
Definition: piecewise_constant.h:47
const LocalDoF & local_dof(size_t dof) const
returns the DoFs storage
Definition: piecewise_constant.h:82
PiecewiseConstantLSFS()
Constructor.
Definition: piecewise_constant.h:54
ReferenceObjectID roid() const
returns the Reference object id of the corresponding grid object
Definition: piecewise_constant.h:69
size_t num_dof(ReferenceObjectID type) const
Definition: piecewise_constant.h:75
bool continuous() const
returns if space constructs continuous functions
Definition: piecewise_constant.h:95
number shape(const size_t i, const MathVector< dim > &x) const
evaluates the shape function
Definition: piecewise_constant.h:101
size_t num_dof() const
Definition: piecewise_constant.h:72
bool position(size_t i, MathVector< dim > &pos) const
local position of DoF i
Definition: piecewise_constant.h:85
bool exact_position_available() const
returns if the local dof position are exact
Definition: piecewise_constant.h:91
size_t num_sh() const
Definition: piecewise_constant.h:98
void grad(MathVector< dim > &g, const size_t i, const MathVector< dim > &x) const
evaluates the shape function
Definition: piecewise_constant.h:107
static const int dim
Dimension, where shape functions are defined.
Definition: piecewise_constant.h:50
MathVector< dim > bary
Barycenter.
Definition: piecewise_constant.h:114
LocalDoF m_vLocalDoF
association to elements
Definition: piecewise_constant.h:115
static TClass & get()
returns a singleton based on the identifier
Definition: provider.h:54
static const int dim
double number
Definition: types.h:124
void VecSet(vector_t &vInOut, typename vector_t::value_type s)
Set each vector component to scalar (componentwise)
Definition: math_vector_functions_common_impl.hpp:539
the ug namespace
ReferenceObjectID
these ids are used to identify the shape of a geometric object.
Definition: grid_base_objects.h:74