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_lattice.h"
22
using namespace rheolef;
23
// -----------------------------------------------------------------------
24
// basis evaluated on nodal lattice of Lagrange basis
25
// -----------------------------------------------------------------------
27
basis_on_nodal_basis::_initialize (reference_element hat_K) const
29
reference_element::variant_type K_type = hat_K.variant();
30
std::vector<point> hat_node;
31
_nb.hat_node (hat_K, hat_node);
32
size_type nq = hat_node.size();
33
_val[K_type].resize (nq);
34
for (size_type q = 0; q < nq; q++) {
35
_b.eval (hat_K, hat_node[q], _val[K_type][q]);
37
_initialized [K_type] = true;
39
basis_on_nodal_basis::const_iterator
40
basis_on_nodal_basis::begin (reference_element hat_K, size_type q) const
42
reference_element::variant_type K_type = hat_K.variant();
43
if (!_initialized [K_type]) _initialize (hat_K);
44
return _val[K_type][q].begin();
46
basis_on_nodal_basis::const_iterator
47
basis_on_nodal_basis::end (reference_element hat_K, size_type q) const
49
reference_element::variant_type K_type = hat_K.variant();
50
if (!_initialized [K_type]) _initialize (hat_K);
51
return _val[K_type][q].end();
53
// -----------------------------------------------------------------------
54
// basis evaluated on lattice of quadrature formulae
55
// -----------------------------------------------------------------------
57
basis_on_quadrature::_initialize (reference_element hat_K) const
59
reference_element::variant_type K_type = hat_K.variant();
60
quadrature::const_iterator first = _p_quad -> begin(hat_K);
61
size_type nq = _p_quad -> size (hat_K);
62
_val[K_type].resize (nq);
63
for (size_type q = 0; q < nq; first++, q++) {
64
const point& hat_node_q = (*first).x;
65
_b.eval (hat_K, hat_node_q, _val[K_type][q]);
67
_initialized [K_type] = true;
69
basis_on_quadrature::const_iterator
70
basis_on_quadrature::begin (reference_element hat_K, size_type q) const
72
reference_element::variant_type K_type = hat_K.variant();
73
if (!_initialized [K_type]) _initialize (hat_K);
74
return _val[K_type][q].begin();
76
basis_on_quadrature::const_iterator
77
basis_on_quadrature::end (reference_element hat_K, size_type q) const
79
reference_element::variant_type K_type = hat_K.variant();
80
if (!_initialized [K_type]) _initialize (hat_K);
81
return _val[K_type][q].end();
85
basis_on_quadrature::_grad_initialize (reference_element hat_K) const
87
reference_element::variant_type K_type = hat_K.variant();
88
quadrature::const_iterator first = _p_quad -> begin(hat_K);
89
size_type nq = _p_quad -> size (hat_K);
90
_grad_val[K_type].resize (nq);
91
for (size_type q = 0; q < nq; first++, q++) {
92
const point& hat_node_q = (*first).x;
93
_b.grad_eval (hat_K, hat_node_q, _grad_val[K_type][q]);
95
_grad_initialized [K_type] = true;
97
basis_on_quadrature::const_iterator_grad
98
basis_on_quadrature::begin_grad (reference_element hat_K, size_type q) const
100
reference_element::variant_type K_type = hat_K.variant();
101
if (!_grad_initialized [K_type]) _grad_initialize (hat_K);
102
return _grad_val[K_type][q].begin();
104
basis_on_quadrature::const_iterator_grad
105
basis_on_quadrature::end_grad (reference_element hat_K, size_type q) const
107
reference_element::variant_type K_type = hat_K.variant();
108
if (!_grad_initialized [K_type]) _grad_initialize (hat_K);
109
return _grad_val[K_type][q].end();
112
basis_on_quadrature::_hessian_initialize (reference_element hat_K) const
114
reference_element::variant_type K_type = hat_K.variant();
115
quadrature::const_iterator first = _p_quad -> begin(hat_K);
116
size_type nq = _p_quad -> size (hat_K);
117
_hessian_val[K_type].resize (nq);
118
for (size_type q = 0; q < nq; first++, q++) {
119
const point& hat_node_q = (*first).x;
120
_b.hessian_eval (hat_K, hat_node_q, _hessian_val[K_type][q]);
122
_hessian_initialized [K_type] = true;
124
basis_on_quadrature::const_iterator_hessian
125
basis_on_quadrature::begin_hessian (reference_element hat_K, size_type q) const
127
reference_element::variant_type K_type = hat_K.variant();
128
if (!_hessian_initialized [K_type]) _hessian_initialize (hat_K);
129
return _hessian_val[K_type][q].begin();
131
basis_on_quadrature::const_iterator_hessian
132
basis_on_quadrature::end_hessian (reference_element hat_K, size_type q) const
134
reference_element::variant_type K_type = hat_K.variant();
135
if (!_hessian_initialized [K_type]) _hessian_initialize (hat_K);
136
return _hessian_val[K_type][q].end();