Skip to content

Commit d4c1b0c

Browse files
19helloFisherd99
authored andcommitted
update bfgs_trad method (deepmodeling#5662)
* update bfgs method * update bfgs method and modify the parameter force in relax_step to be passed by reference * update bfgs method and modify the parameter force in relax_step to be passed by reference * Introduce the differences between the two BFGS methods * Add & in relax_step input parameters
1 parent 281959f commit d4c1b0c

File tree

4 files changed

+11
-7
lines changed

4 files changed

+11
-7
lines changed

docs/advanced/input_files/input-main.md

Lines changed: 1 addition & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -1375,6 +1375,7 @@ These variables are used to control the geometry relaxation.
13751375
- **Description**: The methods to do geometry optimization.
13761376
- cg: using the conjugate gradient (CG) algorithm. Note that there are two implementations of the conjugate gradient (CG) method, see [relax_new](#relax_new).
13771377
- bfgs: using the Broyden–Fletcher–Goldfarb–Shanno (BFGS) algorithm.
1378+
- bfgs_trad: using the traditional Broyden–Fletcher–Goldfarb–Shanno (BFGS) algorithm.
13781379
- cg_bfgs: using the CG method for the initial steps, and switching to BFGS method when the force convergence is smaller than [relax_cg_thr](#relax_cg_thr).
13791380
- sd: using the steepest descent (SD) algorithm.
13801381
- fire: the Fast Inertial Relaxation Engine method (FIRE), a kind of molecular-dynamics-based relaxation algorithm, is implemented in the molecular dynamics (MD) module. The algorithm can be used by setting [calculation](#calculation) to `md` and [md_type](#md_type) to `fire`. Also ionic velocities should be set in this case. See [fire](../md.md#fire) for more details.

docs/advanced/opt.md

Lines changed: 3 additions & 1 deletion
Original file line numberDiff line numberDiff line change
@@ -22,7 +22,9 @@ In the nested procedure mentioned above, we used CG method to perform cell relax
2222

2323
The [BFGS method](https://en.wikipedia.org/wiki/Broyden%E2%80%93Fletcher%E2%80%93Goldfarb%E2%80%93Shanno_algorithm) is a quasi-Newton method for solving nonlinear optimization problem. It belongs to the class of quasi-Newton method where the Hessian matrix is approximated during the optimization process. If the initial point is not far from the extrema, BFGS tends to work better than gradient-based methods.
2424

25-
In ABACUS, we implemented the BFGS method for doing fixed-cell structural relaxation.
25+
There is an alternative traditional BFGS method, which can be called by using the keyword 'bfgs_trad'. The bfgs_trad method is a quasi-Newton method that substitute an approximate matrix B for the Hessian matrix. The main difference between 'bfgs' and 'bfgs_trad' is that 'bfgs' updates the inverse of matrix B while 'bfgs_trad' updates matrix B and obtains the inverse of B by solving the matrix eigenvalues and taking the reciprocal of the eigenvalues. Both methods are mathematically equivalent, but in some cases, 'bfgs_trad' performs better.
26+
27+
In ABACUS, we implemented the BFGS method for doing fixed-cell structural relaxation. Users can choose which implementation of BFGS to call by adding the 'bfgs_trad' or 'bfgs' parameter.
2628

2729
### SD method
2830

source/module_relax/relax_old/bfgs.cpp

Lines changed: 5 additions & 4 deletions
Original file line numberDiff line numberDiff line change
@@ -7,7 +7,7 @@
77
//! initialize H0、H、pos0、force0、force
88
void BFGS::allocate(const int _size)
99
{
10-
alpha=70;
10+
alpha=70;//default value in ase is 70
1111
maxstep=PARAM.inp.relax_bfgs_rmax;
1212
size=_size;
1313
sign =true;
@@ -29,8 +29,9 @@ void BFGS::allocate(const int _size)
2929
steplength = std::vector<double>(size, 0.0);
3030
}
3131

32-
void BFGS::relax_step(ModuleBase::matrix _force,
33-
UnitCell& ucell)
32+
33+
void BFGS::relax_step(const ModuleBase::matrix& _force,UnitCell& ucell)
34+
3435
{
3536
GetPos(ucell,pos);
3637
GetPostaud(ucell,pos_taud);
@@ -380,7 +381,7 @@ void BFGS::IsRestrain(std::vector<std::vector<double>>& dpos)
380381
* ModuleBase::Ry_to_eV / 0.529177<PARAM.inp.force_thr_ev;
381382
}
382383

383-
void BFGS::CalculateLargestGrad(ModuleBase::matrix& _force,UnitCell& ucell)
384+
void BFGS::CalculateLargestGrad(const ModuleBase::matrix& _force,UnitCell& ucell)
384385
{
385386
std::vector<double> grad= std::vector<double>(3*size, 0.0);
386387
int iat = 0;

source/module_relax/relax_old/bfgs.h

Lines changed: 2 additions & 2 deletions
Original file line numberDiff line numberDiff line change
@@ -49,14 +49,14 @@ class BFGS
4949
* @param _size
5050
*/
5151
void allocate(const int _size);//initialize parameters
52-
void relax_step(ModuleBase::matrix _force,UnitCell& ucell);//
52+
void relax_step(const ModuleBase::matrix& _force,UnitCell& ucell);//
5353
void PrepareStep(std::vector<std::vector<double>>& force,std::vector<std::vector<double>>& pos,std::vector<std::vector<double>>& H,std::vector<double>& pos0,std::vector<double>& force0,std::vector<double>& steplength,std::vector<std::vector<double>>& dpos,UnitCell& ucell);
5454
void IsRestrain(std::vector<std::vector<double>>& dpos);
5555

5656
private:
5757
bool sign;
5858

59-
void CalculateLargestGrad(ModuleBase::matrix& _force,UnitCell& ucell);
59+
void CalculateLargestGrad(const ModuleBase::matrix& _force,UnitCell& ucell);
6060
void GetPos(UnitCell& ucell,std::vector<std::vector<double>>& pos);
6161
void GetPostaud(UnitCell& ucell,std::vector<std::vector<double>>& pos_taud);
6262
void Update(std::vector<double>& pos, std::vector<double>& force,std::vector<std::vector<double>>& H,UnitCell& ucell);

0 commit comments

Comments
 (0)