Loading [MathJax]/extensions/tex2jax.js
ug4
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Modules Pages
grid_function_user_data_explicit.h
Go to the documentation of this file.
1/*
2 * Copyright (c) 2014-2015: G-CSC, Goethe University Frankfurt
3 * Author: Arne Nägel
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__CONST_GRID_FUNCTION_USER_DATA__
34#define __H__UG__LIB_DISC__FUNCTION_SPACE__CONST_GRID_FUNCTION_USER_DATA__
35
36
37#include <map>
38#include <string>
39
40// ug4 headers
41
42#include "common/common.h"
43
45
52
53
54namespace ug{
55
66template <typename TImpl, typename TData, typename TGridFunction>
68: public StdUserData<StdExplicitGridFunctionData<TImpl,TData, TGridFunction>, TData, TGridFunction::dim>
69{
70public:
72 static const int dim = TGridFunction::dim;
73
74protected:
76 size_t m_fct;
78
79protected:
81 TImpl& getImpl() {return static_cast<TImpl&>(*this);}
82
84 const TImpl& getImpl() const {return static_cast<const TImpl&>(*this);}
85
86public:
89 (
91 )
92 : m_spGridFct(spGridFct)
93 {}
94
96 virtual bool constant() const {return false;}
97
99
100 virtual bool requires_grid_fct() const {return true;}
101
103 virtual bool continuous() const {return getImpl().continuous(); }
104
105 template <int refDim>
107 const MathVector<dim> vCornerCoords[], bool bDeriv = false)
108 {
109 const int si = this->subset();
110
111 for(size_t s = 0; s < this->num_series(); ++s)
112 {
113 getImpl().template evaluate<refDim>(this->values(s), this->ips(s), this->time(s), si,
114 elem, vCornerCoords,
115 this->template local_ips<refDim>(s), this->num_ip(s),
116 u);
117 }
118 }
119
120 template <int refDim>
122 const MathVector<dim> vCornerCoords[], bool bDeriv = false)
123 {
124 const int si = this->subset();
125
126 for(size_t s = 0; s < this->num_series(); ++s)
127 {
128 getImpl().template evaluate<refDim>(this->values(s), this->ips(s), this->time(s), si,
129 elem, vCornerCoords,
130 this->template local_ips<refDim>(s), this->num_ip(s),
131 &(u->solution(this->time_point(s))));
132 }
133 }
134
135 virtual void compute(LocalVector* u, GridObject* elem,
136 const MathVector<dim> vCornerCoords[], bool bDeriv = false)
137 {
138 UG_ASSERT(elem->base_object_id() == this->dim_local_ips(),
139 "local ip dimension and reference element dimension mismatch.");
140
141 switch(this->dim_local_ips())
142 {
143 case 1: eval<1>(u,elem,vCornerCoords,bDeriv); break;
144 case 2: eval<2>(u,elem,vCornerCoords,bDeriv); break;
145 case 3: eval<3>(u,elem,vCornerCoords,bDeriv); break;
146 default: UG_THROW("StdExplicitGridFunctionData: Dimension not supported.");
147 }
148 }
149
151 const MathVector<dim> vCornerCoords[], bool bDeriv = false){
152
153 UG_ASSERT(elem->base_object_id() == this->dim_local_ips(),
154 "local ip dimension and reference element dimension mismatch.");
155
156 switch(this->dim_local_ips())
157 {
158 case 1: eval<1>(u,elem,vCornerCoords,bDeriv); break;
159 case 2: eval<2>(u,elem,vCornerCoords,bDeriv); break;
160 case 3: eval<3>(u,elem,vCornerCoords,bDeriv); break;
161 default: UG_THROW("StdExplicitGridFunctionData: Dimension not supported.");
162 }
163 }
164
165 virtual void operator() (TData& value,
166 const MathVector<dim>& globIP,
167 number time, int si) const
168 {
169 UG_THROW("StdExplicitGridFunctionData: Solution, element and local ips required "
170 "for evaluation, but not passed. Cannot evaluate.");
171 }
172
173 virtual void operator() (TData vValue[],
174 const MathVector<dim> vGlobIP[],
175 number time, int si, const size_t nip) const
176 {
177 UG_THROW("StdExplicitGridFunctionData: Solution, element and local ips required "
178 "for evaluation, but not passed. Cannot evaluate.");
179 }
180
181 template <int refDim>
182 inline void evaluate(TData vValue[],
183 const MathVector<dim> vGlobIP[],
184 number time, int si,
185 GridObject* elem,
186 const MathVector<dim> vCornerCoords[],
187 const MathVector<refDim> vLocIP[],
188 const size_t nip,
189 LocalVector* u,
190 const MathMatrix<refDim, dim>* vJT = NULL) const
191 {
192 // forward
193 getImpl().template evaluate<refDim> (vValue, vGlobIP, time, si,
194 elem, vCornerCoords, vLocIP, nip, u, vJT);
195 }
196};
197
209template<typename TGridFunction>
211: public StdExplicitGridFunctionData<ExplicitGridFunctionValue<TGridFunction>, number, TGridFunction>
212{
214
216
217 size_t m_fct;
219
220public:
221 // world dimension of grid function
222 static const int dim = TGridFunction::dim;
223
224 // constructor
226 (
228 const char* cmp
229 )
230 : base_type(gf)
231 {
232 // get function id of name
233 m_fct = m_spGridFct->fct_id_by_name(cmp);
234
235 // check that function exists
236 if(m_fct >= m_spGridFct->num_fct())
237 UG_THROW("ExplicitGridFunctionValue: Function space does not contain"
238 " a function with name " << cmp << ".");
239
240 // local finite element id
241 m_lfeID = m_spGridFct->local_finite_element_id(m_fct);
242 }
243
244 template <int refDim>
245 inline void evaluate(number vValue[],
246 const MathVector<dim> vGlobIP[],
247 number time, int si,
248 GridObject* elem,
249 const MathVector<dim> vCornerCoords[],
250 const MathVector<refDim> vLocIP[],
251 const size_t nip,
252 LocalVector* u,
253 const MathMatrix<refDim, dim>* vJT = NULL) const
254 {
255 // reference object id
256 const ReferenceObjectID roid = elem->reference_object_id();
257
258 // get trial space
259 try
260 {
261 const LocalShapeFunctionSet<refDim>& rTrialSpace =
262 LocalFiniteElementProvider::get<refDim>(roid, m_lfeID);
263
264 // memory for shapes
265 std::vector<number> vShape;
266
267 // get multiindices of element
268 std::vector<DoFIndex> ind;
269 m_spGridFct->dof_indices(elem, m_fct, ind);
270
271 // loop ips
272 for(size_t ip = 0; ip < nip; ++ip)
273 {
274 // evaluate at shapes at ip
275 rTrialSpace.shapes(vShape, vLocIP[ip]);
276
277 // compute solution at integration point
278 vValue[ip] = 0.0;
279 for(size_t sh = 0; sh < vShape.size(); ++sh)
280 {
281 const number valSH = DoFRef(*m_spGridFct, ind[sh]);
282 vValue[ip] += valSH * vShape[sh];
283 }
284 }
285
286 }
287 UG_CATCH_THROW("ExplicitGridFunctionValue: Shape Function Set missing for"
288 " Reference Object: "<<roid<<", Trial Space: "
289 <<m_lfeID<<", refDim="<<refDim);
290
291 }
292
294
295};
296
308template<typename TGridFunction>
310: public StdExplicitGridFunctionData<ExplicitGridFunctionVector<TGridFunction>, MathVector<TGridFunction::dim>, TGridFunction>
311{
312public:
313 // world dimension of grid function
314 static const int dim = TGridFunction::dim;
315
316private:
318
320
321 size_t m_vfct[dim];
323
324public:
325
326 // constructor
328 (
330 const char* cmps
331 )
332 : base_type(gf)
333 {
334 try
335 {
336 // get strings
337 std::string fctString (cmps);
338
339 // tokenize strings and select functions
340 std::vector<std::string> tokens;
341 TokenizeString(fctString, tokens, ',');
342
343 for(size_t i = 0; i < tokens.size(); ++i)
345
346 if((int)tokens.size() != dim)
347 UG_THROW("ExplicitGridFunctionVector: Needed " << dim << " components "
348 "in symbolic function names, but given: " << cmps << '.');
349
350 // get function id of name
351 for(int i = 0; i < dim; ++i)
352 {
353 m_vfct[i] = m_spGridFct->fct_id_by_name(tokens[i].c_str());
354 m_vlfeID[i] = m_spGridFct->local_finite_element_id(m_vfct[i]);
355 }
356
357 }
358 UG_CATCH_THROW("ExplicitGridFunctionVector: Cannot find some symbolic function name in '" << cmps << "'.");
359 }
360
361 template <int refDim>
362 inline void evaluate(MathVector<dim> vValue[],
363 const MathVector<dim> vGlobIP[],
364 number time, int si,
365 GridObject* elem,
366 const MathVector<dim> vCornerCoords[],
367 const MathVector<refDim> vLocIP[],
368 const size_t nip,
369 LocalVector* u,
370 const MathMatrix<refDim, dim>* vJT = NULL) const
371 {
372 // reference object id
373 const ReferenceObjectID roid = elem->reference_object_id();
374
375 // memory for shapes
376 std::vector<number> vShape;
377
378 // multiindices of element
379 std::vector<DoFIndex> ind;
380
381 for (int d = 0; d < dim; ++d)
382 {
383 try
384 {
385 // get trial space
386 const LocalShapeFunctionSet<refDim>& rTrialSpace =
387 LocalFiniteElementProvider::get<refDim>(roid, m_vlfeID[d]);
388
389 m_spGridFct->dof_indices(elem, m_vfct[d], ind);
390
391 // loop ips
392 for(size_t ip = 0; ip < nip; ++ip)
393 {
394 // evaluate at shapes at ip
395 rTrialSpace.shapes(vShape, vLocIP[ip]);
396
397 // compute solution at integration point
398 vValue[ip][d] = 0.0;
399 for(size_t sh = 0; sh < vShape.size(); ++sh)
400 {
401 const number valSH = DoFRef(*m_spGridFct, ind[sh]);
402 vValue[ip][d] += valSH * vShape[sh];
403 }
404 }
405 }
406 UG_CATCH_THROW("ExplicitGridFunctionVector: Shape Function Set missing for"
407 " Reference Object: " << roid << ", Trial Space: "
408 << m_vlfeID << ", refDim=" <<refDim);
409 }
410 }
411
412 bool continuous() const
413 {
414 for(int i = 0; i < dim; ++i)
416 return false;
417 return true;
418 }
419
420};
421
433template <typename TGridFunction>
435: public StdExplicitGridFunctionData<ExplicitGridFunctionGradient<TGridFunction>, MathVector<TGridFunction::dim>, TGridFunction >
436{
438
440
441 size_t m_fct;
443 std::map<std::string, double> m_diffCoeffMap;
444
445public:
446 // world dimension of grid function
447 static const int dim = TGridFunction::dim;
448
451 (
453 const char* cmp
454 )
456 {
457 // get function id of name
458 m_fct = m_spGridFct->fct_id_by_name(cmp);
459
460 // check that function exists
461 if(m_fct >= m_spGridFct->num_fct())
462 UG_THROW("ExplicitGridFunctionGradient: Function space does not contain"
463 " a function with name " << cmp << ".");
464
465 // local finite element id
466 m_lfeID = m_spGridFct->local_finite_element_id(m_fct);
467 }
468
469 // evaluate gradient
470 template <int refDim>
472 const MathVector<dim> vGlobIP[],
473 number time, int si,
474 GridObject* elem,
475 const MathVector<dim> vCornerCoords[],
476 const MathVector<refDim> vLocIP[],
477 const size_t nip,
478 LocalVector* u,
479 const MathMatrix<refDim, dim>* vJT = NULL) const
480 {
481 // reference object id
482 const ReferenceObjectID roid = elem->reference_object_id();
483
484 // get reference element mapping by reference object id
485 std::vector<MathMatrix<refDim, dim> > vJTTmp(nip);
486 if(vJT == NULL){
487 try{
489 = ReferenceMappingProvider::get<refDim, dim>(roid, vCornerCoords);
490
491 // compute transformation matrices
492 mapping.jacobian_transposed(&(vJTTmp[0]), vLocIP, nip);
493
494 // store tmp Gradient
495 vJT = &(vJTTmp[0]);
496 }UG_CATCH_THROW("ExplicitGridFunctionGradient: failed.");
497 }
498
499 // scale with coeficient
500 const int subsetInd = m_spGridFct->domain()->subset_handler()->get_subset_index(elem);
501 const char* subsetName = m_spGridFct->domain()->subset_handler()->get_subset_name(subsetInd);
502
503 double diffCoeff = get_subset_coeff(std::string(subsetName));
504
505 // get trial space
506 try{
507 const LocalShapeFunctionSet<refDim>& rTrialSpace =
508 LocalFiniteElementProvider::get<refDim>(roid, m_lfeID);
509
510 // storage for shape function at ip
511 std::vector<MathVector<refDim> > vLocGrad;
512 MathVector<refDim> locGrad;
513
514 // Reference Mapping
516
517 // loop ips
518 for(size_t ip = 0; ip < nip; ++ip)
519 {
520 // evaluate at shapes at ip
521 rTrialSpace.grads(vLocGrad, vLocIP[ip]);
522
523 // get multiindices of element
524 std::vector<DoFIndex > ind;
525 m_spGridFct->dof_indices(elem, m_fct, ind);
526
527 // compute grad at ip
528 VecSet(locGrad, 0.0);
529 for(size_t sh = 0; sh < vLocGrad.size(); ++sh)
530 {
531 const number valSH = DoFRef( *m_spGridFct, ind[sh]);
532 VecScaleAppend(locGrad, valSH, vLocGrad[sh]);
533 }
534
535 Inverse(JTInv, vJT[ip]);
536 MatVecMult(vValue[ip], JTInv, locGrad);
537
538 // scale with diff coeff (TODO: matrix)
539 VecScale(vValue[ip], vValue[ip], diffCoeff);
540 }
541 }
542 UG_CATCH_THROW("ExplicitGridFunctionGradient: Shape Function Set missing for"
543 " Reference Object: "<<roid<<", Trial Space: "
544 <<m_lfeID<<", refDim="<<refDim);
545 }
546
547 bool continuous() const { return false; }
548
549 void add_subset_coeff(const std::string &key, double val)
550 {
551 m_diffCoeffMap.insert(std::pair<std::string, double>(key, val));
552 }
553
554 double get_subset_coeff(const std::string &key) const
555 {
556 std::map<std::string, double>::const_iterator it = m_diffCoeffMap.find(key);
557 double val = (it != m_diffCoeffMap.end()) ? it->second : 1.0;
558 return val;
559 }
560};
561
562
563} // end namespace ug
564
565#endif /* __H__UG__LIB_DISC__FUNCTION_SPACE__GRID_FUNCTION_USER_DATA__ */
parameterString s
Definition smart_pointer.h:108
const TData * values(size_t s) const
returns all values for a series
Definition user_data.h:516
const TData & value(size_t s, size_t ip) const
returns the value at ip
Definition user_data.h:512
size_t num_ip(size_t s) const
returns the number of integration points
Definition user_data.h:327
size_t num_series() const
explicitly forward some functions
Definition user_data.h:324
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_user_data_explicit.h:436
void evaluate(MathVector< dim > vValue[], const MathVector< dim > vGlobIP[], number time, int si, GridObject *elem, const MathVector< dim > vCornerCoords[], const MathVector< refDim > vLocIP[], const size_t nip, LocalVector *u, const MathMatrix< refDim, dim > *vJT=NULL) const
Definition grid_function_user_data_explicit.h:471
bool continuous() const
returns if provided data is continuous over geometric object boundaries
Definition grid_function_user_data_explicit.h:547
void add_subset_coeff(const std::string &key, double val)
Definition grid_function_user_data_explicit.h:549
LFEID m_lfeID
type of the shape functions
Definition grid_function_user_data_explicit.h:442
double get_subset_coeff(const std::string &key) const
Definition grid_function_user_data_explicit.h:554
static const int dim
Definition grid_function_user_data_explicit.h:447
ExplicitGridFunctionGradient(SmartPtr< TGridFunction > gf, const char *cmp)
constructor
Definition grid_function_user_data_explicit.h:451
std::map< std::string, double > m_diffCoeffMap
Definition grid_function_user_data_explicit.h:443
size_t m_fct
component of the function
Definition grid_function_user_data_explicit.h:441
StdExplicitGridFunctionData< ExplicitGridFunctionGradient< TGridFunction >, MathVector< TGridFunction::dim >, TGridFunction > base_type
Definition grid_function_user_data_explicit.h:437
Definition grid_function_user_data_explicit.h:212
bool continuous() const
returns if provided data is continuous over geometric object boundaries
Definition grid_function_user_data_explicit.h:293
StdExplicitGridFunctionData< ExplicitGridFunctionValue< TGridFunction >, number, TGridFunction > base_type
Definition grid_function_user_data_explicit.h:213
void evaluate(number vValue[], const MathVector< dim > vGlobIP[], number time, int si, GridObject *elem, const MathVector< dim > vCornerCoords[], const MathVector< refDim > vLocIP[], const size_t nip, LocalVector *u, const MathMatrix< refDim, dim > *vJT=NULL) const
Definition grid_function_user_data_explicit.h:245
ExplicitGridFunctionValue(SmartPtr< TGridFunction > gf, const char *cmp)
Definition grid_function_user_data_explicit.h:226
LFEID m_lfeID
type of the shape functions
Definition grid_function_user_data_explicit.h:218
static const int dim
Definition grid_function_user_data_explicit.h:222
size_t m_fct
component of the function
Definition grid_function_user_data_explicit.h:217
Definition grid_function_user_data_explicit.h:311
bool continuous() const
returns if provided data is continuous over geometric object boundaries
Definition grid_function_user_data_explicit.h:412
ExplicitGridFunctionVector(SmartPtr< TGridFunction > gf, const char *cmps)
Definition grid_function_user_data_explicit.h:328
LFEID m_vlfeID[dim]
types of the shape functions
Definition grid_function_user_data_explicit.h:322
static const int dim
Definition grid_function_user_data_explicit.h:314
StdExplicitGridFunctionData< ExplicitGridFunctionVector< TGridFunction >, MathVector< TGridFunction::dim >, TGridFunction > base_type
Definition grid_function_user_data_explicit.h:317
void evaluate(MathVector< dim > vValue[], const MathVector< dim > vGlobIP[], number time, int si, GridObject *elem, const MathVector< dim > vCornerCoords[], const MathVector< refDim > vLocIP[], const size_t nip, LocalVector *u, const MathMatrix< refDim, dim > *vJT=NULL) const
Definition grid_function_user_data_explicit.h:362
size_t m_vfct[dim]
components of the function
Definition grid_function_user_data_explicit.h:321
The base class for all geometric objects, such as vertices, edges, faces, volumes,...
Definition grid_base_objects.h:157
virtual int base_object_id() const =0
virtual ReferenceObjectID reference_object_id() const =0
int dim_local_ips() const
returns current local ip dimension
Definition user_data.h:372
int subset() const
returns the subset of evaluation
Definition user_data.h:273
number time() const
get the current evaluation time
Definition user_data.h:285
const MathVector< dim > * ips(size_t s) const
returns global ips
Definition user_data.h:398
const MathVector< dim > & ip(size_t s, size_t ip) const
returns global ip
Definition user_data.h:401
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
Definition local_algebra.h:198
time series of local vectors
Definition solution_time_series.h:167
const LocalVector & solution(size_t i) const
returns the local vector for the i'th time point
Definition solution_time_series.h:182
A class for fixed size, dense matrices.
Definition math_matrix.h:63
a mathematical Vector with N entries.
Definition math_vector.h:97
Definition grid_function_user_data_explicit.h:69
size_t m_fct
component of function
Definition grid_function_user_data_explicit.h:76
void evaluate(TData vValue[], const MathVector< dim > vGlobIP[], number time, int si, GridObject *elem, const MathVector< dim > vCornerCoords[], const MathVector< refDim > vLocIP[], const size_t nip, LocalVector *u, const MathMatrix< refDim, dim > *vJT=NULL) const
Definition grid_function_user_data_explicit.h:182
virtual void compute(LocalVectorTimeSeries *u, GridObject *elem, const MathVector< dim > vCornerCoords[], bool bDeriv=false)
compute values (and derivatives iff compDeriv == true, but only for the 'current' time point)
Definition grid_function_user_data_explicit.h:150
StdExplicitGridFunctionData(SmartPtr< TGridFunction > spGridFct)
common constructor
Definition grid_function_user_data_explicit.h:89
LFEID m_lfeID
Definition grid_function_user_data_explicit.h:77
virtual void operator()(TData &value, const MathVector< dim > &globIP, number time, int si) const
returns value for a global position
Definition grid_function_user_data_explicit.h:165
virtual bool requires_grid_fct() const
returns if grid function is needed for evaluation
Definition grid_function_user_data_explicit.h:100
void eval(LocalVector *u, GridObject *elem, const MathVector< dim > vCornerCoords[], bool bDeriv=false)
Definition grid_function_user_data_explicit.h:106
virtual void compute(LocalVector *u, GridObject *elem, const MathVector< dim > vCornerCoords[], bool bDeriv=false)
compute values (and derivatives iff compDeriv == true)
Definition grid_function_user_data_explicit.h:135
virtual bool continuous() const
returns if provided data is continuous over geometric object boundaries
Definition grid_function_user_data_explicit.h:103
TImpl & getImpl()
access to implementation
Definition grid_function_user_data_explicit.h:81
const TImpl & getImpl() const
const access to implementation
Definition grid_function_user_data_explicit.h:84
static const int dim
world dimension of grid function
Definition grid_function_user_data_explicit.h:72
SmartPtr< TGridFunction > m_spGridFct
grid function
Definition grid_function_user_data_explicit.h:75
void eval(LocalVectorTimeSeries *u, GridObject *elem, const MathVector< dim > vCornerCoords[], bool bDeriv=false)
Definition grid_function_user_data_explicit.h:121
virtual bool constant() const
returns if data is constant
Definition grid_function_user_data_explicit.h:96
Definition std_user_data.h:63
MathMatrix< N, M, T >::value_type Inverse(MathMatrix< N, M, T > &mOut, const MathMatrix< M, N, T > &m)
Inverse of a matrix.
Definition math_matrix_functions_common_impl.hpp:560
void RemoveWhitespaceFromString(std::string &str)
removes all white space from a string, also within the string
Definition string_util.cpp:50
#define UG_ASSERT(expr, msg)
Definition assert.h:70
#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
void VecScale(vector_t &vOut, const vector_t &v, typename vector_t::value_type s)
scales a MathVector<N>
Definition math_vector_functions_common_impl.hpp:252
the ug namespace
ReferenceObjectID
these ids are used to identify the shape of a geometric object.
Definition grid_base_objects.h:74
void TokenizeString(const string &str, vector< string > &vToken, const char delimiter)
Definition string_util.cpp:56
number & DoFRef(TMatrix &mat, const DoFIndex &iInd, const DoFIndex &jInd)
Definition multi_index.h:276