(git:e7e05ae)
qs_cdft_opt_types.F
Go to the documentation of this file.
1 !--------------------------------------------------------------------------------------------------!
2 ! CP2K: A general program to perform molecular dynamics simulations !
3 ! Copyright 2000-2024 CP2K developers group <https://cp2k.org> !
4 ! !
5 ! SPDX-License-Identifier: GPL-2.0-or-later !
6 !--------------------------------------------------------------------------------------------------!
7 
8 ! **************************************************************************************************
9 !> \brief Control parameters for optimizers that work with CDFT constraints
10 !> \par History
11 !> separated from scf_control_types [03.2018]
12 !> \author Nico Holmberg [03.2018]
13 ! **************************************************************************************************
15 
16  USE input_constants, ONLY: &
21  section_vals_type,&
23  USE kinds, ONLY: dp
24 #include "./base/base_uses.f90"
25 
26  IMPLICIT NONE
27 
28  PRIVATE
29 
30  CHARACTER(len=*), PARAMETER, PRIVATE :: moduleN = 'qs_cdft_opt_types'
31 
32  ! Public data types
33 
34  PUBLIC :: cdft_opt_type
35 
36  ! Public subroutines
37 
38  PUBLIC :: cdft_opt_type_create, &
43 
44 ! **************************************************************************************************
45 !> \brief contains the parameters needed by CDFT specific optimizers
46 !> \param build_jacobian logical which determines if the inverse Jacobian should be computed
47 !> \param jacobian_step the step size for calculating the finite difference Jacobian
48 !> \param newton_step the step size used by the Newton optimizer with values between 0 and 1
49 !> \param newton_step_save permanent copy of the above
50 !> \param jacobian_type the finite difference scheme to compute the Jacobian
51 !> \param broyden_type the variant of Broyden's method to use
52 !> \param jacobian_freq control parameters defining how often the Jacobian is built
53 !> \param ijacobian counter to track how many SCF iterations/energy evaluations have passed since
54 !> the last Jacobian rebuild
55 !> \param broyden_update logical which determines if a Broyden update is needed
56 !> \param max_ls the maximum number of backtracking line search steps to perform
57 !> \param continue_ls continue line search until max steps are reached or until the gradient
58 !> no longer decreases
59 !> \param factor_ls line search parameter used in generating a new step size
60 !> \par History
61 !> created [03.2018]
62 !> \author Nico Holmberg [03.2018]
63 ! **************************************************************************************************
64 
65  TYPE cdft_opt_type
66  LOGICAL :: build_jacobian
67  LOGICAL :: broyden_update
68  LOGICAL :: continue_ls
69  LOGICAL :: jacobian_restart
70  REAL(KIND=dp) :: newton_step
71  REAL(KIND=dp) :: newton_step_save
72  REAL(KIND=dp) :: factor_ls
73  REAL(KIND=dp), DIMENSION(:), &
74  ALLOCATABLE :: jacobian_step
75  REAL(KIND=dp), DIMENSION(:), &
76  POINTER :: jacobian_vector
77  INTEGER :: jacobian_type
78  INTEGER :: broyden_type
79  INTEGER :: jacobian_freq(2)
80  INTEGER :: ijacobian(2)
81  INTEGER :: max_ls
82  END TYPE cdft_opt_type
83 
84 CONTAINS
85 
86 ! **************************************************************************************************
87 !> \brief allocates and initializes the CDFT optimizer control object with default values
88 !> \param cdft_opt_control the object to initialize
89 !> \par History
90 !> 03.2018 created [Nico Holmberg]
91 !> \author Nico Holmberg
92 ! **************************************************************************************************
93  SUBROUTINE cdft_opt_type_create(cdft_opt_control)
94 
95  TYPE(cdft_opt_type), POINTER :: cdft_opt_control
96 
97  CHARACTER(LEN=*), PARAMETER :: routinen = 'cdft_opt_type_create'
98 
99  INTEGER :: handle
100 
101  CALL timeset(routinen, handle)
102 
103  cpassert(.NOT. ASSOCIATED(cdft_opt_control))
104  ALLOCATE (cdft_opt_control)
105 
106  ! Load the default values
107 
108  cdft_opt_control%jacobian_type = -1
109  cdft_opt_control%broyden_type = -1
110  cdft_opt_control%jacobian_freq(:) = 1
111  cdft_opt_control%newton_step = 1.0_dp
112  cdft_opt_control%newton_step_save = 1.0_dp
113  cdft_opt_control%factor_ls = 0.5_dp
114  cdft_opt_control%ijacobian(:) = 0
115  cdft_opt_control%max_ls = 0
116  cdft_opt_control%build_jacobian = .false.
117  cdft_opt_control%broyden_update = .false.
118  cdft_opt_control%continue_ls = .false.
119  cdft_opt_control%jacobian_restart = .false.
120  NULLIFY (cdft_opt_control%jacobian_vector)
121 
122  CALL timestop(handle)
123 
124  END SUBROUTINE cdft_opt_type_create
125 
126 ! **************************************************************************************************
127 !> \brief releases the CDFT optimizer control object
128 !> \param cdft_opt_control the object to release
129 !> \par History
130 !> 03.2018 created [Nico Holmberg]
131 !> \author Nico Holmberg
132 ! **************************************************************************************************
133  SUBROUTINE cdft_opt_type_release(cdft_opt_control)
134 
135  TYPE(cdft_opt_type), POINTER :: cdft_opt_control
136 
137  IF (ASSOCIATED(cdft_opt_control)) THEN
138  IF (ASSOCIATED(cdft_opt_control%jacobian_vector)) &
139  DEALLOCATE (cdft_opt_control%jacobian_vector)
140  IF (ALLOCATED(cdft_opt_control%jacobian_step)) &
141  DEALLOCATE (cdft_opt_control%jacobian_step)
142 
143  DEALLOCATE (cdft_opt_control)
144  END IF
145 
146  NULLIFY (cdft_opt_control)
147 
148  END SUBROUTINE cdft_opt_type_release
149 
150 ! **************************************************************************************************
151 !> \brief reads the parameters of the CDFT optimizer type
152 !> \param cdft_opt_control the object that will contain the values read
153 !> \param inp_section the input section that contains the values that are read
154 !> \par History
155 !> 03.2018 created [Nico Holmberg]
156 !> \author Nico Holmberg
157 ! **************************************************************************************************
158  SUBROUTINE cdft_opt_type_read(cdft_opt_control, inp_section)
159 
160  TYPE(cdft_opt_type), POINTER :: cdft_opt_control
161  TYPE(section_vals_type), POINTER :: inp_section
162 
163  CHARACTER(LEN=*), PARAMETER :: routinen = 'cdft_opt_type_read'
164 
165  INTEGER :: handle
166  INTEGER, DIMENSION(:), POINTER :: tmplist
167  LOGICAL :: exists
168  REAL(kind=dp), DIMENSION(:), POINTER :: rtmplist
169  TYPE(section_vals_type), POINTER :: cdft_opt_section
170 
171  CALL timeset(routinen, handle)
172 
173  cpassert(ASSOCIATED(cdft_opt_control))
174  cdft_opt_section => section_vals_get_subs_vals(inp_section, "CDFT_OPT")
175 
176  CALL section_vals_val_get(cdft_opt_section, "MAX_LS", &
177  i_val=cdft_opt_control%max_ls)
178  CALL section_vals_val_get(cdft_opt_section, "JACOBIAN_TYPE", &
179  i_val=cdft_opt_control%jacobian_type)
180  CALL section_vals_val_get(cdft_opt_section, "JACOBIAN_STEP", &
181  r_vals=rtmplist)
182  ALLOCATE (cdft_opt_control%jacobian_step(SIZE(rtmplist)))
183  cdft_opt_control%jacobian_step(:) = rtmplist
184  CALL section_vals_val_get(cdft_opt_section, "BROYDEN_TYPE", &
185  i_val=cdft_opt_control%broyden_type)
186  CALL section_vals_val_get(cdft_opt_section, "CONTINUE_LS", &
187  l_val=cdft_opt_control%continue_ls)
188  CALL section_vals_val_get(cdft_opt_section, "FACTOR_LS", &
189  r_val=cdft_opt_control%factor_ls)
190  IF (cdft_opt_control%factor_ls .LE. 0.0_dp .OR. &
191  cdft_opt_control%factor_ls .GE. 1.0_dp) &
192  CALL cp_abort(__location__, &
193  "Keyword FACTOR_LS must be between 0.0 and 1.0.")
194  CALL section_vals_val_get(cdft_opt_section, "JACOBIAN_FREQ", explicit=exists)
195  IF (exists) THEN
196  CALL section_vals_val_get(cdft_opt_section, "JACOBIAN_FREQ", &
197  i_vals=tmplist)
198  IF (SIZE(tmplist) /= 2) &
199  CALL cp_abort(__location__, &
200  "Keyword JACOBIAN_FREQ takes exactly two input values.")
201  IF (any(tmplist .LT. 0)) &
202  CALL cp_abort(__location__, &
203  "Keyword JACOBIAN_FREQ takes only positive values.")
204  IF (all(tmplist .EQ. 0)) &
205  CALL cp_abort(__location__, &
206  "Both values to keyword JACOBIAN_FREQ cannot be zero.")
207  cdft_opt_control%jacobian_freq(:) = tmplist(1:2)
208  END IF
209  CALL section_vals_val_get(cdft_opt_section, "JACOBIAN_RESTART", &
210  l_val=cdft_opt_control%jacobian_restart)
211  IF (cdft_opt_control%jacobian_restart) THEN
212  CALL section_vals_val_get(cdft_opt_section, "JACOBIAN_VECTOR", &
213  r_vals=rtmplist)
214  ALLOCATE (cdft_opt_control%jacobian_vector(SIZE(rtmplist)))
215  cdft_opt_control%jacobian_vector = rtmplist
216  END IF
217 
218  CALL timestop(handle)
219 
220  END SUBROUTINE cdft_opt_type_read
221 
222 ! **************************************************************************************************
223 !> \brief writes information about the CDFT optimizer object
224 !> \param cdft_opt_control the CDFT optimizer object
225 !> \param optimizer the type of optimizer to use
226 !> \param output_unit the output unit handle
227 !> \par History
228 !> 03.2018 created [Nico Holmberg]
229 !> \author Nico Holmberg
230 ! **************************************************************************************************
231  SUBROUTINE cdft_opt_type_write(cdft_opt_control, optimizer, output_unit)
232  TYPE(cdft_opt_type), POINTER :: cdft_opt_control
233  INTEGER :: optimizer, output_unit
234 
235  cpassert(ASSOCIATED(cdft_opt_control))
236 
237  SELECT CASE (optimizer)
238  CASE DEFAULT
239  ! Do nothing
241  WRITE (output_unit, '(T3,A)') "Optimization with Broyden's method"
242  SELECT CASE (cdft_opt_control%broyden_type)
243  CASE (broyden_type_1)
244  WRITE (output_unit, '(A)') " variant : 1st method"
246  WRITE (output_unit, '(A)') " variant : 1st method with explicit initial Jacobian"
247  CASE (broyden_type_1_ls)
248  WRITE (output_unit, '(A)') " variant : 1st method with backtracking line search"
250  WRITE (output_unit, '(A)') &
251  " variant : 1st method with explicit initial Jacobian"
252  WRITE (output_unit, '(A)') &
253  " and backtracking line search"
254  CASE (broyden_type_2)
255  WRITE (output_unit, '(A)') " variant : 2nd method"
257  WRITE (output_unit, '(A)') " variant : 2nd method with explicit initial Jacobian"
258  CASE (broyden_type_2_ls)
259  WRITE (output_unit, '(A)') " variant : 2nd method with backtracking line search"
261  WRITE (output_unit, '(A)') &
262  " variant : 2nd method with explicit initial Jacobian"
263  WRITE (output_unit, '(A)') &
264  " and backtracking line search"
265  END SELECT
267  WRITE (output_unit, '(T3,A)') "Optimization with Newton's method"
269  WRITE (output_unit, '(T3,A)') "Optimization with Newton's method using backtracking line search"
270  END SELECT
271  SELECT CASE (optimizer)
272  CASE DEFAULT
273  ! Do nothing
275  IF (cdft_opt_control%jacobian_freq(2) > 0) THEN
276  WRITE (output_unit, '(T6,A,I4,A)') &
277  "The Jacobian is restarted every ", cdft_opt_control%jacobian_freq(2), " energy evaluation"
278  IF (cdft_opt_control%jacobian_freq(1) > 0) &
279  WRITE (output_unit, '(T29,A,I4,A)') &
280  "or every ", cdft_opt_control%jacobian_freq(1), " CDFT SCF iteration"
281  ELSE
282  WRITE (output_unit, '(T6,A,I4,A)') &
283  "The Jacobian is restarted every ", cdft_opt_control%jacobian_freq(1), " CDFT SCF iteration"
284  END IF
285  WRITE (output_unit, '(T3,A,F8.4)') &
286  "Optimizer step size: ", cdft_opt_control%newton_step_save
287  END SELECT
288 
289  END SUBROUTINE cdft_opt_type_write
290 
291 ! **************************************************************************************************
292 !> \brief copies settings between two CDFT optimizer control objects retaining both
293 !> \param new the object where to copy the settings
294 !> \param old the object from where to copy the settings
295 !> \par History
296 !> 03.2018 created [Nico Holmberg]
297 !> \author Nico Holmberg
298 ! **************************************************************************************************
299  SUBROUTINE cdft_opt_type_copy(new, old)
300 
301  TYPE(cdft_opt_type), POINTER :: new, old
302 
303  CHARACTER(LEN=*), PARAMETER :: routinen = 'cdft_opt_type_copy'
304 
305  INTEGER :: handle
306 
307  ! Do nothing if cdft_opt_type is not allocated
308  ! this happens if CDFT is performed with an optimizer other than Broyden/Newton
309  IF (.NOT. ASSOCIATED(old)) RETURN
310 
311  CALL timeset(routinen, handle)
312 
313  IF (.NOT. ASSOCIATED(new)) CALL cdft_opt_type_create(new)
314  new%max_ls = old%max_ls
315  new%continue_ls = old%continue_ls
316  new%factor_ls = old%factor_ls
317  new%jacobian_type = old%jacobian_type
318  new%jacobian_freq(:) = old%jacobian_freq(:)
319  new%newton_step = old%newton_step
320  new%newton_step_save = old%newton_step_save
321  new%ijacobian(:) = old%ijacobian(:)
322  new%build_jacobian = old%build_jacobian
323  new%broyden_type = old%broyden_type
324  new%broyden_update = old%broyden_update
325  IF (ALLOCATED(new%jacobian_step)) DEALLOCATE (new%jacobian_step)
326  ALLOCATE (new%jacobian_step(SIZE(old%jacobian_step)))
327  new%jacobian_step(:) = old%jacobian_step
328  IF (old%jacobian_restart) THEN
329  ! Transfer restart vector for inverse Jacobian matrix
330  ! (qs_calculate_inverse_jacobian handles deallocation of transferred vector)
331  new%jacobian_restart = .true.
332  ALLOCATE (new%jacobian_vector(SIZE(old%jacobian_vector)))
333  new%jacobian_vector = old%jacobian_vector
334  DEALLOCATE (old%jacobian_vector)
335  old%jacobian_restart = .false.
336  END IF
337 
338  CALL timestop(handle)
339 
340  END SUBROUTINE cdft_opt_type_copy
341 
342 END MODULE qs_cdft_opt_types
collects all constants needed in input so that they can be used without circular dependencies
integer, parameter, public broyden_type_2_explicit_ls
integer, parameter, public broyden_type_1_explicit
integer, parameter, public broyden_type_2_ls
integer, parameter, public broyden_type_1
integer, parameter, public outer_scf_optimizer_broyden
integer, parameter, public broyden_type_1_explicit_ls
integer, parameter, public broyden_type_2_explicit
integer, parameter, public broyden_type_2
integer, parameter, public outer_scf_optimizer_newton_ls
integer, parameter, public outer_scf_optimizer_newton
integer, parameter, public broyden_type_1_ls
objects that represent the structure of input sections and the data contained in an input section
recursive type(section_vals_type) function, pointer, public section_vals_get_subs_vals(section_vals, subsection_name, i_rep_section, can_return_null)
returns the values of the requested subsection
subroutine, public section_vals_val_get(section_vals, keyword_name, i_rep_section, i_rep_val, n_rep_val, val, l_val, i_val, r_val, c_val, l_vals, i_vals, r_vals, c_vals, explicit)
returns the requested value
Defines the basic variable types.
Definition: kinds.F:23
integer, parameter, public dp
Definition: kinds.F:34
Control parameters for optimizers that work with CDFT constraints.
subroutine, public cdft_opt_type_create(cdft_opt_control)
allocates and initializes the CDFT optimizer control object with default values
subroutine, public cdft_opt_type_release(cdft_opt_control)
releases the CDFT optimizer control object
subroutine, public cdft_opt_type_copy(new, old)
copies settings between two CDFT optimizer control objects retaining both
subroutine, public cdft_opt_type_read(cdft_opt_control, inp_section)
reads the parameters of the CDFT optimizer type
subroutine, public cdft_opt_type_write(cdft_opt_control, optimizer, output_unit)
writes information about the CDFT optimizer object