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
y
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
AugmentedLagrangianForm.hpp
Go to the documentation of this file.
1
#pragma once
2
3
#include <
polyfem/solver/forms/Form.hpp
>
4
5
namespace
polyfem::solver
6
{
8
class
AugmentedLagrangianForm
:
public
Form
9
{
10
11
public
:
12
AugmentedLagrangianForm
() {}
13
14
virtual
~AugmentedLagrangianForm
() {}
15
16
virtual
void
update_lagrangian
(
const
Eigen::VectorXd &
x
,
const
double
k_al) = 0;
17
18
virtual
double
compute_error
(
const
Eigen::VectorXd &
x
)
const
= 0;
19
20
inline
void
set_initial_weight
(
const
double
k_al) {
k_al_
= k_al; }
21
22
inline
double
lagrangian_weight
()
const
{
return
k_al_
; }
23
24
inline
const
StiffnessMatrix
&
constraint_matrix
()
const
{
return
A_
; }
25
inline
const
Eigen::MatrixXd &
constraint_value
()
const
{
return
b_
; }
26
27
inline
const
StiffnessMatrix
&
constraint_projection_matrix
()
const
{
return
A_proj_
; }
28
inline
const
Eigen::MatrixXd &
constraint_projection_vector
()
const
{
return
b_proj_
; }
29
30
inline
bool
has_projection
()
const
{
return
A_proj_
.rows() > 0; }
31
32
virtual
bool
can_project
()
const
{
return
false
; }
33
virtual
void
project_gradient
(Eigen::VectorXd &grad)
const
{ assert(
false
); }
34
virtual
void
project_hessian
(
StiffnessMatrix
&hessian)
const
{ assert(
false
); }
35
38
void
set_scale
(
const
double
scale)
override
{
k_scale_
= scale; }
39
40
protected
:
41
inline
double
L_weight
()
const
{
return
1 /
k_scale_
; }
42
inline
double
A_weight
()
const
{
return
k_al_
/
k_scale_
; }
43
44
double
k_al_
;
45
46
Eigen::VectorXd
lagr_mults_
;
47
48
StiffnessMatrix
A_
;
49
Eigen::MatrixXd
b_
;
50
51
StiffnessMatrix
A_proj_
;
52
Eigen::MatrixXd
b_proj_
;
53
private
:
54
double
k_scale_
= 1;
55
};
8
class
AugmentedLagrangianForm
:
public
Form
{
…
};
56
}
// namespace polyfem::solver
Form.hpp
x
int x
Definition
SplineBasis3d.cpp:55
polyfem::solver::AugmentedLagrangianForm
Form of the augmented lagrangian.
Definition
AugmentedLagrangianForm.hpp:9
polyfem::solver::AugmentedLagrangianForm::k_al_
double k_al_
penalty parameter
Definition
AugmentedLagrangianForm.hpp:44
polyfem::solver::AugmentedLagrangianForm::project_hessian
virtual void project_hessian(StiffnessMatrix &hessian) const
Definition
AugmentedLagrangianForm.hpp:34
polyfem::solver::AugmentedLagrangianForm::project_gradient
virtual void project_gradient(Eigen::VectorXd &grad) const
Definition
AugmentedLagrangianForm.hpp:33
polyfem::solver::AugmentedLagrangianForm::can_project
virtual bool can_project() const
Definition
AugmentedLagrangianForm.hpp:32
polyfem::solver::AugmentedLagrangianForm::compute_error
virtual double compute_error(const Eigen::VectorXd &x) const =0
polyfem::solver::AugmentedLagrangianForm::k_scale_
double k_scale_
Definition
AugmentedLagrangianForm.hpp:54
polyfem::solver::AugmentedLagrangianForm::update_lagrangian
virtual void update_lagrangian(const Eigen::VectorXd &x, const double k_al)=0
polyfem::solver::AugmentedLagrangianForm::AugmentedLagrangianForm
AugmentedLagrangianForm()
Definition
AugmentedLagrangianForm.hpp:12
polyfem::solver::AugmentedLagrangianForm::constraint_projection_matrix
const StiffnessMatrix & constraint_projection_matrix() const
Definition
AugmentedLagrangianForm.hpp:27
polyfem::solver::AugmentedLagrangianForm::lagr_mults_
Eigen::VectorXd lagr_mults_
vector of lagrange multipliers
Definition
AugmentedLagrangianForm.hpp:46
polyfem::solver::AugmentedLagrangianForm::b_
Eigen::MatrixXd b_
Constraints value.
Definition
AugmentedLagrangianForm.hpp:49
polyfem::solver::AugmentedLagrangianForm::constraint_value
const Eigen::MatrixXd & constraint_value() const
Definition
AugmentedLagrangianForm.hpp:25
polyfem::solver::AugmentedLagrangianForm::A_proj_
StiffnessMatrix A_proj_
Constraints projection matrix.
Definition
AugmentedLagrangianForm.hpp:51
polyfem::solver::AugmentedLagrangianForm::set_initial_weight
void set_initial_weight(const double k_al)
Definition
AugmentedLagrangianForm.hpp:20
polyfem::solver::AugmentedLagrangianForm::set_scale
void set_scale(const double scale) override
sets the scale for the form
Definition
AugmentedLagrangianForm.hpp:38
polyfem::solver::AugmentedLagrangianForm::L_weight
double L_weight() const
Definition
AugmentedLagrangianForm.hpp:41
polyfem::solver::AugmentedLagrangianForm::b_proj_
Eigen::MatrixXd b_proj_
Constraints projection value.
Definition
AugmentedLagrangianForm.hpp:52
polyfem::solver::AugmentedLagrangianForm::has_projection
bool has_projection() const
Definition
AugmentedLagrangianForm.hpp:30
polyfem::solver::AugmentedLagrangianForm::A_weight
double A_weight() const
Definition
AugmentedLagrangianForm.hpp:42
polyfem::solver::AugmentedLagrangianForm::constraint_projection_vector
const Eigen::MatrixXd & constraint_projection_vector() const
Definition
AugmentedLagrangianForm.hpp:28
polyfem::solver::AugmentedLagrangianForm::constraint_matrix
const StiffnessMatrix & constraint_matrix() const
Definition
AugmentedLagrangianForm.hpp:24
polyfem::solver::AugmentedLagrangianForm::~AugmentedLagrangianForm
virtual ~AugmentedLagrangianForm()
Definition
AugmentedLagrangianForm.hpp:14
polyfem::solver::AugmentedLagrangianForm::lagrangian_weight
double lagrangian_weight() const
Definition
AugmentedLagrangianForm.hpp:22
polyfem::solver::AugmentedLagrangianForm::A_
StiffnessMatrix A_
Constraints matrix.
Definition
AugmentedLagrangianForm.hpp:48
polyfem::solver::Form
Definition
Form.hpp:11
polyfem::solver
Definition
OptState.hpp:16
polyfem::StiffnessMatrix
Eigen::SparseMatrix< double, Eigen::ColMajor > StiffnessMatrix
Definition
Types.hpp:22
src
polyfem
solver
forms
lagrangian
AugmentedLagrangianForm.hpp
Generated by
1.9.8