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