ug4
grid_function_global_user_data.h
Go to the documentation of this file.
1 /*
2  * Copyright (c) 2013-2015: G-CSC, Goethe University Frankfurt
3  * Author: Sebastian Reiter, Andreas Vogel
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__FUNCTION_SPACE__GRID_FUNCTION_GLOBAL_USER_DATA__
34 #define __H__UG__LIB_DISC__FUNCTION_SPACE__GRID_FUNCTION_GLOBAL_USER_DATA__
35 
36 #include "common/common.h"
37 
39 
43 #include "lib_disc/domain_util.h"
48 
49 #include <math.h> /* fabs */
50 
51 namespace ug{
52 
53 
54 
55 template <typename TGridFunction, int elemDim = TGridFunction::dim>
57  : public StdGlobPosData<GlobalGridFunctionNumberData<TGridFunction, elemDim>, number, TGridFunction::dim>
58 {
59  public:
61  static const int dim = TGridFunction::dim;
62  typedef typename TGridFunction::template dim_traits<elemDim>::grid_base_object element_t;
63 
64  private:
67 
69  size_t m_fct;
70 
73 
76 
77  public:
80  : m_spGridFct(spGridFct),
81  m_tree(*spGridFct->domain()->grid(), spGridFct->domain()->position_attachment())
82  {
83  //this->set_functions(cmp);
84 
85  // get function id of name
86  m_fct = spGridFct->fct_id_by_name(cmp);
87 
88  // check that function exists
89  if(m_fct >= spGridFct->num_fct())
90  UG_THROW("GridFunctionNumberData: Function space does not contain"
91  " a function with name " << cmp << ".");
92 
93  // local finite element id
94  m_lfeID = spGridFct->local_finite_element_id(m_fct);
95 
96 
97 
98 // const MGSubsetHandler& ssh = *m_spGridFct->domain()->subset_handler();
99 
100  SubsetGroup ssGrp(m_spGridFct->domain()->subset_handler());
101  ssGrp.add_all();
102 
103  std::vector<size_t> subsetsOfGridFunction;
104  std::vector<element_t*> elemsWithGridFunctions;
105 
106  typename TGridFunction::template dim_traits<elemDim>::const_iterator iterEnd, iter;
107 
108 
109  for(size_t si = 0; si < ssGrp.size(); si++){
110  if( spGridFct->is_def_in_subset(m_fct, si) )
111  subsetsOfGridFunction.push_back(si);
112  }
113 
114 
115  for(size_t i = 0; i<subsetsOfGridFunction.size(); i++){
116  size_t si = subsetsOfGridFunction[i];
117  iter = spGridFct->template begin<element_t>(si);
118  iterEnd = spGridFct->template end<element_t>(si);
119 
120  for(;iter!=iterEnd; ++iter){
121  element_t *elem = *iter;
122  elemsWithGridFunctions.push_back(elem);
123  }
124  }
125 
126 // m_tree.create_tree(spGridFct->template begin<element_t>(si),
127 // spGridFct->template end<element_t>(si));
128 
129  m_tree.create_tree(elemsWithGridFunctions.begin(), elemsWithGridFunctions.end());
130 
131  };
132 
134 
135  virtual bool continuous() const
136  {
138  }
139 
141  inline void evaluate(number& value, const MathVector<dim>& x, number time, int si) const
142  {
143  if(!evaluate(value, x))
144  UG_THROW("For function "<<m_fct<<" couldn't find an element containing the specified point: " << x);
145  }
146 
147  inline number evaluate(const MathVector<dim>& x) const
148  {
149  number value;
150  if(!evaluate(value, x))
151  UG_THROW("For function "<<m_fct<<" couldn't find an element containing the specified point: " << x);
152  return value;
153  }
154 
156  inline bool evaluate(number& value, const MathVector<dim>& x) const
157  {
158  element_t* elem = NULL;
159  //try{
160 
161  if(!FindContainingElement(elem, m_tree, x)){
162  return false;
163  }
164 
165  // get corners of element
166  std::vector<MathVector<dim> > vCornerCoords;
167  CollectCornerCoordinates(vCornerCoords, *elem, *m_spGridFct->domain());
168 
169  // reference object id
170  const ReferenceObjectID roid = elem->reference_object_id();
171 
172  // get local position of DoF
174  = ReferenceMappingProvider::get<elemDim, dim>(roid, vCornerCoords);
175  MathVector<elemDim> locPos;
176  VecSet(locPos, 0.5);
177  map.global_to_local(locPos, x);
178 
179  // evaluate at shapes at ip
180  const LocalShapeFunctionSet<elemDim>& rTrialSpace =
181  LocalFiniteElementProvider::get<elemDim>(roid, m_lfeID);
182  std::vector<number> vShape;
183  rTrialSpace.shapes(vShape, locPos);
184 
185  // get multiindices of element
186  std::vector<DoFIndex> ind;
187  m_spGridFct->dof_indices(elem, m_fct, ind);
188 
189  // compute solution at integration point
190  value = 0.0;
191  for(size_t sh = 0; sh < vShape.size(); ++sh)
192  {
193  const number valSH = DoFRef(*m_spGridFct, ind[sh]);
194  value += valSH * vShape[sh];
195  }
196 
197  // point is found
198  return true;
199  //}
200  //UG_CATCH_THROW("GlobalGridFunctionNumberData: Evaluation failed."
201  // << "Point: " << x << ", Element: "
202  // << ElementDebugInfo(*m_spGridFct->domain()->grid(), elem));
203  }
204 
206  inline void evaluate_global(number& value, const MathVector<dim>& x) const
207  {
208  // evaluate at this proc
209  bool bFound = this->evaluate(value, x);
210 
211 #ifdef UG_PARALLEL
212  // share value between all procs
214  int numFound = (bFound ? 1 : 0);
215  numFound = com.allreduce(numFound, PCL_RO_SUM);
216 
217  // get overall value
218  if(!bFound) value = 0.0;
219  number globValue = com.allreduce(value, PCL_RO_SUM) / numFound;
220 
221  if(numFound == 0)
222  UG_THROW("Point "<<x<<" not found on all "<<pcl::NumProcs()<<" procs.");
223 
224  // check correctness for continuous spaces
225  // note: if the point is found more than one it is located on the
226  // boundary of some element. thus, if the space is continuous, those
227  // values should match on all procs.
228  if(bFound)
230  if( fabs(value) > 1e-10 && fabs((globValue - value) / value) > 1e-8)
231  UG_THROW("Global mean "<<globValue<<" != local value "<<value);
232  }
233 
234  // set as global value
235  value = globValue;
236  bFound = true;
237 #endif
238 
239  if(!bFound)
240  UG_THROW("Couldn't find an element containing the specified point: " << x);
241  }
242 
243  // evaluates at given position
244  number evaluate_global(std::vector<number> vPos)
245  {
246  if((int)vPos.size() != dim)
247  UG_THROW("Expected "<<dim<<" components, but given "<<vPos.size());
248 
249  MathVector<dim> x;
250  for(int i = 0; i < dim; i++) x[i] = vPos[i];
251 
252  number value;
254 
255  return value;
256  }
257 };
258 
259 
260 
261 template <typename TGridFunction>
263  : public StdGlobPosData<GlobalGridFunctionGradientData<TGridFunction>, MathVector<TGridFunction::dim>, TGridFunction::dim>
264 {
265  public:
267  static const int dim = TGridFunction::dim;
268  typedef typename TGridFunction::element_type element_t;
269 
270  private:
273 
275  size_t m_fct;
276 
279 
282 
283  public:
286  : m_spGridFct(spGridFct),
287  m_tree(*spGridFct->domain()->grid(), spGridFct->domain()->position_attachment())
288  {
289  //this->set_functions(cmp);
290 
291  // get function id of name
292  m_fct = spGridFct->fct_id_by_name(cmp);
293 
294  // check that function exists
295  if(m_fct >= spGridFct->num_fct())
296  UG_THROW("GridFunctionNumberData: Function space does not contain"
297  " a function with name " << cmp << ".");
298 
299  // local finite element id
300  m_lfeID = spGridFct->local_finite_element_id(m_fct);
301 
302 
303 
304 // const MGSubsetHandler& ssh = *m_spGridFct->domain()->subset_handler();
305 
306  SubsetGroup ssGrp(m_spGridFct->domain()->subset_handler());
307  ssGrp.add_all();
308 
309  std::vector<size_t> subsetsOfGridFunction;
310  std::vector<element_t*> elemsWithGridFunctions;
311 
312  typename TGridFunction::const_element_iterator iterEnd, iter;
313 
314 
315  for(size_t si = 0; si < ssGrp.size(); si++){
316  if( spGridFct->is_def_in_subset(m_fct, si) )
317  subsetsOfGridFunction.push_back(si);
318  }
319 
320 
321  for(size_t i = 0; i<subsetsOfGridFunction.size(); i++){
322  size_t si = subsetsOfGridFunction[i];
323  iter = spGridFct->template begin<element_t>(si);
324  iterEnd = spGridFct->template end<element_t>(si);
325 
326  for(;iter!=iterEnd; ++iter){
327  element_t *elem = *iter;
328  elemsWithGridFunctions.push_back(elem);
329  }
330  }
331 
332 // m_tree.create_tree(spGridFct->template begin<element_t>(si),
333 // spGridFct->template end<element_t>(si));
334 
335  m_tree.create_tree(elemsWithGridFunctions.begin(), elemsWithGridFunctions.end());
336 
337  };
338 
340 
341  virtual bool continuous() const
342  {
343  return false;
344  }
345 
347  inline void evaluate(MathVector<dim>& value, const MathVector<dim>& x, number time, int si) const
348  {
349  if(!evaluate(value, x))
350  UG_THROW("For function "<<m_fct<<" couldn't find an element containing the specified point: " << x);
351  }
352 
354  inline bool evaluate(MathVector<dim>& value, const MathVector<dim>& x) const
355  {
356  static const int refDim = dim;
357 
358  element_t* elem = NULL;
359  try{
360 
361  if(!FindContainingElement(elem, m_tree, x)){
362  return false;
363  }
364 
365  // get corners of element
366  std::vector<MathVector<dim> > vCornerCoords;
367  CollectCornerCoordinates(vCornerCoords, *elem, *m_spGridFct->domain());
368 
369  // reference object id
370  const ReferenceObjectID roid = elem->reference_object_id();
371 
372  // get local position of DoF
374  = ReferenceMappingProvider::get<dim, dim>(roid, vCornerCoords);
375  MathVector<dim> locPos;
376  VecSet(locPos, 0.5);
377  map.global_to_local(locPos, x);
378 
379 
381  try{
383  = ReferenceMappingProvider::get<refDim, dim>(roid, vCornerCoords);
384 
385  // compute transformation matrices
386  mapping.jacobian_transposed(JT, x);
387 
388  }UG_CATCH_THROW("GlobalGridFunctionGradientData: failed.");
389 
390  // evaluate at shapes at ip
391  const LocalShapeFunctionSet<dim>& rTrialSpace =
392  LocalFiniteElementProvider::get<dim>(roid, m_lfeID);
393  std::vector<MathVector<refDim> > vLocGrad;
394  rTrialSpace.grads(vLocGrad, locPos);
395 
396 
397  // Reference Mapping
399  RightInverse (JTInv, JT);
400 
401  // get multiindices of element
402  std::vector<DoFIndex> ind;
403  m_spGridFct->dof_indices(elem, m_fct, ind);
404 
405  // storage for shape function at ip
406  MathVector<refDim> locGrad;
407 
408  // compute grad at ip
409  VecSet(locGrad, 0.0);
410  for(size_t sh = 0; sh < vLocGrad.size(); ++sh)
411  {
412  const number valSH = DoFRef( *m_spGridFct, ind[sh]);
413  VecScaleAppend(locGrad, valSH, vLocGrad[sh]);
414  }
415 
416  // transform to global space
417  MatVecMult(value, JTInv, locGrad);
418 
419  // point is found
420  return true;
421  }
422  UG_CATCH_THROW("GlobalGridFunctionGradientData: Evaluation failed."
423  << "Point: " << x << ", Element: "
424  << ElementDebugInfo(*m_spGridFct->domain()->grid(), elem));
425  }
426 
428  inline void evaluate_global(MathVector<dim>& value, const MathVector<dim>& x) const
429  {
430  // evaluate at this proc
431  bool bFound = this->evaluate(value, x);
432 
433  // \todo: (optinal) check for evaluation on other procs
434 
435  if(!bFound)
436  UG_THROW("Couldn't find an element containing the specified point: " << x);
437  }
438 
439  // evaluates at given position
440  std::vector<number> evaluate_global(std::vector<number> vPos)
441  {
442  if((int)vPos.size() != dim)
443  UG_THROW("Expected "<<dim<<" components, but given "<<vPos.size());
444 
445  MathVector<dim> x;
446  for(int i = 0; i < dim; i++) x[i] = vPos[i];
447 
450 
451  for(int i = 0; i < dim; i++) vPos[i] = value[i];
452  return vPos;
453  }
454 };
455 
456 } // end namespace ug
457 
458 #endif /* __H__UG__LIB_DISC__FUNCTION_SPACE__GRID_FUNCTION_GLOBAL_USER_DATA__ */
Definition: pcl_process_communicator.h:70
void allreduce(const void *sendBuf, void *recBuf, int count, DataType type, ReduceOperation op) const
performs MPI_Allreduce on the processes of the communicator.
Definition: pcl_process_communicator.cpp:318
const TData & value(size_t s, size_t ip) const
returns the value at ip
Definition: user_data.h:486
virtual base class for reference mappings
Definition: reference_mapping_provider.h:53
virtual void jacobian_transposed(MathMatrix< dim, worldDim > &JT, const MathVector< dim > &locPos) const =0
returns transposed of jacobian
Definition: grid_function_global_user_data.h:264
static const int dim
world dimension of grid function
Definition: grid_function_global_user_data.h:267
virtual bool continuous() const
returns if provided data is continuous over geometric object boundaries
Definition: grid_function_global_user_data.h:341
void evaluate_global(MathVector< dim > &value, const MathVector< dim > &x) const
evaluate value on all procs
Definition: grid_function_global_user_data.h:428
size_t m_fct
component of function
Definition: grid_function_global_user_data.h:275
bool evaluate(MathVector< dim > &value, const MathVector< dim > &x) const
evaluates the data at a given point, returns false if point not found
Definition: grid_function_global_user_data.h:354
lg_ntree< dim, dim, element_t > tree_t
Definition: grid_function_global_user_data.h:280
std::vector< number > evaluate_global(std::vector< number > vPos)
Definition: grid_function_global_user_data.h:440
LFEID m_lfeID
local finite element id
Definition: grid_function_global_user_data.h:278
SmartPtr< TGridFunction > m_spGridFct
grid function
Definition: grid_function_global_user_data.h:272
GlobalGridFunctionGradientData(SmartPtr< TGridFunction > spGridFct, const char *cmp)
constructor
Definition: grid_function_global_user_data.h:285
tree_t m_tree
Definition: grid_function_global_user_data.h:281
void evaluate(MathVector< dim > &value, const MathVector< dim > &x, number time, int si) const
to full-fill UserData-Interface
Definition: grid_function_global_user_data.h:347
TGridFunction::element_type element_t
Definition: grid_function_global_user_data.h:268
virtual ~GlobalGridFunctionGradientData()
Definition: grid_function_global_user_data.h:339
Definition: grid_function_global_user_data.h:58
lg_ntree< dim, dim, element_t > tree_t
Definition: grid_function_global_user_data.h:74
static const int dim
world dimension of grid function
Definition: grid_function_global_user_data.h:61
virtual ~GlobalGridFunctionNumberData()
Definition: grid_function_global_user_data.h:133
virtual bool continuous() const
returns if provided data is continuous over geometric object boundaries
Definition: grid_function_global_user_data.h:135
bool evaluate(number &value, const MathVector< dim > &x) const
evaluates the data at a given point, returns false if point not found
Definition: grid_function_global_user_data.h:156
void evaluate(number &value, const MathVector< dim > &x, number time, int si) const
to full-fill UserData-Interface
Definition: grid_function_global_user_data.h:141
number evaluate_global(std::vector< number > vPos)
Definition: grid_function_global_user_data.h:244
GlobalGridFunctionNumberData(SmartPtr< TGridFunction > spGridFct, const char *cmp)
constructor
Definition: grid_function_global_user_data.h:79
void evaluate_global(number &value, const MathVector< dim > &x) const
evaluate value on all procs
Definition: grid_function_global_user_data.h:206
size_t m_fct
component of function
Definition: grid_function_global_user_data.h:69
tree_t m_tree
Definition: grid_function_global_user_data.h:75
SmartPtr< TGridFunction > m_spGridFct
grid function
Definition: grid_function_global_user_data.h:66
number evaluate(const MathVector< dim > &x) const
Definition: grid_function_global_user_data.h:147
TGridFunction::template dim_traits< elemDim >::grid_base_object element_t
Definition: grid_function_global_user_data.h:62
LFEID m_lfeID
local finite element id
Definition: grid_function_global_user_data.h:72
number time() const
get the current evaluation time
Definition: user_data.h:259
Identifier for Local Finite Elements.
Definition: local_finite_element_id.h:98
static bool continuous(const LFEID &id, bool bCreate=true)
returns if a Local Shape Function Set is continuous
Definition: local_finite_element_provider.cpp:749
virtual base class for local shape function sets
Definition: local_shape_function_set.h:70
virtual void grads(grad_type *vGrad, const MathVector< dim > &x) const =0
returns all gradients evaluated at a point
virtual void shapes(shape_type *vShape, const MathVector< dim > &x) const =0
returns all shape functions evaluated at a point
A class for fixed size, dense matrices.
Definition: math_matrix.h:52
Definition: std_glob_pos_data.h:55
Group of subsets.
Definition: subset_group.h:51
size_t size() const
number of subsets in this group
Definition: subset_group.h:122
void add_all()
select all subsets of underlying subset handler
Definition: subset_group.cpp:133
const FunctionIndexMapping & map() const
get function mapping
Definition: user_data.h:76
void create_tree(TIterator elemsBegin, TIterator elemsEnd)
Definition: lg_ntree.h:356
void CollectCornerCoordinates(std::vector< typename TAAPos::ValueType > &vCornerCoordsOut, const TElem &elem, const TAAPos &aaPos, bool clearContainer=true)
returns the corner coordinates of a geometric object
Definition: domain_util_impl.h:75
std::string ElementDebugInfo(const Grid &grid, GridObject *e)
Returns a string containing information on the given element.
Definition: debug_util.cpp:991
MathMatrix< N, M, T >::value_type RightInverse(MathMatrix< N, M, T > &mOut, const MathMatrix< M, N, T > &m)
Right-Inverse of a Matrix.
Definition: math_matrix_functions_common_impl.hpp:680
#define PCL_RO_SUM
Definition: pcl_methods.h:63
int NumProcs()
returns the number of processes
Definition: pcl_base.cpp:91
static const int dim
#define UG_CATCH_THROW(msg)
Definition: error.h:64
#define UG_THROW(msg)
Definition: error.h:57
double number
Definition: types.h:124
void MatVecMult(vector_t_out &vOut, const matrix_t &m, const vector_t_in &v)
Matrix - Vector Multiplication.
Definition: math_matrix_vector_functions_common_impl.hpp:49
void VecScaleAppend(vector_t &vOut, typename vector_t::value_type s1, const vector_t &v1)
Scales a Vector and adds it to a second vector.
Definition: math_vector_functions_common_impl.hpp:126
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
number & DoFRef(TMatrix &mat, const DoFIndex &iInd, const DoFIndex &jInd)
Definition: multi_index.h:276
bool FindContainingElement(typename tree_t::elem_t &elemOut, const tree_t &tree, const typename tree_t::vector_t &point)
Definition: ntree_traverser.h:210