Octopus
hamiltonian_abst.F90
Go to the documentation of this file.
1!! Copyright (C) 2019 N. Tancogne-Dejean, M. Oliveira
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 st program; if not, write to the Free Software
15!! Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA
16!! 02110-1301, USA.
17!!
18#include "global.h"
19
23
25 use batch_oct_m
26 use global_oct_m
27 use mesh_oct_m
29
30 implicit none
31
32 private
33
34 public :: &
36
38 type, abstract :: hamiltonian_abst_t
40 real(real64) :: spectral_middle_point
41 real(real64) :: spectral_half_span
42
43 contains
44 procedure(is_hermitian), deferred :: is_hermitian
45 procedure(hamiltonian_update_span), deferred :: update_span
46 procedure(dhamiltonian_apply), deferred :: dapply
47 procedure(zhamiltonian_apply), deferred :: zapply
48 procedure(dhamiltonian_magnus_apply), deferred :: dmagnus_apply
49 procedure(zhamiltonian_magnus_apply), deferred :: zmagnus_apply
50 end type hamiltonian_abst_t
51
52 abstract interface
53 logical function is_hermitian(hm)
54 import
55 class(hamiltonian_abst_t), intent(in) :: hm
56 end function is_hermitian
57
58 subroutine hamiltonian_update_span(hm, delta, emin, namespace)
59 import
60 class(hamiltonian_abst_t), intent(inout) :: hm
61 real(real64), intent(in) :: delta(:)
62 real(real64), intent(in) :: emin
63 type(namespace_t), intent(in) :: namespace
64 end subroutine hamiltonian_update_span
65
66 subroutine dhamiltonian_apply(hm, namespace, mesh, psib, hpsib, terms, set_bc)
67 import
68 class(hamiltonian_abst_t), intent(in) :: hm
69 type(namespace_t), intent(in) :: namespace
70 class(mesh_t), intent(in) :: mesh
71 class(batch_t), target, intent(inout) :: psib
72 class(batch_t), target, intent(inout) :: hpsib
73 integer, optional, intent(in) :: terms
74 logical, optional, intent(in) :: set_bc
75 end subroutine dhamiltonian_apply
76
77 subroutine zhamiltonian_apply(hm, namespace, mesh, psib, hpsib, terms, set_bc)
78 import
79 class(hamiltonian_abst_t), intent(in) :: hm
80 type(namespace_t), intent(in) :: namespace
81 class(mesh_t), intent(in) :: mesh
82 class(batch_t), target, intent(inout) :: psib
83 class(batch_t), target, intent(inout) :: hpsib
84 integer, optional, intent(in) :: terms
85 logical, optional, intent(in) :: set_bc
86 end subroutine zhamiltonian_apply
87
88 subroutine dhamiltonian_magnus_apply(hm, namespace, mesh, psib, hpsib, vmagnus)
89 import
90 class(hamiltonian_abst_t), intent(in) :: hm
91 type(namespace_t), intent(in) :: namespace
92 class(mesh_t), intent(in) :: mesh
93 class(batch_t), intent(inout) :: psib
94 class(batch_t), intent(inout) :: hpsib
95 real(real64), intent(in) :: vmagnus(:, :, :)
96 end subroutine dhamiltonian_magnus_apply
97
98 subroutine zhamiltonian_magnus_apply(hm, namespace, mesh, psib, hpsib, vmagnus)
99 import
100 class(hamiltonian_abst_t), intent(in) :: hm
101 type(namespace_t), intent(in) :: namespace
102 class(mesh_t), intent(in) :: mesh
103 class(batch_t), intent(inout) :: psib
104 class(batch_t), intent(inout) :: hpsib
105 real(real64), intent(in) :: vmagnus(:, :, :)
106 end subroutine zhamiltonian_magnus_apply
107 end interface
108
109end module hamiltonian_abst_oct_m
110
111
112!! Local Variables:
113!! mode: f90
114!! coding: utf-8
115!! End:
This module implements batches of mesh functions.
Definition: batch.F90:133
This module defines an abstract class for Hamiltonians.
This module defines the meshes, which are used in Octopus.
Definition: mesh.F90:118
The abstract Hamiltonian class defines a skeleton for specific implementations.