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