Loading [MathJax]/extensions/tex2jax.js
PolyFEM
Toggle main menu visibility
Main Page
Related Pages
Namespaces
Namespace List
Namespace Members
All
a
b
c
d
e
f
g
h
i
j
l
m
n
o
p
q
r
s
t
u
v
w
x
y
z
Functions
a
b
c
d
e
f
g
h
i
j
l
m
n
o
p
q
r
s
t
u
v
w
y
z
Variables
a
b
c
d
e
f
g
h
i
l
m
n
o
p
r
s
t
u
v
x
y
z
Typedefs
Enumerations
Classes
Class List
Class Index
Class Hierarchy
Class Members
All
_
a
b
c
d
e
f
g
h
i
j
k
l
m
n
o
p
q
r
s
t
u
v
w
x
y
z
~
Functions
_
a
b
c
d
e
f
g
h
i
k
l
m
n
o
p
q
r
s
t
u
v
w
x
z
~
Variables
_
a
b
c
d
e
f
g
h
i
j
k
l
m
n
o
p
q
r
s
t
u
v
w
x
y
Typedefs
b
c
d
e
f
g
h
l
m
n
o
p
q
s
t
v
Enumerations
Related Symbols
a
c
e
g
i
l
o
p
s
Files
File List
File Members
All
_
a
c
d
e
f
g
h
i
j
l
m
o
p
q
r
s
t
u
v
w
x
y
z
Functions
Variables
_
a
c
d
e
f
g
i
l
m
q
s
t
v
w
x
y
z
Macros
•
All
Classes
Namespaces
Files
Functions
Variables
Typedefs
Enumerations
Enumerator
Friends
Macros
Pages
Loading...
Searching...
No Matches
QuadraticBSpline.cpp
Go to the documentation of this file.
1
#include "
QuadraticBSpline.hpp
"
2
3
#include <cmath>
4
5
namespace
polyfem
6
{
7
namespace
basis
8
{
9
void
QuadraticBSpline::init
(
const
std::array<double, 4> &knots)
10
{
11
knots_
= knots;
12
}
9
void
QuadraticBSpline::init
(
const
std::array<double, 4> &knots) {
…
}
13
14
void
QuadraticBSpline::interpolate
(
const
Eigen::MatrixXd &ts, Eigen::MatrixXd &result)
const
15
{
16
result.resize(ts.size(), 1);
17
18
for
(
long
i = 0; i < ts.size(); ++i)
19
result(i) =
interpolate
(ts(i));
20
}
14
void
QuadraticBSpline::interpolate
(
const
Eigen::MatrixXd &ts, Eigen::MatrixXd &result)
const
{
…
}
21
22
double
QuadraticBSpline::interpolate
(
const
double
t)
const
23
{
24
if
(t >=
knots_
[0] && t <
knots_
[1])
25
return
(t -
knots_
[0]) * (t -
knots_
[0]) / (
knots_
[2] -
knots_
[0]) / (
knots_
[1] -
knots_
[0]);
26
27
if
(t >=
knots_
[1] && t <
knots_
[2])
28
return
(t -
knots_
[0]) / (
knots_
[2] -
knots_
[0]) * (
knots_
[2] - t) / (
knots_
[2] -
knots_
[1]) + (
knots_
[3] - t) / (
knots_
[3] -
knots_
[1]) * (t -
knots_
[1]) / (
knots_
[2] -
knots_
[1]);
29
30
if
(t >=
knots_
[2] && t <=
knots_
[3])
31
{
32
if
(fabs(
knots_
[3] -
knots_
[2]) < 1e-12 && fabs(
knots_
[3] -
knots_
[1]) < 1e-12)
33
return
knots_
[3];
34
if
(fabs(
knots_
[3] -
knots_
[2]) < 1e-12)
35
return
0;
36
37
return
(
knots_
[3] - t) * (
knots_
[3] - t) / (
knots_
[3] -
knots_
[1]) / (
knots_
[3] -
knots_
[2]);
38
}
39
40
return
0;
41
}
22
double
QuadraticBSpline::interpolate
(
const
double
t)
const
{
…
}
42
43
void
QuadraticBSpline::derivative
(
const
Eigen::MatrixXd &ts, Eigen::MatrixXd &result)
const
44
{
45
result.resize(ts.size(), 1);
46
47
for
(
long
i = 0; i < ts.size(); ++i)
48
result(i) =
derivative
(ts(i));
49
}
43
void
QuadraticBSpline::derivative
(
const
Eigen::MatrixXd &ts, Eigen::MatrixXd &result)
const
{
…
}
50
51
double
QuadraticBSpline::derivative
(
const
double
t)
const
52
{
53
if
(t >=
knots_
[0] && t <
knots_
[1])
54
return
2 / (
knots_
[2] -
knots_
[0]) * (t -
knots_
[0]) / (
knots_
[1] -
knots_
[0]);
55
56
if
(t >=
knots_
[1] && t <
knots_
[2])
57
return
((-2 * t + 2 *
knots_
[0]) *
knots_
[1] + (2 * t - 2 *
knots_
[3]) *
knots_
[2] - 2 * t * (
knots_
[0] -
knots_
[3])) / (-
knots_
[2] +
knots_
[0]) / (-
knots_
[2] +
knots_
[1]) / (-
knots_
[3] +
knots_
[1]);
58
59
if
(t >=
knots_
[2] && t <=
knots_
[3])
60
{
61
if
(fabs(
knots_
[3] -
knots_
[2]) < 1e-12 && fabs(
knots_
[3] -
knots_
[1]) < 1e-12)
62
return
2 *
knots_
[3];
63
if
(fabs(
knots_
[3] -
knots_
[2]) < 1e-12)
64
return
-2 *
knots_
[3];
65
66
return
(2 * t - 2 *
knots_
[3]) / (-
knots_
[3] +
knots_
[1]) / (-
knots_
[3] +
knots_
[2]);
67
}
68
69
return
0;
70
}
51
double
QuadraticBSpline::derivative
(
const
double
t)
const
{
…
}
71
}
// namespace basis
72
}
// namespace polyfem
QuadraticBSpline.hpp
polyfem::basis::QuadraticBSpline::init
void init(const std::array< double, 4 > &knots)
Definition
QuadraticBSpline.cpp:9
polyfem::basis::QuadraticBSpline::interpolate
void interpolate(const Eigen::MatrixXd &ts, Eigen::MatrixXd &result) const
Definition
QuadraticBSpline.cpp:14
polyfem::basis::QuadraticBSpline::derivative
void derivative(const Eigen::MatrixXd &ts, Eigen::MatrixXd &result) const
Definition
QuadraticBSpline.cpp:43
polyfem::basis::QuadraticBSpline::knots_
std::array< double, 4 > knots_
Definition
QuadraticBSpline.hpp:29
polyfem
Definition
AMIPSEnergy.cpp:6
src
polyfem
basis
function
QuadraticBSpline.cpp
Generated by
1.9.8