Line data Source code
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 :
40 : ////////////////////////////////////////////////////////////////////////////////
41 : // DataImport
42 : ////////////////////////////////////////////////////////////////////////////////
43 :
44 : template <typename TData, int dim>
45 0 : DataImport<TData,dim>::~DataImport()
46 : {
47 0 : if(data_given()) m_spUserData->unregister_storage_callback(this);
48 0 : }
49 :
50 : template <typename TData, int dim>
51 0 : void DataImport<TData,dim>::set_roid(ReferenceObjectID id)
52 : {
53 0 : if(id == ROID_UNKNOWN)
54 0 : UG_THROW("DataImport::set_roid: Setting unknown ReferenceObjectId.");
55 :
56 0 : m_id = id;
57 0 : }
58 :
59 : template <typename TData, int dim>
60 0 : void DataImport<TData,dim>::check_setup()
61 : {
62 : // if lin defect is not supposed to be computed, we're done
63 0 : if(!this->m_bCompLinDefect) return;
64 :
65 0 : if(m_id == ROID_UNKNOWN)
66 0 : 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 0 : if(m_vLinDefectFunc[m_id] != NULL)
71 : return;
72 :
73 : // fails
74 0 : 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
82 0 : DataImport<TData,dim>::
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 0 : if(id >= NUM_REFERENCE_OBJECTS)
89 0 : UG_THROW("Reference Object id invalid: "<<id);
90 :
91 0 : m_vLinDefectFunc[id] = boost::bind(func, obj, _1, _2, _3);
92 0 : }
93 :
94 : template <typename TData, int dim>
95 : void
96 : DataImport<TData,dim>::
97 : set_fct(ReferenceObjectID id,
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>
109 0 : void DataImport<TData,dim>::clear_fct()
110 : {
111 0 : for(size_t i = 0; i < NUM_REFERENCE_OBJECTS; ++i)
112 0 : m_vLinDefectFunc[i] = NULL;
113 0 : }
114 :
115 :
116 : template <typename TData, int dim>
117 0 : void DataImport<TData,dim>::set_data(SmartPtr<CplUserData<TData, dim> > spData)
118 : {
119 : // remember UserData
120 0 : m_spUserData = spData;
121 :
122 : // remember iexport
123 0 : this->m_spICplUserData = spData;
124 :
125 : // remember dependent data (i.e. is NULL iff no dependent data given)
126 0 : m_spDependentUserData = m_spUserData.template cast_dynamic<DependentUserData<TData, dim> >();
127 0 : }
128 :
129 : template <typename TData, int dim>
130 0 : void DataImport<TData,dim>::cache_data_access()
131 : {
132 : // cache the pointer to the data field.
133 0 : m_vValue = m_spUserData->values(m_seriesID);
134 :
135 : // in addition we cache the number of ips
136 0 : m_numIP = m_spUserData->num_ip(m_seriesID);
137 0 : }
138 :
139 : template <typename TData, int dim>
140 : template <int ldim>
141 0 : 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 0 : if(!data_given()) return;
146 :
147 : // request series if first time requested
148 0 : if(m_seriesID == -1)
149 : {
150 0 : m_seriesID = m_spUserData->template
151 0 : register_local_ip_series<ldim>(vPos,numIP,timePointSpec,bMayChange);
152 :
153 : // register callback, invoked when data field is changed
154 0 : 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 0 : 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 0 : if(!bMayChange)
168 0 : UG_THROW("DataImport: Setting different local ips to non-changable ip series.");
169 :
170 : // set new local ips
171 0 : m_spUserData->template set_local_ips<ldim>(m_seriesID,vPos,numIP);
172 0 : m_spUserData->set_time_point(m_seriesID,timePointSpec);
173 :
174 0 : if(numIP != m_numIP)
175 : {
176 : // cache access to the data
177 : cache_data_access();
178 :
179 : // resize also lin defect array
180 0 : 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>
197 : void DataImport<TData,dim>::set_local_ips(const MathVector<ldim>* vPos, size_t numIP,
198 : bool bMayChange)
199 : {
200 0 : set_local_ips<ldim>(vPos, numIP, -1, bMayChange);
201 : }
202 :
203 : template <typename TData, int dim>
204 : void DataImport<TData,dim>::set_local_ips(const MathVector<dim>* vPos, size_t numIP,
205 : bool bMayChange)
206 : {
207 0 : set_local_ips<dim>(vPos, numIP, -1, bMayChange);
208 : }
209 :
210 : template <typename TData, int dim>
211 : void DataImport<TData,dim>::set_time_point(int timePointSpec)
212 : {
213 : m_spUserData->set_time_point(m_seriesID,timePointSpec);
214 : }
215 :
216 : template <typename TData, int dim>
217 0 : void DataImport<TData,dim>::set_global_ips(const MathVector<dim>* vPos, size_t numIP)
218 : {
219 : // if no data set, skip
220 0 : if(!data_given()) return;
221 :
222 : // set global ips for series ID
223 : UG_ASSERT(m_seriesID >= 0, "Wrong series id.");
224 0 : m_spUserData->set_global_ips(m_seriesID,vPos,numIP);
225 : }
226 :
227 : template <typename TData, int dim>
228 0 : void DataImport<TData,dim>::clear_ips()
229 : {
230 0 : if(data_given()) m_spUserData->unregister_storage_callback(this);
231 0 : m_seriesID = -1;
232 0 : m_vValue = 0;
233 0 : m_numIP = 0;
234 0 : m_vvvLinDefect.resize(num_ip());
235 0 : }
236 :
237 : template <typename TData, int dim>
238 0 : void DataImport<TData,dim>::add_jacobian(LocalMatrix& J, const number scale)
239 : {
240 : UG_ASSERT(m_spDependentUserData.valid(), "No Export set.");
241 :
242 : /// compute the linearization only if the export parameter is 'at current time'
243 0 : if (! m_spUserData->at_current_time (m_seriesID))
244 : return;
245 :
246 : // access jacobian by maps
247 0 : J.access_by_map(this->map(), this->conn_map());
248 :
249 : // loop integration points
250 0 : for(size_t ip = 0; ip < num_ip(); ++ip)
251 : {
252 : // loop all functions
253 0 : for(size_t fct1 = 0; fct1 < this->num_fct(); ++fct1)
254 0 : 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 0 : for(size_t sh1 = 0; sh1 < num_sh(fct1); ++sh1)
262 0 : for(size_t sh2 = 0; sh2 < m_spDependentUserData->num_sh(fct2); ++sh2)
263 : {
264 0 : J(fct1, sh1, fct2, sh2) += scale*(LinDef[sh1]*Deriv[sh2]);
265 : }
266 : }
267 : }
268 : }
269 :
270 : template <typename TData, int dim>
271 0 : void DataImport<TData,dim>::update_dof_sizes(const LocalIndices& ind)
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 0 : m_vvNumDoFPerFct.resize(map.num_fct());
279 0 : for(size_t fct = 0; fct < m_vvNumDoFPerFct.size(); ++fct)
280 0 : m_vvNumDoFPerFct[fct] = ind.num_dof(map[fct]);
281 :
282 : m_vvvLinDefect.clear();
283 0 : resize_defect_array();
284 0 : }
285 :
286 : template <typename TData, int dim>
287 0 : void DataImport<TData,dim>::resize_defect_array()
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 0 : m_vvvLinDefect.resize(num_ip());
295 :
296 : // resize num fct
297 0 : for(size_t ip = oldSize; ip < num_ip(); ++ip)
298 : {
299 : // resize num fct
300 0 : m_vvvLinDefect[ip].resize(m_vvNumDoFPerFct.size());
301 :
302 : // resize dofs
303 0 : for(size_t fct = 0; fct < m_vvNumDoFPerFct.size(); ++fct)
304 0 : m_vvvLinDefect[ip][fct].resize(m_vvNumDoFPerFct[fct]);
305 : }
306 0 : }
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>
330 : inline void DataImport<TData,dim>::check_values() const
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__ */
|