Octopus
conductivity.F90
Go to the documentation of this file.
1!! Copyright (C) 2015 X. Andrade
2!!
3!! This program is free software; you can redistribute it and/or modify
4!! it under the terms of the GNU General Public License as published by
5!! the Free Software Foundation; either version 2, or (at your option)
6!! any later version.
7!!
8!! This program is distributed in the hope that it will be useful,
9!! but WITHOUT ANY WARRANTY; without even the implied warranty of
10!! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
11!! GNU General Public License for more details.
12!!
13!! You should have received a copy of the GNU General Public License
14!! along with this program; if not, write to the Free Software
15!! Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA
16!! 02110-1301, USA.
17!!
18
19#include "global.h"
20
21program conductivity
22 use batch_oct_m
24 use global_oct_m
25 use grid_oct_m
26 use io_oct_m
27 use ions_oct_m
28 use, intrinsic :: iso_fortran_env
29 use math_oct_m
32 use parser_oct_m
34 use space_oct_m
37 use smear_oct_m
38 use unit_oct_m
40
41 implicit none
42
43 integer :: iunit, ierr, ii, jj, iter, read_iter, ntime, nvel, ivel
44 integer :: istart, iend, energy_steps, out_file
45 real(real64), allocatable :: time(:), velocities(:, :)
46 real(real64), allocatable :: total_current(:, :), ftcurr(:, :, :), curr(:, :)
47 real(real64), allocatable :: heat_current(:,:), ftheatcurr(:,:,:), heatcurr(:,:)
48 complex(real64), allocatable :: invdielectric(:, :)
49 type(ions_t), pointer :: ions
50 type(space_t) :: space
51 type(spectrum_t) :: spectrum
52 type(block_t) :: blk
53 type(batch_t) :: currb, ftcurrb, heatcurrb, ftheatcurrb
54 real(real64) :: ww, curtime, deltat, v0
55 character(len=MAX_PATH_LEN) :: ref_filename
56 integer :: ref_file, time_steps_ref, kk
57 real(real64), allocatable :: velcm(:), vel0(:), current(:), current_ref(:, :)
58 real(real64) :: dt_ref, tt, start_time
59 integer :: ifreq, max_freq
60 integer :: skip, smearing
61 logical :: from_forces
62 character(len=120) :: header
63 real(real64) :: excess_charge, qtot, javerage
64
65 ! Initialize stuff
67
68 call getopt_init(ierr)
69 call getopt_end()
70
71 call parser_init()
72
73 call messages_init()
74
75 call messages_experimental('oct-conductivity')
76
77 call io_init()
79
81
83
84 call spectrum_init(spectrum, global_namespace, default_energy_step = 0.0001_real64, default_max_energy = m_one)
85
86 !%Variable ConductivitySpectrumTimeStepFactor
87 !%Type integer
88 !%Default 1
89 !%Section Utilities::oct-conductivity_spectrum
90 !%Description
91 !% In the calculation of the conductivity, it is not necessary
92 !% to read the velocity at every time step. This variable controls
93 !% the integer factor between the simulation time step and the
94 !% time step used to calculate the conductivity.
95 !%End
96
97 call messages_obsolete_variable(global_namespace, 'PropagationSpectrumTimeStepFactor', 'ConductivitySpectrumTimeStepFactor')
98 call parse_variable(global_namespace, 'ConductivitySpectrumTimeStepFactor', 1, skip)
99 if (skip <= 0) call messages_input_error(global_namespace, 'ConductivitySpectrumTimeStepFactor')
100
101 !%Variable ConductivityFromForces
102 !%Type logical
103 !%Default no
104 !%Section Utilities::oct-conductivity_spectrum
105 !%Description
106 !% (Experimental) If enabled, Octopus will attempt to calculate the conductivity from the forces instead of the current.
107 !%End
108 call parse_variable(global_namespace, 'ConductivityFromForces', .false., from_forces)
109 if (from_forces) call messages_experimental('ConductivityFromForces')
110
111 max_freq = spectrum_nenergy_steps(spectrum)
112
113 if (spectrum%end_time < m_zero) spectrum%end_time = huge(spectrum%end_time)
115 ions => ions_t(global_namespace)
116
117 !We need the total charge
118 call parse_variable(global_namespace, 'ExcessCharge', m_zero, excess_charge)
119 qtot = -(ions%val_charge() + excess_charge)
120
121 if (from_forces) then
122
123 call messages_write('Info: Reading coordinates from td.general/coordinates')
124 call messages_info()
125
126 ! Opens the coordinates files.
127 iunit = io_open('td.general/coordinates', global_namespace, action='read')
128
129 call io_skip_header(iunit)
130 call spectrum_count_time_steps(global_namespace, iunit, ntime, deltat)
131 call io_close(iunit)
132
133 nvel = ions%natoms*ions%space%dim
134
135 safe_allocate(time(1:ntime))
136 safe_allocate(velocities(1:nvel, 1:ntime))
137
138 ! Opens the coordinates files.
139 iunit = io_open('td.general/coordinates', global_namespace, action='read', status='old', die=.false.)
140
141 call io_skip_header(iunit)
142
143 ntime = 1
144 iter = 1
145
146 do
147 read(unit = iunit, iostat = ierr, fmt = *) read_iter, curtime, &
148 ((ions%pos(jj, ii), jj = 1, 3), ii = 1, ions%natoms), &
149 ((ions%vel(jj, ii), jj = 1, 3), ii = 1, ions%natoms)
150
151 curtime = units_to_atomic(units_out%time, curtime)
152
153 if (ierr /= 0 .or. curtime >= spectrum%end_time) then
154 iter = iter - 1 ! last iteration is not valid
155 ntime = ntime - 1
156 exit
157 end if
158
159 assert(iter == read_iter + 1)
160
161 if (curtime >= spectrum%start_time .and. mod(iter, skip) == 0) then
162
163 time(ntime) = curtime
164 ivel = 1
165 do ii = 1, ions%natoms
166 do jj = 1, ions%space%dim
167 velocities(ivel, ntime) = units_to_atomic(units_out%velocity, ions%vel(jj, ii))
168 ivel = ivel + 1
169 end do
170 end do
171
172 ntime = ntime + 1
173 end if
174
175 iter = iter + 1 !counts number of timesteps (with time larger than zero up to SpecEndTime)
176 end do
177
178 call io_close(iunit)
179
180 call messages_write(' done.')
181 call messages_info()
182
183 else !from_forces
184
185 call messages_write('Info: Reading total current from td.general/total_current')
186 call messages_info()
187
188 iunit = io_open('td.general/total_current', global_namespace, action='read', die=.false.)
189
190 if (iunit /= -1) then
191
192 call io_skip_header(iunit)
193 call spectrum_count_time_steps(global_namespace, iunit, ntime, deltat)
194 deltat = units_to_atomic(units_out%time, deltat)
195 ntime = ntime + 1
196
197 call io_skip_header(iunit)
198
199 safe_allocate(total_current(1:space%dim, 1:ntime))
200 safe_allocate(heat_current(1:space%dim, 1:ntime))
201 safe_allocate(time(1:ntime))
202
203 do iter = 1, ntime
204 read(iunit, *) read_iter, time(iter), (total_current(ii, iter), ii = 1, space%dim)
205 time(iter) = units_to_atomic(units_out%time, time(iter))
206 do ii = 1, space%dim
207 total_current(ii, iter) = units_to_atomic(units_out%length/units_out%time, total_current(ii, iter))
208 end do
209 end do
210
211 call io_close(iunit)
212
213 if (parse_is_defined(global_namespace, 'TransientAbsorptionReference')) then
214 call parse_variable(global_namespace, 'TransientAbsorptionReference', '.', ref_filename)
215 ref_file = io_open(trim(ref_filename)//'/total_current', action='read', status='old')
216 call io_skip_header(ref_file)
217 call spectrum_count_time_steps(global_namespace, ref_file, time_steps_ref, dt_ref)
218 dt_ref = units_to_atomic(units_out%time, dt_ref)
219 if (time_steps_ref < ntime) then
220 message(1) = "The reference calculation does not contain enought time steps"
221 call messages_fatal(1)
222 end if
223
224 if (.not. is_close(dt_ref, time(2)-time(1))) then
225 message(1) = "The time step of the reference calculation is different from the current calculation"
226 call messages_fatal(1)
227 end if
228
229 !We remove the reference
230 time_steps_ref = time_steps_ref + 1
231 safe_allocate(current_ref(1:space%dim, 1:time_steps_ref))
232 call io_skip_header(ref_file)
233 do ii = 1, time_steps_ref
234 read(ref_file, *) jj, tt, (current_ref(kk, ii), kk = 1, space%dim)
235 end do
236 call io_close(ref_file)
237 do iter = 1, ntime
238 do kk = 1, space%dim
239 total_current(kk, iter) = total_current(kk, iter) - current_ref(kk, iter)
240 end do
241 end do
242 safe_deallocate_a(current_ref)
243
244 start_time = spectrum%start_time
245 call parse_variable(global_namespace, 'GaugeFieldDelay', start_time, spectrum%start_time)
246 end if
247
248 !We substract the average value, which can cause a spurious divergence at \omega=0
249 !Only valid for non-metallic systems
250 call parse_variable(global_namespace, 'SmearingFunction', smear_semiconductor, smearing)
251 if(smearing == smear_semiconductor) then
252 do kk = 1, space%dim
253 javerage = m_zero
254 do iter = 1, ntime
255 javerage = javerage + total_current(kk, iter)
256 end do
257 javerage = javerage/ntime
258 do iter = 1, ntime
259 total_current(kk, iter) = total_current(kk, iter) - javerage
260 end do
261 end do
262 end if
263
264 else
265
266 call messages_write("Cannot find the 'td.general/total_current' file.")
267 call messages_new_line()
268 call messages_write("Optical conductivity cannot be calculated.")
269 call messages_warning()
270
271 ntime = 1
272 safe_allocate(total_current(1:space%dim, 1:ntime))
273 safe_allocate(heat_current(1:space%dim, 1:ntime))
274 safe_allocate(time(1:ntime))
275
276 total_current(1:space%dim, 1:ntime) = m_zero
277
278 end if
279
280 iunit = io_open('td.general/total_heat_current', global_namespace, action='read', status='old', die=.false.)
281
282 if (iunit /= -1) then
283
284 if (ntime == 1) then
285 call io_skip_header(iunit)
286 call spectrum_count_time_steps(global_namespace, iunit, ntime, deltat)
287 ntime = ntime + 1
288 end if
289
290 call io_skip_header(iunit)
291
292 do iter = 1, ntime
293 read(iunit, *) read_iter, time(iter), (heat_current(ii, iter), ii = 1, space%dim)
294 end do
295
296 call io_close(iunit)
297
298 else
299
300 call messages_write("Cannot find the 'td.general/heat_current' file.")
301 call messages_new_line()
302 call messages_write("Thermal conductivity cannot be calculated.")
303 call messages_new_line()
304 call messages_write("Check ConductivityFromForces to possibly compute it from forces.")
305 call messages_warning()
306
307 heat_current(1:space%dim, 1:ntime) = m_zero
308
309 end if
310 end if
311
312
313 safe_allocate(curr(ntime, 1:space%dim))
314 safe_allocate(heatcurr(ntime, 1:space%dim))
315 safe_allocate(vel0(1:space%dim))
316 vel0 = m_zero
317 safe_allocate(velcm(1:space%dim))
318 safe_allocate(current(1:space%dim))
319
320 if (from_forces) iunit = io_open('td.general/current_from_forces', global_namespace, action='write')
321
322 do iter = 1, ntime
323
324 if (from_forces) then
325 if (iter == 1) then
326 vel0 = m_zero
327 ivel = 1
328 do ii = 1, ions%natoms
329 do jj = 1, space%dim
330 vel0(jj) = vel0(jj) + velocities(ivel, iter)/real(ions%natoms, real64)
331 ivel = ivel + 1
332 end do
333 end do
334 end if
335
336 velcm = m_zero
337 current = m_zero
338
339 ivel = 1
340 do ii = 1, ions%natoms
341 do jj = 1, space%dim
342 velcm(jj) = velcm(jj) + velocities(ivel, iter)/real(ions%natoms, real64)
343 current(jj) = current(jj) + &
344 ions%mass(ii)/ions%latt%rcell_volume*(velocities(ivel, iter) - vel0(jj))
345 ivel = ivel + 1
346 end do
347 end do
348
349 curr(iter, 1:space%dim) = vel0(1:space%dim)*qtot/ions%latt%rcell_volume + current(1:space%dim)
350
351 else
352 curr(iter, 1:space%dim) = total_current(1:space%dim, iter)/ions%latt%rcell_volume
353 heatcurr(iter,1:space%dim) = heat_current(1:space%dim, iter)/ions%latt%rcell_volume
354 end if
355
356 if (from_forces) write(iunit,*) iter, iter*deltat, curr(iter, 1:space%dim)
357
358 end do
359
360 safe_deallocate_a(velcm)
361 safe_deallocate_a(current)
362
363
364 if (from_forces) call io_close(iunit)
365
366 ! Find out the iteration numbers corresponding to the time limits.
367 call spectrum_fix_time_limits(spectrum, ntime, deltat, istart, iend, max_freq)
368 istart = max(1, istart)
369 energy_steps = spectrum_nenergy_steps(spectrum)
370
371 safe_allocate(ftcurr(1:energy_steps, 1:space%dim, 1:2))
372
373 call batch_init(currb, 1, 1, space%dim, curr)
374 call spectrum_signal_damp(spectrum%damp, spectrum%damp_factor, istart, iend, spectrum%start_time, deltat, currb)
375
376 call batch_init(ftcurrb, 1, 1, space%dim, ftcurr(:, :, 1))
377 call spectrum_fourier_transform(spectrum%method, spectrum_transform_cos, spectrum%noise, &
378 istart, iend, spectrum%start_time, deltat, currb, spectrum%min_energy, spectrum%max_energy, &
379 spectrum%energy_step, ftcurrb)
380 call ftcurrb%end()
381
382 call batch_init(ftcurrb, 1, 1, space%dim, ftcurr(:, :, 2))
383 call spectrum_fourier_transform(spectrum%method, spectrum_transform_sin, spectrum%noise, &
384 istart, iend, spectrum%start_time, deltat, currb, spectrum%min_energy, spectrum%max_energy, &
385 spectrum%energy_step, ftcurrb)
386
387 call ftcurrb%end()
388 call currb%end()
389
390
391 !and print the spectrum
392 iunit = io_open('td.general/conductivity', global_namespace, action='write')
393
394
395 write(unit = iunit, iostat = ierr, fmt = '(a)') &
396 '###########################################################################################################################'
397 write(unit = iunit, iostat = ierr, fmt = '(8a)') '# HEADER'
398 write(unit = iunit, iostat = ierr, fmt = '(a,a,a)') &
399 '# Energy [', trim(units_abbrev(units_out%energy)), '] Conductivity [a.u.] ReX ImX ReY ImY ReZ ImZ'
400 write(unit = iunit, iostat = ierr, fmt = '(a)') &
401 '###########################################################################################################################'
402
403 v0 = norm2(vel0(1:space%dim))
404 if (.not. from_forces .or. v0 < epsilon(v0)) v0 = m_one
405
406 if (.not. from_forces .and. parse_is_defined(global_namespace, 'GaugeVectorField')) then
407 if (parse_block(global_namespace, 'GaugeVectorField', blk) == 0) then
408 do ii = 1, space%dim
409 call parse_block_float(blk, 0, ii - 1, vel0(ii))
410 end do
411 call parse_block_end(blk)
412 end if
413 vel0(1:space%dim) = vel0(1:space%dim) / p_c
414 v0 = norm2(vel0(1:space%dim))
415 end if
416
417
418 do ifreq = 1, energy_steps
419 ww = spectrum%energy_step*(ifreq - 1) + spectrum%min_energy
420 write(unit = iunit, iostat = ierr, fmt = '(7e20.10)') units_from_atomic(units_out%energy, ww), &
421 transpose(ftcurr(ifreq, 1:space%dim, 1:2)/v0)
422 end do
423
424 call io_close(iunit)
425
426
427 !Compute the inverse dielectric function from the conductivity
428 ! We have \chi = -i \sigma / \omega
429 ! and \epsilon^-1 = 1 + 4 \pi \chi
430 safe_allocate(invdielectric(1:space%dim, 1:energy_steps))
431 do ifreq = 1, energy_steps
432 ww = max((ifreq-1)*spectrum%energy_step + spectrum%min_energy, m_epsilon)
433
434 invdielectric(1:space%dim, ifreq) = (vel0(1:space%dim) + m_four * m_pi * &
435 cmplx(ftcurr(ifreq, 1:space%dim, 2),-ftcurr(ifreq, 1:space%dim, 1), real64) / ww) / v0
436
437 end do
438
439 out_file = io_open('td.general/inverse_dielectric_function_from_current', global_namespace, action='write')
440 select case (space%dim)
441 case (1)
442 write(header, '(3a15)') '# energy', 'Re x', 'Im x'
443 case (2)
444 write(header, '(5a15)') '# energy', 'Re x', 'Im x', 'Re y', 'Im y'
445 case (3)
446 write(header, '(7a15)') '# energy', 'Re x', 'Im x', 'Re y', 'Im y', 'Re z', 'Im z'
447 end select
448 write(out_file,'(a)') trim(header)
449 do ifreq = 1, energy_steps
450 ww = (ifreq-1)*spectrum%energy_step + spectrum%min_energy
451 select case (space%dim)
452 case (1)
453 write(out_file, '(3e15.6)') ww, &
454 real(invdielectric(1, ifreq), real64), aimag(invdielectric(1, ifreq))
455 case (2)
456 write(out_file, '(5e15.6)') ww, &
457 real(invdielectric(1, ifreq), real64), aimag(invdielectric(1, ifreq)), &
458 real(invdielectric(2, ifreq), real64), aimag(invdielectric(2, ifreq))
459 case (3)
460 write(out_file, '(7e15.6)') ww, &
461 real(invdielectric(1, ifreq), real64), aimag(invdielectric(1, ifreq)), &
462 real(invdielectric(2, ifreq), real64), aimag(invdielectric(2, ifreq)), &
463 real(invdielectric(3, ifreq), real64), aimag(invdielectric(3, ifreq))
464 end select
465 end do
466 call io_close(out_file)
467
468 safe_deallocate_a(ftcurr)
469 safe_deallocate_a(invdielectric)
470
471
472!!!!!!!!!!!!!!!!!!!!!
473 if (from_forces) then
474 safe_allocate(ftheatcurr(1:energy_steps, 1:space%dim, 1:2))
475
476 ftheatcurr = m_one
477
478 call batch_init(heatcurrb, 1, 1, space%dim, heatcurr)
479
480 call spectrum_signal_damp(spectrum%damp, spectrum%damp_factor, 1, ntime, m_zero, deltat, heatcurrb)
481
482 call batch_init(ftheatcurrb, 1, 1, space%dim, ftheatcurr(:, :, 1))
483 call spectrum_fourier_transform(spectrum%method, spectrum_transform_cos, spectrum%noise, &
484 1, ntime, m_zero, deltat, heatcurrb, spectrum%min_energy, spectrum%max_energy, spectrum%energy_step, ftheatcurrb)
485 call ftheatcurrb%end()
486
487 call batch_init(ftheatcurrb, 1, 1, space%dim, ftheatcurr(:, :, 2))
488 call spectrum_fourier_transform(spectrum%method, spectrum_transform_sin, spectrum%noise, &
489 1, ntime, m_zero, deltat, heatcurrb, spectrum%min_energy, spectrum%max_energy, spectrum%energy_step, ftheatcurrb)
490 call ftheatcurrb%end()
491
492 call heatcurrb%end()
493
494
495 !and print the spectrum
496 iunit = io_open('td.general/heat_conductivity', global_namespace, action='write')
497
498
499 write(unit = iunit, iostat = ierr, fmt = '(a)') &
500 '###########################################################################################################################'
501 write(unit = iunit, iostat = ierr, fmt = '(8a)') '# HEADER'
502 write(unit = iunit, iostat = ierr, fmt = '(a,a,a)') &
503 '# Energy [', trim(units_abbrev(units_out%energy)), '] Conductivity [a.u.] ReX ImX ReY ImY ReZ ImZ'
504 write(unit = iunit, iostat = ierr, fmt = '(a)') &
505 '###########################################################################################################################'
506
507 v0 = norm2(vel0(1:space%dim))
508 if (.not. from_forces .or. v0 < epsilon(v0)) v0 = m_one
509
510 do ifreq = 1, energy_steps
511 ww = spectrum%energy_step*(ifreq - 1) + spectrum%min_energy
512 write(unit = iunit, iostat = ierr, fmt = '(7e20.10)') units_from_atomic(units_out%energy, ww), &
513 transpose(ftheatcurr(ifreq, 1:space%dim, 1:2)/v0)
514 !print *, ifreq, ftheatcurr(ifreq, 1:3, 1:2)
515 end do
516
517 call io_close(iunit)
518
519 safe_deallocate_a(ftheatcurr)
520 end if
521!!!!!!!!!!!!!!!!!!!!!!!!!!!
522
523
524 safe_deallocate_a(vel0)
525 safe_deallocate_p(ions)
526
527 safe_deallocate_a(total_current)
528 safe_deallocate_a(heat_current)
529 safe_deallocate_a(curr)
530 safe_deallocate_a(heatcurr)
531 safe_deallocate_a(time)
532
534 call io_end()
535 call messages_end()
536
537 call parser_end()
538 call global_end()
539
540end program conductivity
541
542!! Local Variables:
543!! mode: f90
544!! coding: utf-8
545!! End:
program conductivity
initialize a batch with existing memory
Definition: batch.F90:267
This module implements batches of mesh functions.
Definition: batch.F90:133
subroutine, public getopt_init(ierr)
Initializes the getopt machinery. Must be called before attempting to parse the options....
subroutine, public getopt_end
subroutine, public global_end()
Finalise parser varinfo file, and MPI.
Definition: global.F90:382
real(real64), parameter, public m_zero
Definition: global.F90:188
real(real64), parameter, public m_four
Definition: global.F90:192
real(real64), parameter, public m_pi
some mathematical constants
Definition: global.F90:186
type(mpi_comm), parameter, public serial_dummy_comm
Alias MPI_COMM_UNDEFINED for the specific use case of initialising Octopus utilities with no MPI supp...
Definition: global.F90:265
subroutine, public init_octopus_globals(comm)
Initialise Octopus-specific global constants and files. This routine performs no initialisation calls...
Definition: global.F90:354
real(real64), parameter, public m_epsilon
Definition: global.F90:204
real(real64), parameter, public p_c
Electron gyromagnetic ratio, see Phys. Rev. Lett. 130, 071801 (2023)
Definition: global.F90:224
real(real64), parameter, public m_one
Definition: global.F90:189
This module implements the underlying real-space grid.
Definition: grid.F90:117
Definition: io.F90:114
subroutine, public io_init(defaults)
If the argument defaults is present and set to true, then the routine will not try to read anything f...
Definition: io.F90:161
subroutine, public io_close(iunit, grp)
Definition: io.F90:418
subroutine, public io_skip_header(iunit)
Definition: io.F90:597
subroutine, public io_end()
Definition: io.F90:258
integer function, public io_open(file, namespace, action, status, form, position, die, recl, grp)
Definition: io.F90:352
This module is intended to contain "only mathematical" functions and procedures.
Definition: math.F90:115
subroutine, public messages_end()
Definition: messages.F90:278
subroutine, public messages_init(output_dir)
Definition: messages.F90:225
subroutine, public messages_warning(no_lines, all_nodes, namespace)
Definition: messages.F90:538
subroutine, public messages_obsolete_variable(namespace, name, rep)
Definition: messages.F90:1046
subroutine, public messages_new_line()
Definition: messages.F90:1135
character(len=256), dimension(max_lines), public message
to be output by fatal, warning
Definition: messages.F90:161
subroutine, public messages_fatal(no_lines, only_root_writes, namespace)
Definition: messages.F90:415
subroutine, public messages_input_error(namespace, var, details, row, column)
Definition: messages.F90:714
subroutine, public messages_experimental(name, namespace)
Definition: messages.F90:1086
subroutine, public messages_info(no_lines, iunit, debug_only, stress, all_nodes, namespace)
Definition: messages.F90:617
type(namespace_t), public global_namespace
Definition: namespace.F90:132
logical function, public parse_is_defined(namespace, name)
Definition: parser.F90:502
subroutine, public parser_init()
Initialise the Octopus parser.
Definition: parser.F90:450
subroutine, public parser_end()
End the Octopus parser.
Definition: parser.F90:481
integer function, public parse_block(namespace, name, blk, check_varinfo_)
Definition: parser.F90:618
subroutine, public profiling_end(namespace)
Definition: profiling.F90:413
subroutine, public profiling_init(namespace)
Create profiling subdirectory.
Definition: profiling.F90:255
integer, parameter, public smear_semiconductor
Definition: smear.F90:171
subroutine, public spectrum_fix_time_limits(spectrum, time_steps, dt, istart, iend, ntiter)
Definition: spectrum.F90:2507
subroutine, public spectrum_fourier_transform(method, transform, noise, time_start, time_end, t0, time_step, time_function, energy_start, energy_end, energy_step, energy_function)
Computes the sine, cosine, (or "exponential") Fourier transform of the real function given in the tim...
Definition: spectrum.F90:2637
subroutine, public spectrum_init(spectrum, namespace, default_energy_step, default_max_energy)
Definition: spectrum.F90:213
subroutine, public spectrum_signal_damp(damp_type, damp_factor, time_start, time_end, t0, time_step, time_function)
Definition: spectrum.F90:2551
integer, parameter, public spectrum_transform_cos
Definition: spectrum.F90:171
integer, parameter, public spectrum_transform_sin
Definition: spectrum.F90:171
subroutine, public spectrum_count_time_steps(namespace, iunit, time_steps, dt)
Definition: spectrum.F90:2385
pure integer function, public spectrum_nenergy_steps(spectrum)
Definition: spectrum.F90:2946
brief This module defines the class unit_t which is used by the unit_systems_oct_m module.
Definition: unit.F90:132
character(len=20) pure function, public units_abbrev(this)
Definition: unit.F90:223
This module defines the unit system, used for input and output.
type(unit_system_t), public units_out
subroutine, public unit_system_init(namespace)