ug4
data_import_impl.h
Go to the documentation of this file.
1 /*
2  * Copyright (c) 2011-2015: G-CSC, Goethe University Frankfurt
3  * Author: 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__SPATIAL_DISC__DATA_IMPORT_IMPL__
34 #define __H__UG__LIB_DISC__SPATIAL_DISC__DATA_IMPORT_IMPL__
35 
36 #include "data_import.h"
37 
38 namespace ug{
39 
41 // DataImport
43 
44 template <typename TData, int dim>
46 {
47  if(data_given()) m_spUserData->unregister_storage_callback(this);
48 }
49 
50 template <typename TData, int dim>
52 {
53  if(id == ROID_UNKNOWN)
54  UG_THROW("DataImport::set_roid: Setting unknown ReferenceObjectId.");
55 
56  m_id = id;
57 }
58 
59 template <typename TData, int dim>
61 {
62 // if lin defect is not supposed to be computed, we're done
63  if(!this->m_bCompLinDefect) return;
64 
65  if(m_id == ROID_UNKNOWN)
66  UG_THROW("DataImport::check_setup: The reference element "
67  "type has not been set for evaluation.");
68 
69 // Check for evaluation function and choose it if present
70  if(m_vLinDefectFunc[m_id] != NULL)
71  return;
72 
73 // fails
74  UG_THROW("DataImport::check_setup: No evaluation function for computation of "
75  "linearized defect registered for "<<m_id<<", but required. "
76  "(world dim: "<<dim<<", part: "<<this->part()<<")");
77 }
78 
79 template <typename TData, int dim>
80 template <typename TClass>
81 void
83 set_fct(ReferenceObjectID id, TClass* obj,
84  void (TClass::*func)(const LocalVector& u,
85  std::vector<std::vector<TData> > vvvLinDefect[],
86  const size_t nip))
87 {
88  if(id >= NUM_REFERENCE_OBJECTS)
89  UG_THROW("Reference Object id invalid: "<<id);
90 
91  m_vLinDefectFunc[id] = boost::bind(func, obj, _1, _2, _3);
92 }
93 
94 template <typename TData, int dim>
95 void
98  void (*func)(const LocalVector& u,
99  std::vector<std::vector<TData> > vvvLinDefect[],
100  const size_t nip))
101 {
102  if(id >= NUM_REFERENCE_OBJECTS)
103  UG_THROW("Reference Object id invalid: "<<id);
104 
105  m_vLinDefectFunc[id] = func;
106 }
107 
108 template <typename TData, int dim>
110 {
111  for(size_t i = 0; i < NUM_REFERENCE_OBJECTS; ++i)
112  m_vLinDefectFunc[i] = NULL;
113 }
114 
115 
116 template <typename TData, int dim>
118 {
119 // remember UserData
120  m_spUserData = spData;
121 
122 // remember iexport
123  this->m_spICplUserData = spData;
124 
125 // remember dependent data (i.e. is NULL iff no dependent data given)
126  m_spDependentUserData = m_spUserData.template cast_dynamic<DependentUserData<TData, dim> >();
127 }
128 
129 template <typename TData, int dim>
131 {
132  // cache the pointer to the data field.
133  m_vValue = m_spUserData->values(m_seriesID);
134 
135  // in addition we cache the number of ips
136  m_numIP = m_spUserData->num_ip(m_seriesID);
137 }
138 
139 template <typename TData, int dim>
140 template <int ldim>
141 void DataImport<TData,dim>::set_local_ips(const MathVector<ldim>* vPos, size_t numIP, int timePointSpec,
142  bool bMayChange)
143 {
144 // if no data set, skip
145  if(!data_given()) return;
146 
147 // request series if first time requested
148  if(m_seriesID == -1)
149  {
150  m_seriesID = m_spUserData->template
151  register_local_ip_series<ldim>(vPos,numIP,timePointSpec,bMayChange);
152 
153  // register callback, invoked when data field is changed
154  m_spUserData->register_storage_callback(this, &DataImport<TData,dim>::cache_data_access);
155 
156  // cache access to the data
157  cache_data_access();
158 
159  // resize also lin defect array
160  resize_defect_array();
161 
162  // check that num ip is correct
163  UG_ASSERT(m_numIP == numIP, "Different number of ips than requested.");
164  }
165  else
166  {
167  if(!bMayChange)
168  UG_THROW("DataImport: Setting different local ips to non-changable ip series.");
169 
170  // set new local ips
171  m_spUserData->template set_local_ips<ldim>(m_seriesID,vPos,numIP);
172  m_spUserData->set_time_point(m_seriesID,timePointSpec);
173 
174  if(numIP != m_numIP)
175  {
176  // cache access to the data
177  cache_data_access();
178 
179  // resize also lin defect array
180  resize_defect_array();
181  }
182 
183  // check that num ip is correct
184  UG_ASSERT(m_numIP == numIP, "Different number of ips than requested.");
185  }
186 }
187 
188 template <typename TData, int dim>
189 void DataImport<TData,dim>::set_local_ips(const MathVector<dim>* vPos, size_t numIP, int timePointSpec,
190  bool bMayChange)
191 {
192  set_local_ips<dim>(vPos, numIP, timePointSpec, bMayChange);
193 }
194 
195 template <typename TData, int dim>
196 template <int ldim>
198  bool bMayChange)
199 {
200  set_local_ips<ldim>(vPos, numIP, -1, bMayChange);
201 }
202 
203 template <typename TData, int dim>
205  bool bMayChange)
206 {
207  set_local_ips<dim>(vPos, numIP, -1, bMayChange);
208 }
209 
210 template <typename TData, int dim>
212 {
213  m_spUserData->set_time_point(m_seriesID,timePointSpec);
214 }
215 
216 template <typename TData, int dim>
218 {
219 // if no data set, skip
220  if(!data_given()) return;
221 
222 // set global ips for series ID
223  UG_ASSERT(m_seriesID >= 0, "Wrong series id.");
224  m_spUserData->set_global_ips(m_seriesID,vPos,numIP);
225 }
226 
227 template <typename TData, int dim>
229 {
230  if(data_given()) m_spUserData->unregister_storage_callback(this);
231  m_seriesID = -1;
232  m_vValue = 0;
233  m_numIP = 0;
234  m_vvvLinDefect.resize(num_ip());
235 }
236 
237 template <typename TData, int dim>
239 {
240  UG_ASSERT(m_spDependentUserData.valid(), "No Export set.");
241 
243  if (! m_spUserData->at_current_time (m_seriesID))
244  return;
245 
246 // access jacobian by maps
247  J.access_by_map(this->map(), this->conn_map());
248 
249 // loop integration points
250  for(size_t ip = 0; ip < num_ip(); ++ip)
251  {
252  // loop all functions
253  for(size_t fct1 = 0; fct1 < this->num_fct(); ++fct1)
254  for(size_t fct2 = 0; fct2 < m_spDependentUserData->num_fct(); ++fct2)
255  {
256  // get array of linearized defect and derivative
257  const TData* LinDef = lin_defect(ip, fct1);
258  const TData* Deriv = m_spDependentUserData->deriv(m_seriesID, ip, fct2);
259 
260  // loop shapes of functions
261  for(size_t sh1 = 0; sh1 < num_sh(fct1); ++sh1)
262  for(size_t sh2 = 0; sh2 < m_spDependentUserData->num_sh(fct2); ++sh2)
263  {
264  J(fct1, sh1, fct2, sh2) += scale*(LinDef[sh1]*Deriv[sh2]);
265  }
266  }
267  }
268 }
269 
270 template <typename TData, int dim>
272 {
273 // check size
274  const FunctionIndexMapping& map = this->map();
275  UG_ASSERT(map.num_fct() == this->num_fct(), "Number function mismatch.");
276 
277 // cache numFct and their numDoFs
278  m_vvNumDoFPerFct.resize(map.num_fct());
279  for(size_t fct = 0; fct < m_vvNumDoFPerFct.size(); ++fct)
280  m_vvNumDoFPerFct[fct] = ind.num_dof(map[fct]);
281 
282  m_vvvLinDefect.clear();
283  resize_defect_array();
284 }
285 
286 template <typename TData, int dim>
288 {
289 // get old size
290 // NOTE: for all ips up to oldSize the arrays are already resized
291  const size_t oldSize = m_vvvLinDefect.size();
292 
293 // resize ips
294  m_vvvLinDefect.resize(num_ip());
295 
296 // resize num fct
297  for(size_t ip = oldSize; ip < num_ip(); ++ip)
298  {
299  // resize num fct
300  m_vvvLinDefect[ip].resize(m_vvNumDoFPerFct.size());
301 
302  // resize dofs
303  for(size_t fct = 0; fct < m_vvNumDoFPerFct.size(); ++fct)
304  m_vvvLinDefect[ip][fct].resize(m_vvNumDoFPerFct[fct]);
305  }
306 }
307 
308 template <typename TData, int dim>
309 inline void DataImport<TData,dim>::check_ip_fct(size_t ip, size_t fct) const
310 {
311  check_ip(ip);
312  UG_ASSERT(ip < m_vvvLinDefect.size(), "Invalid index.");
313  UG_ASSERT(fct < m_vvvLinDefect[ip].size(), "Invalid index.");
314 }
315 
316 template <typename TData, int dim>
317 inline void DataImport<TData,dim>::check_ip_fct_sh(size_t ip, size_t fct, size_t sh) const
318 {
319  check_ip_fct(ip, fct);
320  UG_ASSERT(sh < m_vvvLinDefect[ip][fct].size(), "Invalid index.");
321 }
322 
323 template <typename TData, int dim>
324 inline void DataImport<TData,dim>::check_ip(size_t ip) const
325 {
326  UG_ASSERT(ip < m_numIP, "Invalid index.");
327 }
328 
329 template <typename TData, int dim>
331 {
332  UG_ASSERT(m_vValue != NULL, "Data Value field not set.");
333 }
334 
335 } // end namespace ug
336 
337 #endif /* __H__UG__LIB_DISC__SPATIAL_DISC__DATA_IMPORT_IMPL__ */
Definition: smart_pointer.h:108
Type based UserData.
Definition: user_data.h:501
Data import.
Definition: data_import.h:180
void add_jacobian(LocalMatrix &J, const number scale)
compute jacobian for derivative w.r.t. non-system owned unknowns
Definition: data_import_impl.h:238
virtual void clear_ips()
removes the positions
Definition: data_import_impl.h:228
void check_ip(size_t ip) const
checks in debug mode the correct index
Definition: data_import_impl.h:324
virtual void update_dof_sizes(const LocalIndices &ind)
resize lin defect arrays
Definition: data_import_impl.h:271
void clear_fct()
clear all evaluation functions
Definition: data_import_impl.h:109
virtual void set_roid(ReferenceObjectID id)
sets the geometric object type
Definition: data_import_impl.h:51
void set_data(SmartPtr< CplUserData< TData, dim > > spData)
set the user data
Definition: data_import_impl.h:117
~DataImport()
Destructor.
Definition: data_import_impl.h:45
void set_global_ips(const MathVector< dim > *vPos, size_t numIP)
sets the global positions
Definition: data_import_impl.h:217
void cache_data_access()
caches data access
Definition: data_import_impl.h:130
void resize_defect_array()
resizes the lin defect arrays for current number of ips.
Definition: data_import_impl.h:287
void check_ip_fct_sh(size_t ip, size_t fct, size_t sh) const
checks in debug mode the correct index
Definition: data_import_impl.h:317
void set_fct(ReferenceObjectID id, TClass *obj, void(TClass::*func)(const LocalVector &u, std::vector< std::vector< TData > > vvvLinDefect[], const size_t nip))
register evaluation of linear defect for a element
Definition: data_import_impl.h:83
virtual void check_setup()
checks if ready for evaluation
Definition: data_import_impl.h:60
void set_local_ips(const MathVector< ldim > *vPos, size_t numIP, int timePointSpec, bool bMayChange=true)
set the local integration points
Definition: data_import_impl.h:141
void check_values() const
checks in debug mode the correct index
Definition: data_import_impl.h:330
void check_ip_fct(size_t ip, size_t fct) const
checks in debug mode the correct index
Definition: data_import_impl.h:309
void set_time_point(int timePointSpec)
set the time point specification
Definition: data_import_impl.h:211
describes a mapping between two local index sets
Definition: function_group.h:186
size_t num_fct() const
returns the number of indices that are mapped
Definition: function_group.h:195
Definition: local_algebra.h:50
size_t num_dof(size_t fct) const
number of dofs for accessible function
Definition: local_algebra.h:117
Definition: local_algebra.h:422
void access_by_map(const FunctionIndexMapping &funcMap)
access only part of the functions using mapping (restrict functions)
Definition: local_algebra.h:538
Definition: local_algebra.h:198
a mathematical Vector with N entries.
Definition: math_vector.h:97
static const int dim
#define UG_ASSERT(expr, msg)
Definition: assert.h:70
#define UG_THROW(msg)
Definition: error.h:57
double number
Definition: types.h:124
the ug namespace
ReferenceObjectID
these ids are used to identify the shape of a geometric object.
Definition: grid_base_objects.h:74
@ NUM_REFERENCE_OBJECTS
Definition: grid_base_objects.h:85
@ ROID_UNKNOWN
Definition: grid_base_objects.h:75
function func(x, y, z, t, si)
bool resize(size_t newRows, size_t newCols)