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
53 use kdotp_oct_m
55 use pulpo_oct_m
59 use system_oct_m
60 use td_oct_m
61 use test_oct_m
64 use unocc_oct_m
66 use vdw_oct_m
67
68 implicit none
69
70 private
71 public :: &
72 run
73
74 integer, parameter :: LR = 1, fd = 2
75
76contains
77
78 ! ---------------------------------------------------------
80 integer function get_resp_method(namespace)
81 type(namespace_t), intent(in) :: namespace
82
83 push_sub(get_resp_method)
84
85 !%Variable ResponseMethod
86 !%Type integer
87 !%Default sternheimer
88 !%Section Linear Response
89 !%Description
90 !% Some response properties can be calculated either via
91 !% Sternheimer linear response or by using finite
92 !% differences. You can use this variable to select how you want
93 !% them to be calculated, it applies to <tt>em_resp</tt> and <tt>vib_modes</tt>
94 !% calculation modes. By default, the Sternheimer linear-response
95 !% technique is used.
96 !%Option sternheimer 1
97 !% The linear response is obtained by solving a self-consistent
98 !% Sternheimer equation for the variation of the orbitals. This
99 !% is the recommended method.
100 !%Option finite_differences 2
101 !% Properties are calculated as a finite-differences derivative of
102 !% the energy obtained by several ground-state calculations. This
103 !% method, slow and limited only to static response, is kept
104 !% mainly because it is simple and useful for testing purposes.
105 !%End
106
107 call parse_variable(namespace, 'ResponseMethod', lr, get_resp_method)
108
109 if (.not. varinfo_valid_option('ResponseMethod', get_resp_method)) then
110 call messages_input_error(namespace, 'ResponseMethod')
111 end if
112
113 pop_sub(get_resp_method)
114 end function get_resp_method
115
116 ! ---------------------------------------------------------
120 subroutine run(namespace, calc_mode_id)
121 type(namespace_t), intent(in) :: namespace
122 integer, intent(in) :: calc_mode_id
123
124 type(partner_list_t) :: partners
125 class(system_t), pointer :: systems
126 type(system_factory_t) :: system_factory
127 type(interactions_factory_t) :: interactions_factory
128 logical :: from_scratch
129 integer :: iunit_out
130 type(partner_iterator_t) :: iter
131 class(interaction_partner_t), pointer :: partner
132
133 push_sub(run)
134
135 call messages_print_with_emphasis(msg="Calculation Mode", namespace=namespace)
136 call messages_print_var_option("CalculationMode", calc_mode_id, namespace=namespace)
137 call messages_print_with_emphasis(namespace=namespace)
138
139 call calc_mode_init()
140
141 if (calc_mode_id == option__calculationmode__recipe) then
142 call pulpo_print()
143 pop_sub(run)
144 return
145 end if
146
147 call restart_module_init(namespace)
148
149 call unit_system_init(namespace)
150
151 call accel_init(mpi_world, namespace)
152
153 ! initialize FFTs
154 call fft_all_init(namespace)
155
156 if (calc_mode_id == option__calculationmode__test) then
157 call test_run(namespace)
158 call fft_all_end()
160 pop_sub(run)
161 return
162 end if
163
164 ! Create systems
165 if (parse_is_defined(namespace, "Systems")) then
166 ! We are running in multi-system mode
167 systems => system_factory%create(namespace, system_multisystem)
168 else
169 ! Fall back to old behaviour
170 systems => electrons_t(namespace, generate_epot = calc_mode_id /= option__calculationmode__dummy)
171 end if
172
173 ! initialize everything that needs parallelization
174 call systems%init_parallelization(mpi_world)
175
176 ! Create list of partners
177 select type (systems)
178 class is (multisystem_basic_t)
179 ! Systems are also partners
180 partners = systems%list
181 ! Add external potentials to partners list
182 call load_external_potentials(partners, namespace)
183
184 call load_external_waves(partners, namespace)
185
186 ! Add lasers to the partner list
187 call load_lasers(partners, namespace)
188
189 type is (electrons_t)
190 call partners%add(systems)
191 end select
192
193 ! Initialize algorithms (currently only propagators are supported)
194 select type (systems)
195 class is (multisystem_basic_t)
196 select case (calc_mode_id)
197 case (option__calculationmode__gs)
198 call systems%new_algorithm(minimizer_factory_t(systems%namespace))
199 case (option__calculationmode__td)
200 call systems%new_algorithm(propagator_factory_t(systems%namespace))
201 end select
202 end select
203
204 ![create_interactions] !doxygen marker. Dont delete
205 ! Create and initialize interactions
206 !
207 ! This function is called recursively for all subsystems of systems.
208 ! If systems is a multisystem_basic_t container, the partners list contains all subsystems.
209 call systems%create_interactions(interactions_factory, partners)
210 ![create_interactions]
211
212 select type (systems)
214 ! Write the interaction graph as a DOT graph for debug
215 if ((debug%interaction_graph .or. debug%interaction_graph_full) .and. mpi_grp_is_root(mpi_world)) then
216 iunit_out = io_open('debug/interaction_graph.dot', systems%namespace, action='write')
217 write(iunit_out, '(a)') 'digraph {'
218 call systems%write_interaction_graph(iunit_out, debug%interaction_graph_full)
219 write(iunit_out, '(a)') '}'
220 call io_close(iunit_out)
221 end if
222 end select
223
224 if (.not. systems%process_is_slave()) then
225 call messages_write('Info: Octopus initialization completed.', new_line = .true.)
226 call messages_write('Info: Starting calculation mode.')
227 call messages_info(namespace=namespace)
228
229 !%Variable FromScratch
230 !%Type logical
231 !%Default false
232 !%Section Execution
233 !%Description
234 !% When this variable is set to true, <tt>Octopus</tt> will perform a
235 !% calculation from the beginning, without looking for restart
236 !% information.
237 !% NOTE: If available, mesh partitioning information will be used for
238 !% initializing the calculation regardless of the set value for this variable.
239 !%End
240 call parse_variable(namespace, 'FromScratch', .false., from_scratch)
241
242 call profiling_in("CALC_MODE")
243
244 select case (calc_mode_id)
245 case (option__calculationmode__gs)
246 select type (systems)
247 class is (multisystem_basic_t)
248 call multisystem_run(systems, from_scratch)
249 type is (electrons_t)
250 call ground_state_run(systems, from_scratch)
251 end select
252 case (option__calculationmode__unocc)
253 call unocc_run(systems, from_scratch)
254 case (option__calculationmode__td)
255 select type (systems)
256 class is (multisystem_basic_t)
257 call multisystem_run(systems, from_scratch)
258 type is (electrons_t)
259 call time_dependent_run(systems, from_scratch)
260 end select
261 case (option__calculationmode__go)
262 call geom_opt_run(systems, from_scratch)
263 case (option__calculationmode__opt_control)
264 call opt_control_run(systems)
265 case (option__calculationmode__em_resp)
266 select case (get_resp_method(namespace))
267 case (fd)
268 call static_pol_run(systems, from_scratch)
269 case (lr)
270 call em_resp_run(systems, from_scratch)
271 end select
272 case (option__calculationmode__casida)
273 call casida_run(systems, from_scratch)
274 case (option__calculationmode__vdw)
275 call vdw_run(systems, from_scratch)
276 case (option__calculationmode__vib_modes)
277 select case (get_resp_method(namespace))
278 case (fd)
279 call phonons_run(systems)
280 case (lr)
281 call phonons_lr_run(systems, from_scratch)
282 end select
283 case (option__calculationmode__one_shot)
284 message(1) = "CalculationMode = one_shot is obsolete. Please use gs with MaximumIter = 0."
285 call messages_fatal(1, namespace=namespace)
286 case (option__calculationmode__kdotp)
287 call kdotp_lr_run(systems, from_scratch)
288 case (option__calculationmode__dummy)
289 case (option__calculationmode__invert_ks)
290 call invert_ks_run(systems)
291 case (option__calculationmode__recipe)
292 assert(.false.) !this is handled before, if we get here, it is an error
293 end select
294
295 call profiling_out("CALC_MODE")
296 end if
297
298 select type (systems)
299 class is (multisystem_basic_t)
300 !Deallocate the external potentials
301 call iter%start(partners)
302 do while (iter%has_next())
303 select type(ptr => iter%get_next())
304 class is(external_potential_t)
305 partner => ptr
306 safe_deallocate_p(partner)
307 class is(external_waves_t)
308 partner => ptr
309 safe_deallocate_p(partner)
310 class is(lasers_t)
311 partner => ptr
312 safe_deallocate_p(partner)
313 end select
314 end do
315 end select
316
317 ! Finalize systems
318 safe_deallocate_p(systems)
319
320 call fft_all_end()
321
323
325
326 pop_sub(run)
327
328 contains
329
330 subroutine calc_mode_init()
331
332 push_sub(calc_mode_init)
333
334 select case (calc_mode_id)
335 case (option__calculationmode__gs, option__calculationmode__go, option__calculationmode__unocc)
337 case (option__calculationmode__td)
338 call td_run_init()
339 case (option__calculationmode__casida)
340 call casida_run_init()
341 end select
342
343 pop_sub(calc_mode_init)
344 end subroutine calc_mode_init
345
346 end subroutine run
347
348end module run_oct_m
349
350!! Local Variables:
351!! mode: f90
352!! coding: utf-8
353!! End:
subroutine, public accel_init(base_grp, namespace)
Definition: accel.F90:420
subroutine, public accel_end(namespace)
Definition: accel.F90:1056
This module implements the Casida equations for excited states.
Definition: casida.F90:118
subroutine, public casida_run_init()
Definition: casida.F90:284
subroutine, public casida_run(system, from_scratch)
Definition: casida.F90:313
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:277
subroutine, public fft_all_end()
delete all plans
Definition: fft.F90:390
subroutine, public geom_opt_run(system, from_scratch)
Definition: geom_opt.F90:206
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:184
subroutine, public load_lasers(partners, namespace)
Definition: lasers.F90:1224
subroutine, public messages_print_with_emphasis(msg, iunit, namespace)
Definition: messages.F90:930
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:420
subroutine, public messages_input_error(namespace, var, details, row, column)
Definition: messages.F90:723
subroutine, public messages_info(no_lines, iunit, debug_only, stress, all_nodes, namespace)
Definition: messages.F90:624
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:430
type(mpi_grp_t), public mpi_world
Definition: mpi.F90:266
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.
subroutine, public pulpo_print()
Definition: pulpo.F90:129
subroutine, public restart_module_init(namespace)
Definition: restart.F90:308
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:174
integer, parameter fd
Definition: run.F90:167
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:214
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:116
subroutine, public test_run(namespace)
Definition: test.F90:206
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:157
subroutine, public vdw_run(system, from_scratch)
Definition: vdw.F90:152
subroutine calc_mode_init()
Definition: run.F90:424
This class defines the factory for minimizers.
Container class for lists of system_oct_m::system_t.
This class defines the factory for propagators.
int true(void)