DFT-FE 1.1.0-pre
Density Functional Theory With Finite-Elements
Loading...
Searching...
No Matches
geoOptIon.h
Go to the documentation of this file.
1// ---------------------------------------------------------------------
2//
3// Copyright (c) 2017-2025 The Regents of the University of Michigan and DFT-FE
4// authors.
5//
6// This file is part of the DFT-FE code.
7//
8// The DFT-FE code is free software; you can use it, redistribute
9// it, and/or modify it under the terms of the GNU Lesser General
10// Public License as published by the Free Software Foundation; either
11// version 2.1 of the License, or (at your option) any later version.
12// The full text of the license can be found in the file LICENSE at
13// the top level of the DFT-FE distribution.
14//
15// ---------------------------------------------------------------------
16
17#ifndef geoOptIon_H_
18#define geoOptIon_H_
19#include "constants.h"
21#include "nonLinearSolver.h"
22#include "dftBase.h"
23#include "dftfeWrapper.h"
24
25namespace dftfe
26{
27 /**
28 * @brief problem class for atomic force relaxation solver.
29 *
30 * @author Sambit Das
31 */
32
34 {
35 public:
36 /** @brief Constructor.
37 *
38 * @param _dftPtr pointer to dftClass
39 * @param mpi_comm_parent parent mpi_communicator
40 */
42 const MPI_Comm &mpi_comm_parent,
43 const bool restart = false);
44
45 /**
46 * @brief initializes the data member d_relaxationFlags.
47 *
48 */
49 void
50 init(const std::string &restartPath);
51
52 /**
53 * @brief calls the atomic force relaxation solver.
54 *
55 * Currently we have option of one solver: Polak–Ribière nonlinear CG solver
56 * with secant based line search. In future releases, we will have more
57 * options like BFGS solver.
58 *
59 * @return int total geometry update calls
60 */
61 int
62 run();
63
64 /**
65 * @brief Obtain number of unknowns (total number of force components to be relaxed).
66 *
67 * @return int Number of unknowns.
68 */
69 unsigned int
71
72 /**
73 * @brief Compute function gradient (aka forces).
74 *
75 * @param gradient STL vector for gradient values.
76 */
77 void
78 gradient(std::vector<double> &gradient);
79
80 /**
81 * @brief Update atomic positions.
82 *
83 * @param solution displacement of the atoms with respect to their current position.
84 * The size of the solution vector is equal to the number of unknowns.
85 */
86 void
87 update(const std::vector<double> &solution,
88 const bool computeForces = true,
89 const bool useSingleAtomSolutionsInitialGuess = false);
90
91 /**
92 * @brief create checkpoint file for current domain bounding vectors and atomic coordinates.
93 *
94 */
95 void
97
98 /**
99 * @brief check for convergence.
100 *
101 */
102 bool
103 isConverged() const;
104
105 const MPI_Comm &
107
108 /// not implemented
109 void
110 value(std::vector<double> &functionValue);
111
112 /// not implemented
113 void
114 precondition(std::vector<double> &s, const std::vector<double> &gradient);
115
116 /// not implemented
117 void
118 solution(std::vector<double> &solution);
119
120 /// not implemented
121 std::vector<unsigned int>
123
124 private:
125 /// storage for relaxation flags and external force components for each
126 /// global atom. each atom has three flags corresponding to three components
127 /// (0- no relax, 1- relax) and three external force components
128 std::vector<unsigned int> d_relaxationFlags;
129 std::vector<double> d_externalForceOnAtom;
130 std::vector<std::vector<double>> d_atomLocationsInitial;
131 std::string d_restartPath;
137 /// maximum force component to be relaxed
139
140 /// total number of calls to update()
142
143 /// pointer to dft class
145 std::unique_ptr<nonLinearSolver> d_nonLinearSolverPtr;
146
147 /// parallel communication objects
148 const MPI_Comm mpi_communicator;
149 const unsigned int n_mpi_processes;
150 const unsigned int this_mpi_process;
151
152 /// conditional stream object
153 dealii::ConditionalOStream pcout;
154 };
155
156} // namespace dftfe
157#endif
abstract base class for dft
Definition dftBase.h:34
bool d_solverRestart
Definition geoOptIon.h:134
void precondition(std::vector< double > &s, const std::vector< double > &gradient)
not implemented
std::vector< double > d_externalForceOnAtom
Definition geoOptIon.h:129
std::vector< std::vector< double > > d_atomLocationsInitial
Definition geoOptIon.h:130
bool isConverged() const
check for convergence.
bool d_isScfRestart
Definition geoOptIon.h:135
int d_totalUpdateCalls
total number of calls to update()
Definition geoOptIon.h:141
void gradient(std::vector< double > &gradient)
Compute function gradient (aka forces).
void solution(std::vector< double > &solution)
not implemented
const MPI_Comm & getMPICommunicator()
get MPI communicator.
std::unique_ptr< nonLinearSolver > d_nonLinearSolverPtr
Definition geoOptIon.h:145
int d_solver
Definition geoOptIon.h:136
const MPI_Comm mpi_communicator
parallel communication objects
Definition geoOptIon.h:148
const unsigned int this_mpi_process
Definition geoOptIon.h:150
bool d_isRestart
Definition geoOptIon.h:133
int run()
calls the atomic force relaxation solver.
dealii::ConditionalOStream pcout
conditional stream object
Definition geoOptIon.h:153
unsigned int getNumberUnknowns() const
Obtain number of unknowns (total number of force components to be relaxed).
std::string d_solverRestartPath
Definition geoOptIon.h:132
std::vector< unsigned int > getUnknownCountFlag() const
not implemented
dftBase * d_dftPtr
pointer to dft class
Definition geoOptIon.h:144
double d_maximumAtomForceToBeRelaxed
maximum force component to be relaxed
Definition geoOptIon.h:138
const unsigned int n_mpi_processes
Definition geoOptIon.h:149
void init(const std::string &restartPath)
initializes the data member d_relaxationFlags.
void update(const std::vector< double > &solution, const bool computeForces=true, const bool useSingleAtomSolutionsInitialGuess=false)
Update atomic positions.
geoOptIon(dftBase *dftPtr, const MPI_Comm &mpi_comm_parent, const bool restart=false)
Constructor.
std::string d_restartPath
Definition geoOptIon.h:131
std::vector< unsigned int > d_relaxationFlags
Definition geoOptIon.h:128
void save()
create checkpoint file for current domain bounding vectors and atomic coordinates.
void value(std::vector< double > &functionValue)
not implemented
nonlinearSolverProblem()
Constructor.
Definition pseudoPotentialToDftfeConverter.cc:34