![]() |
MFEM v4.8.0
Finite element discretization library
|
MMA (Method of Moving Asymptotes) solves an optimization problem of the form: More...
#include <mma.hpp>
Public Member Functions | |
MMA (int nVar, int nCon, real_t *xval, int iterationNumber=0) | |
Serial MMA. | |
MMA (const int nVar, int nCon, Vector &xval, int iterationNumber=0) | |
MMA (MPI_Comm comm_, int nVar, int nCon, real_t *xval, int iterationNumber=0) | |
MMA (MPI_Comm comm_, const int &nVar, const int &nCon, const Vector &xval, int iterationNumber=0) | |
~MMA () | |
Destructor. | |
void | Update (const Vector &dfdx, const Vector &gx, const Vector &dgdx, const Vector &xmin, const Vector &xmax, Vector &xval) |
void | Update (const Vector &dfdx, const Vector &xmin, const Vector &xmax, Vector &xval) |
Unconstrained. | |
void | SetIteration (int iterationNumber) |
int | GetIteration () |
void | SetPrintLevel (int print_lvl) |
Protected Attributes | |
::std::unique_ptr< real_t[]> | a |
::std::unique_ptr< real_t[]> | b |
::std::unique_ptr< real_t[]> | c |
::std::unique_ptr< real_t[]> | d |
real_t | a0 |
real_t | machineEpsilon |
real_t | epsimin |
real_t | z |
real_t | zet |
int | nCon |
int | nVar |
int | iter = 0 |
int | print_level = 1 |
::std::unique_ptr< real_t[]> | low |
::std::unique_ptr< real_t[]> | upp |
::std::unique_ptr< real_t[]> | x |
::std::unique_ptr< real_t[]> | y |
::std::unique_ptr< real_t[]> | xsi |
::std::unique_ptr< real_t[]> | eta |
::std::unique_ptr< real_t[]> | lam |
::std::unique_ptr< real_t[]> | mu |
::std::unique_ptr< real_t[]> | s |
Friends | |
class | MMASubSvanberg |
MMA (Method of Moving Asymptotes) solves an optimization problem of the form:
Find x that minimizes the objective function F(x), subject to C(x)_i <= 0, for all i = 1, ... m x_lo <= x <= x_hi.
The objective functions are replaced by convex functions chosen based on gradient information, and solved using a dual method. The unique optimal solution of this subproblem is returned as the next iteration point. Optimality is determined by the KKT conditions.
The "Update" function in MMA advances the optimization and must be called in every optimization iteration. Current and previous iteration points construct the "moving asymptotes". The design variables, objective function, constraints are passed to an approximating subproblem. The design variables are updated and returned. Its implementation closely follows the original formulation of 'Svanberg, K. (2007). MMA and GCMMA-two methods for nonlinear optimization. vol, 1, 1-15.'
When used in parallel, all Vectors are assumed to be true dof vectors, and the operators are expected to be defined for tdof vectors.
mfem::MMA::MMA | ( | int | nVar, |
int | nCon, | ||
real_t * | xval, | ||
int | iterationNumber = 0 ) |
mfem::MMA::MMA | ( | const int | nVar, |
int | nCon, | ||
Vector & | xval, | ||
int | iterationNumber = 0 ) |
mfem::MMA::MMA | ( | MPI_Comm | comm_, |
int | nVar, | ||
int | nCon, | ||
real_t * | xval, | ||
int | iterationNumber = 0 ) |
mfem::MMA::MMA | ( | MPI_Comm | comm_, |
const int & | nVar, | ||
const int & | nCon, | ||
const Vector & | xval, | ||
int | iterationNumber = 0 ) |
void mfem::MMA::Update | ( | const Vector & | dfdx, |
const Vector & | gx, | ||
const Vector & | dgdx, | ||
const Vector & | xmin, | ||
const Vector & | xmax, | ||
Vector & | xval ) |
Update the optimization parameters dfdx[nVar] - gradients of the objective gx[nCon] - values of the constraints dgdx[nCon*nVar] - gradients of the constraints ordered constraint by constraint, e.g. {dg0dx0, dg0dx1, ... ,} {dg1dx0, dg1dx1, ... ,} xmin[nVar] - lower bounds xmax[nVar] - upper bounds xval[nVar] - input/output for optimization parameters