ug4
Loading...
Searching...
No Matches
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
38namespace ug{
39
41// DataImport
43
44template <typename TData, int dim>
46{
47 if(data_given()) m_spUserData->unregister_storage_callback(this);
48}
49
50template <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
59template <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
79template <typename TData, int dim>
80template <typename TClass>
81void
83set_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] = std::bind(func, obj,
92 std::placeholders::_1, std::placeholders::_2, std::placeholders::_3);
93}
94
95template <typename TData, int dim>
96void
99 void (*func)(const LocalVector& u,
100 std::vector<std::vector<TData> > vvvLinDefect[],
101 const size_t nip))
102{
103 if(id >= NUM_REFERENCE_OBJECTS)
104 UG_THROW("Reference Object id invalid: "<<id);
105
106 m_vLinDefectFunc[id] = func;
107}
108
109template <typename TData, int dim>
111{
112 for(size_t i = 0; i < NUM_REFERENCE_OBJECTS; ++i)
113 m_vLinDefectFunc[i] = NULL;
114}
115
116
117template <typename TData, int dim>
119{
120// remember UserData
121 m_spUserData = spData;
122
123// remember iexport
124 this->m_spICplUserData = spData;
125
126// remember dependent data (i.e. is NULL iff no dependent data given)
127 m_spDependentUserData = m_spUserData.template cast_dynamic<DependentUserData<TData, dim> >();
128}
129
130template <typename TData, int dim>
132{
133 // cache the pointer to the data field.
134 m_vValue = m_spUserData->values(m_seriesID);
135
136 // in addition we cache the number of ips
137 m_numIP = m_spUserData->num_ip(m_seriesID);
138}
139
140template <typename TData, int dim>
141template <int ldim>
142void DataImport<TData,dim>::set_local_ips(const MathVector<ldim>* vPos, size_t numIP, int timePointSpec,
143 bool bMayChange)
144{
145// if no data set, skip
146 if(!data_given()) return;
147
148// request series if first time requested
149 if(m_seriesID == -1)
150 {
151 m_seriesID = m_spUserData->template
152 register_local_ip_series<ldim>(vPos,numIP,timePointSpec,bMayChange);
153
154 // register callback, invoked when data field is changed
155 m_spUserData->register_storage_callback(this, &DataImport<TData,dim>::cache_data_access);
156
157 // cache access to the data
158 cache_data_access();
159
160 // resize also lin defect array
161 resize_defect_array();
162
163 // check that num ip is correct
164 UG_ASSERT(m_numIP == numIP, "Different number of ips than requested.");
165 }
166 else
167 {
168 if(!bMayChange)
169 UG_THROW("DataImport: Setting different local ips to non-changable ip series.");
170
171 // set new local ips
172 m_spUserData->template set_local_ips<ldim>(m_seriesID,vPos,numIP);
173 m_spUserData->set_time_point(m_seriesID,timePointSpec);
174
175 if(numIP != m_numIP)
176 {
177 // cache access to the data
178 cache_data_access();
179
180 // resize also lin defect array
181 resize_defect_array();
182 }
183
184 // check that num ip is correct
185 UG_ASSERT(m_numIP == numIP, "Different number of ips than requested.");
186 }
187}
188
189template <typename TData, int dim>
190void DataImport<TData,dim>::set_local_ips(const MathVector<dim>* vPos, size_t numIP, int timePointSpec,
191 bool bMayChange)
192{
193 set_local_ips<dim>(vPos, numIP, timePointSpec, bMayChange);
195
196template <typename TData, int dim>
197template <int ldim>
199 bool bMayChange)
200{
201 set_local_ips<ldim>(vPos, numIP, -1, bMayChange);
202}
203
204template <typename TData, int dim>
206 bool bMayChange)
207{
208 set_local_ips<dim>(vPos, numIP, -1, bMayChange);
209}
210
211template <typename TData, int dim>
213{
214 m_spUserData->set_time_point(m_seriesID,timePointSpec);
215}
216
217template <typename TData, int dim>
219{
220// if no data set, skip
221 if(!data_given()) return;
222
223// set global ips for series ID
224 UG_ASSERT(m_seriesID >= 0, "Wrong series id.");
225 m_spUserData->set_global_ips(m_seriesID,vPos,numIP);
226}
227
228template <typename TData, int dim>
230{
231 if(data_given()) m_spUserData->unregister_storage_callback(this);
232 m_seriesID = -1;
233 m_vValue = 0;
234 m_numIP = 0;
235 m_vvvLinDefect.resize(num_ip());
236}
237
238template <typename TData, int dim>
240{
241 UG_ASSERT(m_spDependentUserData.valid(), "No Export set.");
242
244 if (! m_spUserData->at_current_time (m_seriesID))
245 return;
246
247// access jacobian by maps
248 J.access_by_map(this->map(), this->conn_map());
249
250// loop integration points
251 for(size_t ip = 0; ip < num_ip(); ++ip)
253 // loop all functions
254 for(size_t fct1 = 0; fct1 < this->num_fct(); ++fct1)
255 for(size_t fct2 = 0; fct2 < m_spDependentUserData->num_fct(); ++fct2)
256 {
257 // get array of linearized defect and derivative
258 const TData* LinDef = lin_defect(ip, fct1);
259 const TData* Deriv = m_spDependentUserData->deriv(m_seriesID, ip, fct2);
260
261 // loop shapes of functions
262 for(size_t sh1 = 0; sh1 < num_sh(fct1); ++sh1)
263 for(size_t sh2 = 0; sh2 < m_spDependentUserData->num_sh(fct2); ++sh2)
264 {
265 J(fct1, sh1, fct2, sh2) += scale*(LinDef[sh1]*Deriv[sh2]);
266 }
267 }
269}
270
271template <typename TData, int dim>
273{
274// check size
275 const FunctionIndexMapping& map = this->map();
276 UG_ASSERT(map.num_fct() == this->num_fct(), "Number function mismatch.");
277
278// cache numFct and their numDoFs
279 m_vvNumDoFPerFct.resize(map.num_fct());
280 for(size_t fct = 0; fct < m_vvNumDoFPerFct.size(); ++fct)
281 m_vvNumDoFPerFct[fct] = ind.num_dof(map[fct]);
282
283 m_vvvLinDefect.clear();
284 resize_defect_array();
285}
286
287template <typename TData, int dim>
289{
290// get old size
291// NOTE: for all ips up to oldSize the arrays are already resized
292 const size_t oldSize = m_vvvLinDefect.size();
293
294// resize ips
295 m_vvvLinDefect.resize(num_ip());
296
297// resize num fct
298 for(size_t ip = oldSize; ip < num_ip(); ++ip)
299 {
300 // resize num fct
301 m_vvvLinDefect[ip].resize(m_vvNumDoFPerFct.size());
302
303 // resize dofs
304 for(size_t fct = 0; fct < m_vvNumDoFPerFct.size(); ++fct)
305 m_vvvLinDefect[ip][fct].resize(m_vvNumDoFPerFct[fct]);
306 }
307}
308
309template <typename TData, int dim>
310inline void DataImport<TData,dim>::check_ip_fct(size_t ip, size_t fct) const
311{
312 check_ip(ip);
313 UG_ASSERT(ip < m_vvvLinDefect.size(), "Invalid index.");
314 UG_ASSERT(fct < m_vvvLinDefect[ip].size(), "Invalid index.");
315}
316
317template <typename TData, int dim>
318inline void DataImport<TData,dim>::check_ip_fct_sh(size_t ip, size_t fct, size_t sh) const
319{
320 check_ip_fct(ip, fct);
321 UG_ASSERT(sh < m_vvvLinDefect[ip][fct].size(), "Invalid index.");
322}
323
324template <typename TData, int dim>
325inline void DataImport<TData,dim>::check_ip(size_t ip) const
326{
327 UG_ASSERT(ip < m_numIP, "Invalid index.");
329
330template <typename TData, int dim>
332{
333 UG_ASSERT(m_vValue != NULL, "Data Value field not set.");
334}
336} // end namespace ug
337
338#endif /* __H__UG__LIB_DISC__SPATIAL_DISC__DATA_IMPORT_IMPL__ */
Definition smart_pointer.h:107
Type based UserData.
Definition user_data.h:502
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:239
virtual void clear_ips()
removes the positions
Definition data_import_impl.h:229
void check_ip(size_t ip) const
checks in debug mode the correct index
Definition data_import_impl.h:325
virtual void update_dof_sizes(const LocalIndices &ind)
resize lin defect arrays
Definition data_import_impl.h:272
void clear_fct()
clear all evaluation functions
Definition data_import_impl.h:110
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:118
~DataImport()
Destructor.
Definition data_import_impl.h:45
LinDefectFunc m_vLinDefectFunc[NUM_REFERENCE_OBJECTS]
function pointers for all elem types
Definition data_import.h:381
void set_global_ips(const MathVector< dim > *vPos, size_t numIP)
sets the global positions
Definition data_import_impl.h:218
void cache_data_access()
caches data access
Definition data_import_impl.h:131
void resize_defect_array()
resizes the lin defect arrays for current number of ips.
Definition data_import_impl.h:288
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:318
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:142
void check_values() const
checks in debug mode the correct index
Definition data_import_impl.h:331
void check_ip_fct(size_t ip, size_t fct) const
checks in debug mode the correct index
Definition data_import_impl.h:310
void set_time_point(int timePointSpec)
set the time point specification
Definition data_import_impl.h:212
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
#define UG_ASSERT(expr, msg)
Definition assert.h:70
#define UG_THROW(msg)
Definition error.h:57
double number
Definition types.h:124
Definition smart_pointer.h:817
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)