Octopus
gauge_field.F90
Go to the documentation of this file.
1!! Copyright (C) 2008 X. Andrade
2!! Copyright (C) 2020 N. Tancogne-Dejean
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
24 use debug_oct_m
25 use global_oct_m
26 use grid_oct_m
29 use, intrinsic :: iso_fortran_env
31 use mesh_oct_m
34 use mpi_oct_m
36 use parser_oct_m
40 use space_oct_m
43 use unit_oct_m
46
47 implicit none
48
49 private
50
51 public :: &
70
71 type, extends(interaction_partner_t) :: gauge_field_t
72 private
73 type(space_t), public :: space
74 real(real64), allocatable :: vecpot(:)
75 real(real64), allocatable :: vecpot_vel(:)
76 real(real64), allocatable :: vecpot_acc(:)
77 real(real64), allocatable :: vecpot_kick(:)
78 real(real64), allocatable :: force(:)
79 real(real64) :: wp2
80 logical :: with_gauge_field = .false.
81 integer :: dynamics
82 real(real64) :: kicktime
83
84 real(real64) :: volume
85
86 contains
87 procedure :: init_interaction_as_partner => gauge_field_init_interaction_as_partner
88 procedure :: copy_quantities_to_interaction => gauge_field_copy_quantities_to_interaction
90 end type gauge_field_t
91
92 interface gauge_field_t
93 module procedure gauge_field_init
94 end interface gauge_field_t
95
96
97contains
98
99 ! ---------------------------------------------------------
100 function gauge_field_init(namespace, volume) result(this)
101 class(gauge_field_t), pointer :: this
102 type(namespace_t), intent(in) :: namespace
103 real(real64), intent(in) :: volume
104
105 integer :: ii
106 type(block_t) :: blk
107
108 push_sub(gauge_field_init)
109
110 safe_allocate(this)
111
112 this%namespace = namespace_t("GaugeField", parent=namespace)
113
114 ! Temporary initialization of the space.
115 ! This will be inherited from system_t once this is turned into a system
116 this%space = space_t(namespace)
117
118 this%with_gauge_field = .false.
119
120 this%volume = volume
121
122 safe_allocate(this%vecpot(1:this%space%dim))
123 safe_allocate(this%vecpot_vel(1:this%space%dim))
124 safe_allocate(this%vecpot_acc(1:this%space%dim))
125 safe_allocate(this%vecpot_kick(1:this%space%dim))
126 safe_allocate(this%force(1:this%space%dim))
127 this%vecpot = m_zero
128 this%vecpot_vel = m_zero
129 this%vecpot_acc = m_zero
130 this%vecpot_kick = m_zero
131 this%force = m_zero
132
133 !%Variable GaugeFieldDynamics
134 !%Type integer
135 !%Default polarization
136 !%Section Hamiltonian
137 !%Description
138 !% This variable select the dynamics of the gauge field used to
139 !% apply a finite electric field to periodic systems in
140 !% time-dependent runs.
141 !%Option none 0
142 !% The gauge field does not have dynamics. The induced polarization field is zero.
143 !%Option polarization 1
144 !% The gauge field follows the dynamic described in
145 !% Bertsch et al, Phys. Rev. B 62 7998 (2000).
146 !%End
147
148 call parse_variable(namespace, 'GaugeFieldDynamics', option__gaugefielddynamics__polarization, this%dynamics)
149
150 !%Variable GaugeFieldPropagate
151 !%Type logical
152 !%Default no
153 !%Section Hamiltonian
154 !%Description
155 !% Propagate the gauge field with initial condition set by GaugeVectorField or zero if not specified
156 !%End
157
158 call parse_variable(namespace, 'GaugeFieldPropagate', .false., this%with_gauge_field)
159
160 !%Variable GaugeVectorField
161 !%Type block
162 !%Section Hamiltonian
163 !%Description
164 !% The gauge vector field is used to include a uniform (but time-dependent)
165 !% external electric field in a time-dependent run for
166 !% a periodic system. An optional second row specifies the initial
167 !% value for the time derivative of the gauge field (which is set
168 !% to zero by default). By default this field is not included.
169 !% If <tt>KPointsUseSymmetries = yes</tt>, then <tt>SymmetryBreakDir</tt>
170 !% must be set in the same direction.
171 !% This is used with utility <tt>oct-dielectric_function</tt>
172 !% according to GF Bertsch, J-I Iwata, A Rubio, and K Yabana,
173 !% <i>Phys. Rev. B</i> <b>62</b>, 7998-8002 (2000).
174 !%End
175 ! Read the initial gauge vector field
176
177 if (parse_block(namespace, 'GaugeVectorField', blk) == 0) then
178
179 this%with_gauge_field = .true.
181 do ii = 1, this%space%dim
182 call parse_block_float(blk, 0, ii - 1, this%vecpot_kick(ii))
183 end do
184
185 call parse_block_end(blk)
186 if (.not. this%space%is_periodic()) then
187 message(1) = "GaugeVectorField is intended for periodic systems."
188 call messages_warning(1, namespace=namespace)
189 end if
190
191 end if
192
193 !%Variable GaugeFieldDelay
194 !%Type float
195 !%Default 0.
196 !%Section Hamiltonian
197 !%Description
198 !% The application of the gauge field acts as a probe of the system. For dynamical
199 !% systems one can apply this probe with a delay relative to the start of the simulation.
200 !%End
201
202 call parse_variable(namespace, 'GaugeFieldDelay', m_zero, this%kicktime)
203
204 if (abs(this%kicktime) <= m_epsilon) then
205 this%vecpot(1:this%space%dim) = this%vecpot_kick(1:this%space%dim)
206 end if
207
208 pop_sub(gauge_field_init)
209 end function gauge_field_init
210
211 ! ---------------------------------------------------------
212 subroutine gauge_field_check_symmetries(this, kpoints)
213 type(gauge_field_t), intent(in) :: this
214 type(kpoints_t), intent(in) :: kpoints
215
216 integer :: iop
217
219
220 if (kpoints%use_symmetries) then
221 do iop = 1, symmetries_number(kpoints%symm)
222 if (iop == symmetries_identity_index(kpoints%symm)) cycle
223 if (.not. symm_op_invariant_cart(kpoints%symm%ops(iop), this%vecpot_kick, 1e-5_real64)) then
224 message(1) = "The GaugeVectorField breaks (at least) one of the symmetries used to reduce the k-points."
225 message(2) = "Set SymmetryBreakDir equal to GaugeVectorField."
226 call messages_fatal(2, namespace=this%namespace)
227 end if
228 end do
229 end if
230
232 end subroutine gauge_field_check_symmetries
233
234
235 ! ---------------------------------------------------------
236 subroutine gauge_field_end(this)
237 type(gauge_field_t), intent(inout) :: this
238
239 push_sub(gauge_field_end)
240
241 this%with_gauge_field = .false.
242 safe_deallocate_a(this%vecpot)
243 safe_deallocate_a(this%vecpot_vel)
244 safe_deallocate_a(this%vecpot_acc)
245 safe_deallocate_a(this%vecpot_kick)
246 safe_deallocate_a(this%force)
247
248 pop_sub(gauge_field_end)
249 end subroutine gauge_field_end
250
251
252 ! ---------------------------------------------------------
253 logical pure function gauge_field_is_propagated(this) result(is_propagated)
254 type(gauge_field_t), intent(in) :: this
255
256 is_propagated = this%with_gauge_field
257 end function gauge_field_is_propagated
258
259 ! ---------------------------------------------------------
260 logical pure function gauge_field_is_used(this) result(is_used)
261 type(gauge_field_t), intent(in) :: this
262
263 is_used = this%with_gauge_field .or. (norm2(this%vecpot_kick(1:this%space%dim)) > m_epsilon)
264 end function gauge_field_is_used
265
266
267
268 ! ---------------------------------------------------------
269 subroutine gauge_field_set_vec_pot(this, vec_pot)
270 type(gauge_field_t), intent(inout) :: this
271 real(real64), intent(in) :: vec_pot(:)
272
274 this%vecpot(1:this%space%dim) = vec_pot(1:this%space%dim)
275
277 end subroutine gauge_field_set_vec_pot
278
279
280 ! ---------------------------------------------------------
281 subroutine gauge_field_set_vec_pot_vel(this, vec_pot_vel)
282 type(gauge_field_t), intent(inout) :: this
283 real(real64), intent(in) :: vec_pot_vel(:)
284
286 this%vecpot_vel(1:this%space%dim) = vec_pot_vel(1:this%space%dim)
287
289 end subroutine gauge_field_set_vec_pot_vel
290
291
292 ! ---------------------------------------------------------
293 subroutine gauge_field_get_vec_pot(this, vec_pot)
294 type(gauge_field_t), intent(in) :: this
295 real(real64), intent(out) :: vec_pot(:)
296
298 vec_pot(1:this%space%dim) = this%vecpot(1:this%space%dim)
299
301 end subroutine gauge_field_get_vec_pot
302
303
304 ! ---------------------------------------------------------
305 subroutine gauge_field_get_vec_pot_vel(this, vec_pot_vel)
306 type(gauge_field_t), intent(in) :: this
307 real(real64), intent(out) :: vec_pot_vel(:)
308
310 vec_pot_vel(1:this%space%dim) = this%vecpot_vel(1:this%space%dim)
311
313 end subroutine gauge_field_get_vec_pot_vel
314
315
316 ! ---------------------------------------------------------
317 subroutine gauge_field_get_vec_pot_acc(this, vec_pot_acc)
318 type(gauge_field_t), intent(in) :: this
319 real(real64), intent(out) :: vec_pot_acc(:)
320
322 vec_pot_acc(1:this%space%dim) = this%vecpot_acc(1:this%space%dim)
323
325 end subroutine gauge_field_get_vec_pot_acc
326
327 ! ---------------------------------------------------------
328 subroutine gauge_field_propagate(this, dt, time)
329 type(gauge_field_t), intent(inout) :: this
330 real(real64), intent(in) :: dt
331 real(real64), intent(in) :: time
332
333 logical, save :: warning_shown = .false.
334 integer :: idim
335
336 push_sub(gauge_field_propagate)
337
338 this%vecpot_acc(1:this%space%dim) = this%force(1:this%space%dim)
339
340 ! apply kick, in case kicktime=0 the kick has already been applied
341 if (this%kicktime > m_zero .and. time-dt <= this%kicktime .and. time > this%kicktime) then
342 this%vecpot(1:this%space%dim) = this%vecpot(1:this%space%dim) + this%vecpot_kick(1:this%space%dim)
343 call messages_write(' ---------------- Applying gauge kick ----------------')
344 call messages_info(namespace=this%namespace)
345 end if
347 this%vecpot(1:this%space%dim) = this%vecpot(1:this%space%dim) + dt * this%vecpot_vel(1:this%space%dim) + &
348 m_half * dt**2 * this%force(1:this%space%dim)
349
350 !In the case of a kick, the induced field could not be higher than the initial kick
351 do idim = 1, this%space%dim
352 if (.not. warning_shown .and. abs(this%vecpot_kick(idim)) > m_epsilon .and. &
353 abs(this%vecpot(idim))> abs(this%vecpot_kick(idim))*1.01 .and. .not. this%kicktime > m_zero) then
354
355 warning_shown = .true.
356
357 write(message(1),'(a)') 'It seems that the gauge-field might be diverging. You should probably check'
358 write(message(2),'(a)') 'the simulation parameters, in particular the number of k-points.'
359 call messages_warning(2, namespace=this%namespace)
360 end if
361 end do
363 end subroutine gauge_field_propagate
364
365 ! ---------------------------------------------------------
366 subroutine gauge_field_init_vec_pot(this, qtot)
367 type(gauge_field_t), intent(inout) :: this
368 real(real64), intent(in) :: qtot
369
371
372 this%wp2 = m_four*m_pi*qtot/this%volume
373
374 write (message(1), '(a,f12.6,a)') "Info: Electron-gas plasmon frequency", &
375 units_from_atomic(units_out%energy, sqrt(this%wp2)), " ["//trim(units_abbrev(units_out%energy))//"]"
376 call messages_info(1, namespace=this%namespace)
377
379 end subroutine gauge_field_init_vec_pot
380
381 ! ---------------------------------------------------------
382 real(real64) function gauge_field_get_energy(this) result(energy)
383 type(gauge_field_t), intent(in) :: this
384
385 push_sub(gauge_field_get_energy)
387 if(allocated(this%vecpot_vel)) then
388 energy = this%volume / (8.0_real64 * m_pi * p_c**2) * sum(this%vecpot_vel(1:this%space%dim)**2)
389 else
390 energy = m_zero
391 end if
392
394 end function gauge_field_get_energy
395
396
397 ! ---------------------------------------------------------
398 subroutine gauge_field_dump(restart, gfield, ierr)
399 type(restart_t), intent(in) :: restart
400 type(gauge_field_t), intent(in) :: gfield
401 integer, intent(out) :: ierr
402
403 integer :: err
404 real(real64), allocatable :: vecpot(:,:)
405
406 push_sub(gauge_field_dump)
407
408 ierr = 0
409
410 if (restart_skip(restart)) then
411 pop_sub(gauge_field_dump)
412 return
413 end if
414
415 if (debug%info) then
416 message(1) = "Debug: Writing gauge field restart."
417 call messages_info(1, namespace=gfield%namespace)
418 end if
419
420 safe_allocate(vecpot(1:gfield%space%dim, 1:2))
421 vecpot = m_zero
422 call gauge_field_get_vec_pot(gfield, vecpot(:, 1))
423 call gauge_field_get_vec_pot_vel(gfield, vecpot(:, 2))
424
425 call drestart_write_binary(restart, "gauge_field", 2*gfield%space%dim, vecpot, err)
426 safe_deallocate_a(vecpot)
427 if (err /= 0) ierr = ierr + 1
428
429 if (debug%info) then
430 message(1) = "Debug: Writing gauge field restart done."
431 call messages_info(1, namespace=gfield%namespace)
432 end if
433
434 pop_sub(gauge_field_dump)
435 end subroutine gauge_field_dump
436
437
438 ! ---------------------------------------------------------
439 subroutine gauge_field_load(restart, gfield, ierr)
440 type(restart_t), intent(in) :: restart
441 type(gauge_field_t), intent(inout) :: gfield
442 integer, intent(out) :: ierr
443
444 integer :: err
445 real(real64), allocatable :: vecpot(:,:)
446
447 push_sub(gauge_field_load)
448
449 ierr = 0
450
451 if (restart_skip(restart)) then
452 ierr = -1
453 pop_sub(gauge_field_load)
454 return
455 end if
456
457 if (debug%info) then
458 message(1) = "Debug: Reading gauge field restart."
459 call messages_info(1, namespace=gfield%namespace)
460 end if
461
462 safe_allocate(vecpot(1:gfield%space%dim, 1:2))
463 call drestart_read_binary(restart, "gauge_field", 2*gfield%space%dim, vecpot, err)
464 if (err /= 0) ierr = ierr + 1
465
466 call gauge_field_set_vec_pot(gfield, vecpot(:,1))
467 call gauge_field_set_vec_pot_vel(gfield, vecpot(:,2))
468 safe_deallocate_a(vecpot)
469
470 if (debug%info) then
471 message(1) = "Debug: Reading gauge field restart done."
472 call messages_info(1, namespace=gfield%namespace)
473 end if
474
476 end subroutine gauge_field_load
477
478 ! ---------------------------------------------------------
479
480 subroutine gauge_field_get_force(this, gr, spin_channels, current, lrc_alpha)
481 type(gauge_field_t), intent(inout) :: this
482 type(grid_t), intent(in) :: gr
483 integer, intent(in) :: spin_channels
484 real(real64), intent(in) :: current(:,:,:)
485 real(real64), optional, intent(in) :: lrc_alpha
486
487 integer :: idir, ispin
488 real(real64) :: lrc_alpha_
489
490 push_sub(gauge_field_get_force)
492 lrc_alpha_ = optional_default(lrc_alpha, m_zero)
493
494 select case (this%dynamics)
495 case (option__gaugefielddynamics__none)
496 this%force(1:this%space%dim) = m_zero
497
498 case (option__gaugefielddynamics__polarization)
499 do idir = 1, this%space%periodic_dim
500 this%force(idir) = m_zero
501 do ispin = 1, spin_channels
502 if(lrc_alpha_ > m_epsilon) then
503 this%force(idir) = this%force(idir) + &
504 lrc_alpha*p_c/this%volume*dmf_integrate(gr, current(:, idir, ispin))
505 else
506 this%force(idir) = this%force(idir) - &
507 m_four*m_pi*p_c/this%volume*dmf_integrate(gr, current(:, idir, ispin))
508 endif
509 end do
510 end do
511
512 case default
513 assert(.false.)
514 end select
515
516 pop_sub(gauge_field_get_force)
517 end subroutine gauge_field_get_force
518
519 ! ---------------------------------------------------------
520
521 subroutine gauge_field_do_algorithmic_operation(this, operation, dt, time)
522 class(gauge_field_t), intent(inout) :: this
523 type(algorithmic_operation_t), intent(in) :: operation
524 real(real64), intent(in) :: dt
525 real(real64), intent(in) :: time
526
528
529 select case (operation%id)
530 case (verlet_start)
531 !Does nothing at the moment
532 case (verlet_finish)
533 !Does nothing at the moment
534 case (verlet_update_pos)
535 !This is inside the gauge_field_propagate routine at the moment
536 case (verlet_compute_acc)
537 call gauge_field_propagate(this, dt, time)
538
539 case (verlet_compute_vel)
540 this%vecpot_vel(1:this%space%dim) = this%vecpot_vel(1:this%space%dim) + &
541 m_half * dt * (this%vecpot_acc(1:this%space%dim) + this%force(1:this%space%dim))
542 case default
543 message(1) = "Unsupported TD operation."
544 call messages_fatal(1, namespace=this%namespace)
545 end select
546
549
550
551 ! ---------------------------------------------------------
552 subroutine gauge_field_output_write(this, out_gauge, iter)
553 type(gauge_field_t), intent(in) :: this
554 type(c_ptr), intent(inout) :: out_gauge
555 integer, intent(in) :: iter
556
557 integer :: idir
558 character(len=50) :: aux
559 real(real64) :: temp(this%space%dim)
560
561 if (.not. mpi_grp_is_root(mpi_world)) return ! only first node outputs
562
564
565 if (iter == 0) then
566 call write_iter_clear(out_gauge)
567 call write_iter_string(out_gauge,'################################################################################')
568 call write_iter_nl(out_gauge)
569 call write_iter_string(out_gauge,'# HEADER')
570 call write_iter_nl(out_gauge)
571
572 ! first line: column names
573 call write_iter_header_start(out_gauge)
574
575 do idir = 1, this%space%dim
576 write(aux, '(a2,i1,a1)') 'A(', idir, ')'
577 call write_iter_header(out_gauge, aux)
578 end do
579 do idir = 1, this%space%dim
580 write(aux, '(a6,i1,a1)') 'dA/dt(', idir, ')'
581 call write_iter_header(out_gauge, aux)
582 end do
583 do idir = 1, this%space%dim
584 write(aux, '(a10,i1,a1)') 'd^2A/dt^2(', idir, ')'
585 call write_iter_header(out_gauge, aux)
586 end do
587 call write_iter_nl(out_gauge)
588
589 ! second line: units
590 !call write_iter_string(out_gauge, '#[Iter n.]')
591 !call write_iter_header(out_gauge, '[' // trim(units_abbrev(units_out%time)) // ']')
592 !call write_iter_string(out_gauge, &
593 ! 'A Vector potential in ' // trim(units_abbrev(units_out%length)) &
594 ! 'A dot in ' // trim(units_abbrev(units_out%length)) &
595 ! 'A dot dot in ' // trim(units_abbrev(units_out%length))
596 !call write_iter_nl(out_gauge)
597
598 call write_iter_string(out_gauge,'################################################################################')
599 call write_iter_nl(out_gauge)
600
601 end if
602
603 call write_iter_start(out_gauge)
604
605 do idir = 1, this%space%dim
606 temp(idir) = units_from_atomic(units_out%energy, this%vecpot(idir))
607 end do
608 call write_iter_double(out_gauge, temp, this%space%dim)
609
610 do idir = 1, this%space%dim
611 temp(idir) = units_from_atomic(units_out%energy / units_out%time, this%vecpot_vel(idir))
612 end do
613 call write_iter_double(out_gauge, temp, this%space%dim)
615 do idir = 1, this%space%dim
616 temp(idir) = units_from_atomic(units_out%energy / units_out%time**2, this%vecpot_acc(idir))
617 end do
618 call write_iter_double(out_gauge, temp, this%space%dim)
619
620 call write_iter_nl(out_gauge)
621
623 end subroutine gauge_field_output_write
624
625 ! ---------------------------------------------------------
626 subroutine gauge_field_finalize(this)
627 type(gauge_field_t), intent(inout) :: this
628
629 push_sub(gauge_field_finalize)
630
631 call gauge_field_end(this)
632
633 pop_sub(gauge_field_finalize)
634 end subroutine gauge_field_finalize
635
636 ! ---------------------------------------------------------
637 subroutine gauge_field_init_interaction_as_partner(partner, interaction)
638 class(gauge_field_t), intent(in) :: partner
639 class(interaction_surrogate_t), intent(inout) :: interaction
640
642
643 select type (interaction)
644 class default
645 message(1) = "Unsupported interaction."
646 call messages_fatal(1)
647 end select
648
651
652 ! ---------------------------------------------------------
653 subroutine gauge_field_copy_quantities_to_interaction(partner, interaction)
654 class(gauge_field_t), intent(inout) :: partner
655 class(interaction_surrogate_t), intent(inout) :: interaction
656
658
659 select type (interaction)
660 class default
661 message(1) = "Unsupported interaction."
662 call messages_fatal(1)
663 end select
664
667
668end module gauge_field_oct_m
669
670!! Local Variables:
671!! mode: f90
672!! coding: utf-8
673!! End:
double sqrt(double __x) __attribute__((__nothrow__
This module implements the basic elements defining algorithms.
Definition: algorithm.F90:141
subroutine, public gauge_field_set_vec_pot_vel(this, vec_pot_vel)
subroutine, public gauge_field_output_write(this, out_gauge, iter)
subroutine, public gauge_field_set_vec_pot(this, vec_pot)
subroutine gauge_field_finalize(this)
subroutine, public gauge_field_get_vec_pot_acc(this, vec_pot_acc)
subroutine, public gauge_field_load(restart, gfield, ierr)
subroutine gauge_field_copy_quantities_to_interaction(partner, interaction)
subroutine, public gauge_field_end(this)
subroutine, public gauge_field_get_vec_pot(this, vec_pot)
subroutine, public gauge_field_do_algorithmic_operation(this, operation, dt, time)
subroutine, public gauge_field_check_symmetries(this, kpoints)
subroutine, public gauge_field_get_force(this, gr, spin_channels, current, lrc_alpha)
subroutine, public gauge_field_get_vec_pot_vel(this, vec_pot_vel)
class(gauge_field_t) function, pointer, public gauge_field_init(namespace, volume)
subroutine, public gauge_field_dump(restart, gfield, ierr)
subroutine gauge_field_init_interaction_as_partner(partner, interaction)
logical pure function, public gauge_field_is_propagated(this)
logical pure function, public gauge_field_is_used(this)
subroutine, public gauge_field_init_vec_pot(this, qtot)
real(real64) function, public gauge_field_get_energy(this)
subroutine gauge_field_propagate(this, dt, time)
real(real64), parameter, public m_zero
Definition: global.F90:187
real(real64), parameter, public m_epsilon
Definition: global.F90:203
This module implements the underlying real-space grid.
Definition: grid.F90:117
This module defines classes and functions for interaction partners.
This module defines various routines, operating on mesh functions.
This module defines the meshes, which are used in Octopus.
Definition: mesh.F90:118
subroutine, public messages_warning(no_lines, all_nodes, namespace)
Definition: messages.F90:543
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
integer function, public parse_block(namespace, name, blk, check_varinfo_)
Definition: parser.F90:618
integer pure function, public symmetries_identity_index(this)
Definition: symmetries.F90:602
integer pure function, public symmetries_number(this)
Definition: symmetries.F90:556
brief This module defines the class unit_t which is used by the unit_systems_oct_m module.
Definition: unit.F90:132
This module defines the unit system, used for input and output.
Explicit interfaces to C functions, defined in write_iter_low.cc.
Definition: write_iter.F90:114
abstract class for general interaction partners
int true(void)