2
/// This file is part of Rheolef.
4
/// Copyright (C) 2000-2009 Pierre Saramito <Pierre.Saramito@imag.fr>
6
/// Rheolef is free software; you can redistribute it and/or modify
7
/// it under the terms of the GNU General Public License as published by
8
/// the Free Software Foundation; either version 2 of the License, or
9
/// (at your option) any later version.
11
/// Rheolef is distributed in the hope that it will be useful,
12
/// but WITHOUT ANY WARRANTY; without even the implied warranty of
13
/// MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
14
/// GNU General Public License for more details.
16
/// You should have received a copy of the GNU General Public License
17
/// along with Rheolef; if not, write to the Free Software
18
/// Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
20
/// =========================================================================
21
#include "rheolef/basis_on_pointset.h"
24
// -----------------------------------------------------------------------
25
// basis evaluated on lattice of quadrature formulae
26
// -----------------------------------------------------------------------
29
basis_on_pointset<T>::_initialize (reference_element hat_K) const
31
reference_element::variant_type K_variant = hat_K.variant();
32
if (_mode == quad_mode) {
33
typename quadrature<T>::const_iterator first = _p_quad -> begin(hat_K);
34
size_type nq = _p_quad -> size (hat_K);
35
_val[K_variant].resize (nq);
36
for (size_type q = 0; q < nq; first++, q++) {
37
const point_basic<T>& hat_node_q = (*first).x;
38
_b.eval (hat_K, hat_node_q, _val[K_variant][q]);
40
} else if (_mode == nodal_mode) {
41
std::vector<point_basic<T> > hat_node;
42
_nb.hat_node (hat_K, hat_node);
43
size_type nq = hat_node.size();
44
_val[K_variant].resize (nq);
45
for (size_type q = 0; q < nq; q++) {
46
_b.eval (hat_K, hat_node[q], _val[K_variant][q]);
49
fatal_macro ("unsupported eval in specific mode");
51
_initialized [K_variant] = true;
55
basis_on_pointset<T>::_grad_initialize (reference_element hat_K) const
57
reference_element::variant_type K_variant = hat_K.variant();
58
if (_mode == quad_mode) {
59
typename quadrature<T>::const_iterator first = _p_quad -> begin(hat_K);
60
size_type nq = _p_quad -> size (hat_K);
61
_grad_val[K_variant].resize (nq);
62
for (size_type q = 0; q < nq; first++, q++) {
63
const point_basic<T>& hat_node_q = (*first).x;
64
_b.grad_eval (hat_K, hat_node_q, _grad_val[K_variant][q]);
67
fatal_macro ("unsupported grad eval in non-quadrature mode");
69
_grad_initialized [K_variant] = true;
72
typename basis_on_pointset<T>::size_type
73
basis_on_pointset<T>::size (reference_element hat_K) const
75
if (_mode == quad_mode) {
76
return _p_quad -> size (hat_K);
77
} else if (_mode == nodal_mode) {
78
return _nb.size (hat_K);
80
fatal_macro ("unsupported size(hat_K) in specific mode");
83
// ----------------------------------------------------------------------------
84
// instanciation in library
85
// ----------------------------------------------------------------------------
86
#define _RHEOLEF_instanciation(T) \
87
template class basis_on_pointset<T>;
89
_RHEOLEF_instanciation(Float)
91
} // namespace rheolef