Octopus
run.F90
Go to the documentation of this file.
1!! Copyright (C) 2002-2006 M. Marques, A. Castro, A. Rubio, G. Bertsch
2!! Copyright (C) 2020 M. Oliveira
3!!
4!! This program is free software; you can redistribute it and/or modify
5!! it under the terms of the GNU General Public License as published by
6!! the Free Software Foundation; either version 2, or (at your option)
7!! any later version.
8!!
9!! This program is distributed in the hope that it will be useful,
10!! but WITHOUT ANY WARRANTY; without even the implied warranty of
11!! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
12!! GNU General Public License for more details.
13!!
14!! You should have received a copy of the GNU General Public License
15!! along with this program; if not, write to the Free Software
16!! Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA
17!! 02110-1301, USA.
18!!
19
20#include "global.h"
21
24module run_oct_m
25 use accel_oct_m
26 use casida_oct_m
30 use fft_oct_m
32 use global_oct_m
38 use lasers_oct_m
41 use mpi_oct_m
48 use parser_oct_m
54 use kdotp_oct_m
56 use pulpo_oct_m
60 use system_oct_m
61 use td_oct_m
62 use test_oct_m
65 use unocc_oct_m
67 use vdw_oct_m
68
69 implicit none
70
71 private
72 public :: &
73 run
74
75 integer, parameter :: LR = 1, fd = 2
76
77contains
78
79 ! ---------------------------------------------------------
81 integer function get_resp_method(namespace)
82 type(namespace_t), intent(in) :: namespace
83
84 push_sub(get_resp_method)
85
86 !%Variable ResponseMethod
87 !%Type integer
88 !%Default sternheimer
89 !%Section Linear Response
90 !%Description
91 !% Some response properties can be calculated either via
92 !% Sternheimer linear response or by using finite
93 !% differences. You can use this variable to select how you want
94 !% them to be calculated, it applies to <tt>em_resp</tt> and <tt>vib_modes</tt>
95 !% calculation modes. By default, the Sternheimer linear-response
96 !% technique is used.
97 !%Option sternheimer 1
98 !% The linear response is obtained by solving a self-consistent
99 !% Sternheimer equation for the variation of the orbitals. This
100 !% is the recommended method.
101 !%Option finite_differences 2
102 !% Properties are calculated as a finite-differences derivative of
103 !% the energy obtained by several ground-state calculations. This
104 !% method, slow and limited only to static response, is kept
105 !% mainly because it is simple and useful for testing purposes.
106 !%End
107
108 call parse_variable(namespace, 'ResponseMethod', lr, get_resp_method)
109
110 if (.not. varinfo_valid_option('ResponseMethod', get_resp_method)) then
111 call messages_input_error(namespace, 'ResponseMethod')
112 end if
113
114 pop_sub(get_resp_method)
115 end function get_resp_method
116
117 ! ---------------------------------------------------------
121 subroutine run(namespace, calc_mode_id)
122 type(namespace_t), intent(in) :: namespace
123 integer, intent(in) :: calc_mode_id
124
125 type(partner_list_t) :: partners
126 class(system_t), pointer :: systems
127 type(system_factory_t) :: system_factory
128 type(interactions_factory_t) :: interactions_factory
129 logical :: from_scratch
130 integer :: iunit_out
131 type(partner_iterator_t) :: iter
132 type(system_list_t) :: flat_list
133 class(interaction_partner_t), pointer :: partner
134 real(real64) :: largest_dt, largest_allowed_time
135
136 push_sub(run)
137
138 call messages_print_with_emphasis(msg="Calculation Mode", namespace=namespace)
139 call messages_print_var_option("CalculationMode", calc_mode_id, namespace=namespace)
140 call messages_print_with_emphasis(namespace=namespace)
141
142 call calc_mode_init()
143
144 if (calc_mode_id == option__calculationmode__recipe) then
145 call pulpo_print()
146 pop_sub(run)
147 return
148 end if
149
150 call restart_module_init(namespace)
151
152 call unit_system_init(namespace)
153
154 call accel_init(mpi_world, namespace)
155
156 ! initialize FFTs
157 call fft_all_init(namespace)
158
159 if (calc_mode_id == option__calculationmode__test) then
160 call test_run(namespace)
161 call fft_all_end()
163 pop_sub(run)
164 return
165 end if
166
167 ! Create systems
168 if (parse_is_defined(namespace, "Systems")) then
169 ! We are running in multi-system mode
170 systems => system_factory%create(namespace, system_multisystem)
171 else
172 ! Fall back to old behaviour
173 systems => electrons_t(namespace, generate_epot = calc_mode_id /= option__calculationmode__dummy)
174 end if
175
176 ! initialize everything that needs parallelization
177 call systems%init_parallelization(mpi_world)
178
179 ! Create list of partners
180 select type (systems)
181 class is (multisystem_basic_t)
182 ! Systems are also partners
183 partners = systems%list
184 ! Add external potentials to partners list
185 call load_external_potentials(partners, namespace)
186
187 call load_external_waves(partners, namespace)
188
189 ! Add lasers to the partner list
190 call load_lasers(partners, namespace)
191
192 type is (electrons_t)
193 call partners%add(systems)
194 end select
195
196 ! Initialize algorithms (currently only propagators are supported)
197 select type (systems)
198 class is (multisystem_basic_t)
199 select case (calc_mode_id)
200 case (option__calculationmode__gs)
201 call systems%new_algorithm(minimizer_factory_t(systems%namespace))
202 case (option__calculationmode__td)
203 call systems%new_algorithm(propagator_factory_t(systems%namespace))
204 end select
205 end select
206
207
208 ! Check whether the final propagation time would lead to an infinite loop, and stop the code if necessary.
209 ! See issue: https://gitlab.com/octopus-code/octopus/-/issues/1226
210 !
211 ! The deadlock currently appears when systems have different timesteps, and one system could do a
212 !
213 ! 2 4 6 8 10 propagation time=11
214 ! system1 (dt=2): * * * * * |
215 ! system2 (dt=4): * * |
216 !
217 ! As system2 is not allowed to do to timestep 12, it is stuck at 8. System1 will try to get to 10, but is stuck at a barrier, as
218 ! system2 is behind, and keeps waiting for system2.
219 !
220 ! We can only allow propagations up to the largest time, the fastest system (largest dt)
221 ! can reach within the given propagation time. If one of the slower systems could perform
222 ! an extra step within the propagation time, we need to stop the calculation, before the
223 ! propagation, to prevent wasting CPU time in the deadlock.
224 !
225 ! TODO: Fix the framework, so that this infinite loop can be avoided in the first place.
226
227 if(calc_mode_id == option__calculationmode__td) then
228
229 select type (sys => systems)
230 type is (multisystem_basic_t)
231
232 largest_dt = m_zero
233 largest_allowed_time = m_zero
234
235 call sys%get_flat_list(flat_list)
236
237 call iter%start(flat_list)
238 do while (iter%has_next())
239 select type (subsystem => iter%get_next())
240 class is (system_t)
241 select type (subalgorithm => subsystem%algo)
242 class is (propagator_t)
243 largest_dt = max(largest_dt, subalgorithm%dt)
244 end select
245 end select
246 end do
247
248 select type(prop => systems%algo)
249 class is (propagator_t)
250 largest_allowed_time = floor(prop%final_time / largest_dt)* largest_dt
251 class default
252 assert(.false.)
253 end select
254
255 call iter%start(flat_list)
256 do while (iter%has_next())
257 select type(subsystem => iter%get_next())
258 class is (system_t)
259 select type (subalgorithm => subsystem%algo)
260 class is (propagator_t)
261 if( floor(subalgorithm%final_time/subalgorithm%dt) * subalgorithm%dt > largest_allowed_time ) then
262 write(message(1), *) "Incommensurate propagation time: The calculation would run into a deadlock, as a system with"
263 write(message(2), *) "a smaller timestep would attempt a timestep beyond the last step" &
264 // " of a system with a larger time step."
265 write(message(3), *) "Please, adjust the TDPropagationTime and/or the TDTimeStep variables of the systems."
266 write(message(4), *) ""
267 write(message(5), *) "With the current timesteps, the TDPropagationTime should be ", largest_allowed_time
268 call messages_fatal(5, namespace=namespace)
269 end if
270 end select
271 end select
272 end do
273
274 end select
275
276 end if
277
278 ![create_interactions] !doxygen marker. Dont delete
279 ! Create and initialize interactions
280 !
281 ! This function is called recursively for all subsystems of systems.
282 ! If systems is a multisystem_basic_t container, the partners list contains all subsystems.
283 call systems%create_interactions(interactions_factory, partners)
284 ![create_interactions]
285
286 select type (systems)
287 class is (multisystem_basic_t)
288 ! Write the interaction graph as a DOT graph for debug
289 if ((debug%interaction_graph .or. debug%interaction_graph_full) .and. mpi_grp_is_root(mpi_world)) then
290 iunit_out = io_open('debug/interaction_graph.dot', systems%namespace, action='write')
291 write(iunit_out, '(a)') 'digraph {'
292 call systems%write_interaction_graph(iunit_out, debug%interaction_graph_full)
293 write(iunit_out, '(a)') '}'
294 call io_close(iunit_out)
295 end if
296 end select
297
298 if (.not. systems%process_is_slave()) then
299 call messages_write('Info: Octopus initialization completed.', new_line = .true.)
300 call messages_write('Info: Starting calculation mode.')
301 call messages_info(namespace=namespace)
302
303 !%Variable FromScratch
304 !%Type logical
305 !%Default false
306 !%Section Execution
307 !%Description
308 !% When this variable is set to true, <tt>Octopus</tt> will perform a
309 !% calculation from the beginning, without looking for restart
310 !% information.
311 !% NOTE: If available, mesh partitioning information will be used for
312 !% initializing the calculation regardless of the set value for this variable.
313 !%End
314 call parse_variable(namespace, 'FromScratch', .false., from_scratch)
315
316 call profiling_in("CALC_MODE")
317
318 select case (calc_mode_id)
319 case (option__calculationmode__gs)
320 select type (systems)
321 class is (multisystem_basic_t)
322 call multisystem_run(systems, from_scratch)
323 type is (electrons_t)
324 call ground_state_run(systems, from_scratch)
325 end select
326 case (option__calculationmode__unocc)
327 call unocc_run(systems, from_scratch)
328 case (option__calculationmode__td)
329 select type (systems)
330 class is (multisystem_basic_t)
331 call multisystem_run(systems, from_scratch)
332 type is (electrons_t)
333 call time_dependent_run(systems, from_scratch)
334 end select
335 case (option__calculationmode__go)
336 call geom_opt_run(systems, from_scratch)
337 case (option__calculationmode__opt_control)
338 call opt_control_run(systems)
339 case (option__calculationmode__em_resp)
340 select case (get_resp_method(namespace))
341 case (fd)
342 call static_pol_run(systems, from_scratch)
343 case (lr)
344 call em_resp_run(systems, from_scratch)
345 end select
346 case (option__calculationmode__casida)
347 call casida_run(systems, from_scratch)
348 case (option__calculationmode__vdw)
349 call vdw_run(systems, from_scratch)
350 case (option__calculationmode__vib_modes)
351 select case (get_resp_method(namespace))
352 case (fd)
353 call phonons_run(systems)
354 case (lr)
355 call phonons_lr_run(systems, from_scratch)
356 end select
357 case (option__calculationmode__one_shot)
358 message(1) = "CalculationMode = one_shot is obsolete. Please use gs with MaximumIter = 0."
359 call messages_fatal(1, namespace=namespace)
360 case (option__calculationmode__kdotp)
361 call kdotp_lr_run(systems, from_scratch)
362 case (option__calculationmode__dummy)
363 case (option__calculationmode__invert_ks)
364 call invert_ks_run(systems)
365 case (option__calculationmode__recipe)
366 assert(.false.) !this is handled before, if we get here, it is an error
367 end select
368
369 call profiling_out("CALC_MODE")
370 end if
371
372 select type (systems)
373 class is (multisystem_basic_t)
374 !Deallocate the external potentials
375 call iter%start(partners)
376 do while (iter%has_next())
377 select type(ptr => iter%get_next())
378 class is(external_potential_t)
379 partner => ptr
380 safe_deallocate_p(partner)
381 class is(external_waves_t)
382 partner => ptr
383 safe_deallocate_p(partner)
384 class is(lasers_t)
385 partner => ptr
386 safe_deallocate_p(partner)
387 end select
388 end do
389 end select
390
391 ! Finalize systems
392 safe_deallocate_p(systems)
393
394 call fft_all_end()
395
397
399
400 pop_sub(run)
401
402 contains
403
404 subroutine calc_mode_init()
405
406 push_sub(calc_mode_init)
407
408 select case (calc_mode_id)
409 case (option__calculationmode__gs, option__calculationmode__go, option__calculationmode__unocc)
411 case (option__calculationmode__td)
412 call td_run_init()
413 case (option__calculationmode__casida)
414 call casida_run_init()
415 end select
416
417 pop_sub(calc_mode_init)
418 end subroutine calc_mode_init
419
420 end subroutine run
421
422end module run_oct_m
423
424!! Local Variables:
425!! mode: f90
426!! coding: utf-8
427!! End:
double floor(double __x) __attribute__((__nothrow__
subroutine, public accel_init(base_grp, namespace)
Definition: accel.F90:422
subroutine, public accel_end(namespace)
Definition: accel.F90:1060
This module implements the Casida equations for excited states.
Definition: casida.F90:118
subroutine, public casida_run_init()
Definition: casida.F90:285
subroutine, public casida_run(system, from_scratch)
Definition: casida.F90:314
subroutine, public em_resp_run(system, from_scratch)
Definition: em_resp.F90:232
subroutine, public load_external_potentials(external_potentials, namespace)
subroutine, public load_external_waves(partners, namespace)
Fast Fourier Transform module. This module provides a single interface that works with different FFT ...
Definition: fft.F90:118
subroutine, public fft_all_init(namespace)
initialize the table
Definition: fft.F90:278
subroutine, public fft_all_end()
delete all plans
Definition: fft.F90:391
subroutine, public geom_opt_run(system, from_scratch)
Definition: geom_opt.F90:207
real(real64), parameter, public m_zero
Definition: global.F90:188
subroutine, public ground_state_run_init()
subroutine, public ground_state_run(electrons, from_scratch)
This module defines classes and functions for interaction partners.
subroutine, public invert_ks_run(system)
Definition: invert_ks.F90:148
subroutine, public kdotp_lr_run(system, from_scratch)
Definition: kdotp.F90:185
subroutine, public load_lasers(partners, namespace)
Definition: lasers.F90:1229
subroutine, public messages_print_with_emphasis(msg, iunit, namespace)
Definition: messages.F90:920
character(len=512), private msg
Definition: messages.F90:165
character(len=256), dimension(max_lines), public message
to be output by fatal, warning
Definition: messages.F90:160
subroutine, public messages_fatal(no_lines, only_root_writes, namespace)
Definition: messages.F90:414
subroutine, public messages_input_error(namespace, var, details, row, column)
Definition: messages.F90:713
subroutine, public messages_info(no_lines, iunit, debug_only, stress, all_nodes, namespace)
Definition: messages.F90:616
This module implements the factory for ground state algorithm.
subroutine, public mpi_debug_statistics()
Definition: mpi_debug.F90:205
logical function mpi_grp_is_root(grp)
Is the current MPI process of grpcomm, root.
Definition: mpi.F90:434
type(mpi_grp_t), public mpi_world
Definition: mpi.F90:270
This module handles the communicators for the various parallelization strategies.
Definition: multicomm.F90:145
This module implements the basic mulsisystem class, a container system for other systems.
This module implements the multisystem debug functionality.
subroutine, public multisystem_run(systems, from_scratch)
type(namespace_t), public global_namespace
Definition: namespace.F90:132
This module contains the main procedure ("opt_control_run") that is used when optimal control runs ar...
subroutine, public opt_control_run(system)
logical function, public parse_is_defined(namespace, name)
Definition: parser.F90:502
subroutine, public phonons_run(system)
Definition: phonons_fd.F90:154
subroutine, public phonons_lr_run(system, from_scratch)
Definition: phonons_lr.F90:170
subroutine, public profiling_out(label)
Increment out counter and sum up difference between entry and exit time.
Definition: profiling.F90:623
subroutine, public profiling_in(label, exclude)
Increment in counter and save entry time.
Definition: profiling.F90:552
This module implements the factory for propagators.
This module implements the basic propagator framework.
Definition: propagator.F90:117
subroutine, public pulpo_print()
Definition: pulpo.F90:129
subroutine, public restart_module_init(namespace)
Definition: restart.F90:309
top level module for all calculation modes
Definition: run.F90:117
integer function get_resp_method(namespace)
query input file for the response mode.
Definition: run.F90:175
integer, parameter fd
Definition: run.F90:168
subroutine, public run(namespace, calc_mode_id)
main routine to run all calculations: This routine parses the input file, sets up the systems and int...
Definition: run.F90:215
subroutine, public static_pol_run(system, from_scratch)
Definition: static_pol.F90:157
integer, parameter, public system_multisystem
container system. (multisystem_basic_oct_m::multisystem_basic_t)
This module implements the abstract system type.
Definition: system.F90:118
Definition: td.F90:114
subroutine, public td_run_init()
Definition: td.F90:236
This module implements a unit-test like runmode for Octopus.
Definition: test.F90:115
subroutine, public test_run(namespace)
Components and integration test runner.
Definition: test.F90:210
subroutine, public time_dependent_run(electrons, from_scratch)
This module defines the unit system, used for input and output.
subroutine, public unit_system_init(namespace)
subroutine, public unocc_run(system, from_scratch)
Definition: unocc.F90:159
subroutine, public vdw_run(system, from_scratch)
Definition: vdw.F90:152
subroutine calc_mode_init()
Definition: run.F90:498
This class defines the factory for minimizers.
Container class for lists of system_oct_m::system_t.
This class defines the factory for propagators.
Abstract class implementing propagators.
Definition: propagator.F90:138
Abstract class for systems.
Definition: system.F90:172
int true(void)