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
27 use casida_oct_m
31 use fft_oct_m
33 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_parallel_strategy_init(calc_mode_id)
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 unit_system_init(namespace)
151
152 call accel_init(mpi_world, namespace)
153
154 ! initialize FFTs
155 call fft_all_init(namespace)
156
157 if (calc_mode_id == option__calculationmode__test) then
158 call test_run(namespace)
159 call fft_all_end()
161 pop_sub(run)
162 return
163 end if
164
165 ! Create systems
166 if (parse_is_defined(namespace, "Systems")) then
167 ! We are running in multi-system mode
168 systems => system_factory%create(namespace, system_multisystem)
169 else
170 ! Fall back to old behaviour
171 systems => electrons_t(namespace, generate_epot = calc_mode_id /= option__calculationmode__dummy)
172 end if
173
174 ! initialize everything that needs parallelization
175 call systems%init_parallelization(mpi_world)
177 ! Create list of partners
178 select type (systems)
179 class is (multisystem_basic_t)
180 ! Systems are also partners
181 partners = systems%list
182 ! Add external potentials to partners list
183 call load_external_potentials(partners, namespace)
184
185 call load_external_waves(partners, namespace)
186
187 ! Add lasers to the partner list
188 call load_lasers(partners, namespace)
189
190 type is (electrons_t)
191 call partners%add(systems)
192 end select
193
194 ! Initialize algorithms (currently only propagators are supported)
195 select type (systems)
196 class is (multisystem_basic_t)
197 select case (calc_mode_id)
198 case (option__calculationmode__gs)
199 call systems%new_algorithm(minimizer_factory_t(systems%namespace))
200 case (option__calculationmode__td)
201 call systems%new_algorithm(propagator_factory_t(systems%namespace))
202 end select
203 end select
204
205
206 ! Check whether the final propagation time would lead to an infinite loop, and stop the code if necessary.
207 ! See issue: https://gitlab.com/octopus-code/octopus/-/issues/1226
208 !
209 ! The deadlock currently appears when systems have different timesteps, and one system could do a
210 !
211 ! 2 4 6 8 10 propagation time=11
212 ! system1 (dt=2): * * * * * |
213 ! system2 (dt=4): * * |
214 !
215 ! 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
216 ! system2 is behind, and keeps waiting for system2.
217 !
218 ! We can only allow propagations up to the largest time, the fastest system (largest dt)
219 ! can reach within the given propagation time. If one of the slower systems could perform
220 ! an extra step within the propagation time, we need to stop the calculation, before the
221 ! propagation, to prevent wasting CPU time in the deadlock.
222 !
223 ! TODO: Fix the framework, so that this infinite loop can be avoided in the first place.
224
225 if(calc_mode_id == option__calculationmode__td) then
226
227 select type (sys => systems)
228 type is (multisystem_basic_t)
229
230 largest_dt = m_zero
231 largest_allowed_time = m_zero
232
233 call sys%get_flat_list(flat_list)
234
235 call iter%start(flat_list)
236 do while (iter%has_next())
237 select type (subsystem => iter%get_next())
238 class is (system_t)
239 select type (subalgorithm => subsystem%algo)
240 class is (propagator_t)
241 largest_dt = max(largest_dt, subalgorithm%dt)
242 end select
243 end select
244 end do
245
246 select type(prop => systems%algo)
247 class is (propagator_t)
248 largest_allowed_time = floor(prop%final_time / largest_dt)* largest_dt
249 class default
250 assert(.false.)
251 end select
252
253 call iter%start(flat_list)
254 do while (iter%has_next())
255 select type(subsystem => iter%get_next())
256 class is (system_t)
257 select type (subalgorithm => subsystem%algo)
258 class is (propagator_t)
259 if( floor(subalgorithm%final_time/subalgorithm%dt) * subalgorithm%dt > largest_allowed_time ) then
260 write(message(1), *) "Incommensurate propagation time: The calculation would run into a deadlock, as a system with"
261 write(message(2), *) "a smaller timestep would attempt a timestep beyond the last step" &
262 // " of a system with a larger time step."
263 write(message(3), *) "Please, adjust the TDPropagationTime and/or the TDTimeStep variables of the systems."
264 write(message(4), *) ""
265 write(message(5), *) "With the current timesteps, the TDPropagationTime should be ", largest_allowed_time
266 call messages_fatal(5, namespace=namespace)
267 end if
268 end select
269 end select
270 end do
271
272 end select
273
274 end if
275
276 ![create_interactions] !doxygen marker. Dont delete
277 ! Create and initialize interactions
278 !
279 ! This function is called recursively for all subsystems of systems.
280 ! If systems is a multisystem_basic_t container, the partners list contains all subsystems.
281 call systems%create_interactions(interactions_factory, partners)
282 ![create_interactions]
283
284 select type (systems)
285 class is (multisystem_basic_t)
286 ! Write the interaction graph as a DOT graph for debug
287 if ((debug%interaction_graph .or. debug%interaction_graph_full) .and. mpi_world%is_root()) then
288 iunit_out = io_open('debug/interaction_graph.dot', systems%namespace, action='write')
289 write(iunit_out, '(a)') 'digraph {'
290 call systems%write_interaction_graph(iunit_out, debug%interaction_graph_full)
291 write(iunit_out, '(a)') '}'
292 call io_close(iunit_out)
293 end if
294 end select
295
296 if (.not. systems%process_is_slave()) then
297 call messages_write('Info: Octopus initialization completed.', new_line = .true.)
298 call messages_write('Info: Starting calculation mode.')
299 call messages_info(namespace=namespace)
300
301 !%Variable FromScratch
302 !%Type logical
303 !%Default false
304 !%Section Execution
305 !%Description
306 !% When this variable is set to true, <tt>Octopus</tt> will perform a
307 !% calculation from the beginning, without looking for restart
308 !% information.
309 !% NOTE: If available, mesh partitioning information will be used for
310 !% initializing the calculation regardless of the set value for this variable.
311 !%End
312 call parse_variable(namespace, 'FromScratch', .false., from_scratch)
313
314 call profiling_in("CALC_MODE")
315
316 select case (calc_mode_id)
317 case (option__calculationmode__gs)
318 select type (systems)
319 class is (multisystem_basic_t)
320 call multisystem_run(systems, from_scratch)
321 type is (electrons_t)
322 call systems%ground_state_run(from_scratch)
323 end select
324 case (option__calculationmode__unocc)
325 call unocc_run(systems, from_scratch)
326 case (option__calculationmode__td)
327 select type (systems)
328 class is (multisystem_basic_t)
329 call multisystem_run(systems, from_scratch)
330 type is (electrons_t)
331 call time_dependent_run(systems, from_scratch)
332 end select
333 case (option__calculationmode__go)
334 call geom_opt_run(systems, from_scratch)
335 case (option__calculationmode__opt_control)
336 call opt_control_run(systems)
337 case (option__calculationmode__em_resp)
338 select case (get_resp_method(namespace))
339 case (fd)
340 call static_pol_run(systems, from_scratch)
341 case (lr)
342 call em_resp_run(systems, from_scratch)
343 end select
344 case (option__calculationmode__casida)
345 call casida_run(systems, from_scratch)
346 case (option__calculationmode__vdw)
347 call vdw_run(systems, from_scratch)
348 case (option__calculationmode__vib_modes)
349 select case (get_resp_method(namespace))
350 case (fd)
351 call phonons_run(systems)
352 case (lr)
353 call phonons_lr_run(systems, from_scratch)
354 end select
355 case (option__calculationmode__one_shot)
356 message(1) = "CalculationMode = one_shot is obsolete. Please use gs with MaximumIter = 0."
357 call messages_fatal(1, namespace=namespace)
358 case (option__calculationmode__kdotp)
359 call kdotp_lr_run(systems, from_scratch)
360 case (option__calculationmode__dummy)
361 case (option__calculationmode__invert_ks)
362 call invert_ks_run(systems)
363 case (option__calculationmode__recipe)
364 assert(.false.) !this is handled before, if we get here, it is an error
365 end select
366
367 call profiling_out("CALC_MODE")
368 end if
369
370 select type (systems)
371 class is (multisystem_basic_t)
372 !Deallocate the external potentials
373 call iter%start(partners)
374 do while (iter%has_next())
375 select type(ptr => iter%get_next())
376 class is(external_potential_t)
377 partner => ptr
378 safe_deallocate_p(partner)
379 class is(external_waves_t)
380 partner => ptr
381 safe_deallocate_p(partner)
382 class is(lasers_t)
383 partner => ptr
384 safe_deallocate_p(partner)
385 end select
386 end do
387 end select
388
389 ! Finalize systems
390 safe_deallocate_p(systems)
391
392 call fft_all_end()
393
395
397
398 pop_sub(run)
399
400 end subroutine run
401
403 subroutine calc_mode_parallel_strategy_init(calc_mode_id)
404 integer, intent(in) :: calc_mode_id
405
407
408 select case (calc_mode_id)
409 case (option__calculationmode__gs, option__calculationmode__go, option__calculationmode__unocc)
410 call calc_mode_par%set_parallelization(p_strategy_states, default = .false.)
411#ifdef HAVE_SCALAPACK
412 call calc_mode_par%set_scalapack_compat()
413#endif
414
415 case (option__calculationmode__td)
416 call calc_mode_par%set_parallelization(p_strategy_states, default = .true.)
417
418 case (option__calculationmode__casida)
419 ! Pure 'other' parallelization is a bad idea. Trying to solve the Poisson equation separately on each node
420 ! consumes excessive memory and time (easily more than is available). In principle, the line below would setup
421 ! joint domain/other parallelization, but 'other' parallelization takes precedence, especially since
422 ! multicomm_init does not know the actual problem size and uses a fictitious value of 10000, making it
423 ! impossible to choose joint parallelization wisely, and generally resulting in a choice of only one domain
424 ! group. FIXME! --DAS
425
426 ! With the recent improvements, the 'other' parallelization over electron-hole
427 ! pairs works quite well now. For smaller matrices, a combination
428 ! seems to give the fastest run times. For larger matrices (more than a few
429 ! thousand entries per dimension), CasidaDistributedMatrix is needed for
430 ! the Casida matrix to fit into memory; this takes the cores from the
431 ! 'other' strategy. For very large matrices (more than 100000), it is
432 ! advisable to use only the 'other' strategy because the diagonalization
433 ! uses most of the computation time.
434 ! Thus you may want to enable this or a combination of other and domain to get
435 ! better performance - STO
436 call calc_mode_par%set_parallelization(p_strategy_other, default = .false.) ! enabled, but not default
437 call calc_mode_par%unset_parallelization(p_strategy_kpoints) ! disabled. FIXME: could be implemented.
438
439 end select
440
442
444
445end module run_oct_m
446
447!! Local Variables:
448!! mode: f90
449!! coding: utf-8
450!! End:
double floor(double __x) __attribute__((__nothrow__
subroutine, public accel_init(base_grp, namespace)
Definition: accel.F90:438
subroutine, public accel_end(namespace)
Definition: accel.F90:745
This module handles the calculation mode.
integer, parameter, public p_strategy_kpoints
parallelization in k-points
integer, parameter, public p_strategy_other
something else like e-h pairs
type(calc_mode_par_t), public calc_mode_par
Singleton instance of parallel calculation mode.
integer, parameter, public p_strategy_states
parallelization in states
This module implements the Casida equations for excited states.
Definition: casida.F90:120
subroutine, public casida_run(system, from_scratch)
Definition: casida.F90:292
subroutine, public em_resp_run(system, from_scratch)
Definition: em_resp.F90:234
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:120
subroutine, public fft_all_init(namespace)
initialize the table
Definition: fft.F90:280
subroutine, public fft_all_end()
delete all plans
Definition: fft.F90:391
subroutine, public geom_opt_run(system, from_scratch)
Definition: geom_opt.F90:213
real(real64), parameter, public m_zero
Definition: global.F90:190
This module defines classes and functions for interaction partners.
subroutine, public invert_ks_run(system)
Definition: invert_ks.F90:150
subroutine, public kdotp_lr_run(system, from_scratch)
Definition: kdotp.F90:187
subroutine, public load_lasers(partners, namespace)
Definition: lasers.F90:1230
subroutine, public messages_print_with_emphasis(msg, iunit, namespace)
Definition: messages.F90:904
character(len=512), private msg
Definition: messages.F90:167
character(len=256), dimension(max_lines), public message
to be output by fatal, warning
Definition: messages.F90:162
subroutine, public messages_fatal(no_lines, only_root_writes, namespace)
Definition: messages.F90:416
subroutine, public messages_input_error(namespace, var, details, row, column)
Definition: messages.F90:697
subroutine, public messages_info(no_lines, iunit, debug_only, stress, all_nodes, namespace)
Definition: messages.F90:600
This module implements the factory for ground state algorithm.
subroutine, public mpi_debug_statistics()
Definition: mpi_debug.F90:219
type(mpi_grp_t), public mpi_world
Definition: mpi.F90:272
This module handles the communicators for the various parallelization strategies.
Definition: multicomm.F90:147
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:134
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:505
subroutine, public phonons_run(system)
Definition: phonons_fd.F90:156
subroutine, public phonons_lr_run(system, from_scratch)
Definition: phonons_lr.F90:172
subroutine, public profiling_out(label)
Increment out counter and sum up difference between entry and exit time.
Definition: profiling.F90:625
subroutine, public profiling_in(label, exclude)
Increment in counter and save entry time.
Definition: profiling.F90:554
This module implements the factory for propagators.
This module implements the basic propagator framework.
Definition: propagator.F90:119
subroutine, public pulpo_print()
Definition: pulpo.F90:131
top level module for all calculation modes
Definition: run.F90:119
integer function get_resp_method(namespace)
query input file for the response mode.
Definition: run.F90:177
integer, parameter fd
Definition: run.F90:170
subroutine calc_mode_parallel_strategy_init(calc_mode_id)
Calculation mode initialization (parallelization strategy)
Definition: run.F90:499
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:217
subroutine, public static_pol_run(system, from_scratch)
Definition: static_pol.F90:159
integer, parameter, public system_multisystem
container system. (multisystem_basic_oct_m::multisystem_basic_t)
This module implements the abstract system type.
Definition: system.F90:120
Definition: td.F90:116
This module implements a unit-test like runmode for Octopus.
Definition: test.F90:117
subroutine, public test_run(namespace)
Components and integration test runner.
Definition: test.F90:213
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:161
subroutine, public vdw_run(system, from_scratch)
Definition: vdw.F90:154
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:140
Abstract class for systems.
Definition: system.F90:174
int true(void)