Octopus
wfs_elec.F90
Go to the documentation of this file.
1!! Copyright (C) 2019 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 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
21module wfs_elec_oct_m
22 use batch_oct_m
23 use debug_oct_m
24 use global_oct_m
27
28 implicit none
29
30 private
31 public :: &
32 wfs_elec_t, &
36
44 !
45 type, extends(batch_t) :: wfs_elec_t
46 private
47 integer, public :: ik
51 logical, public :: has_phase
55 contains
56 procedure :: clone_to => wfs_elec_clone_to
57 procedure :: clone_to_array => wfs_elec_clone_to_array
58 procedure :: copy_to => wfs_elec_copy_to
59 procedure :: check_compatibility_with => wfs_elec_check_compatibility_with
60 procedure :: end => wfs_elec_end
61 end type wfs_elec_t
62
63 !--------------------------------------------------------------
64 interface wfs_elec_init
65 module procedure dwfs_elec_init_with_memory_3
66 module procedure zwfs_elec_init_with_memory_3
67 module procedure dwfs_elec_init_with_memory_2
68 module procedure zwfs_elec_init_with_memory_2
69 end interface wfs_elec_init
70
71contains
72
73 !--------------------------------------------------------------
77 subroutine wfs_elec_clone_to(this, dest, pack, copy_data, new_np)
78 class(wfs_elec_t), intent(in) :: this
79 class(batch_t), allocatable, intent(out) :: dest
80 logical, optional, intent(in) :: pack
81 logical, optional, intent(in) :: copy_data
82 integer, optional, intent(in) :: new_np
83
84 push_sub(wfs_elec_clone_to)
85
86 safe_allocate_type(wfs_elec_t, dest)
87
88 select type (dest)
89 class is (wfs_elec_t)
90 call this%copy_to(dest, pack, copy_data, new_np)
91 class default
92 message(1) = "Internal error: imcompatible batches in wfs_elec_clone_to."
93 call messages_fatal(1)
94 end select
95
96 pop_sub(wfs_elec_clone_to)
97 end subroutine wfs_elec_clone_to
98
99 !--------------------------------------------------------------
101 !
102 subroutine wfs_elec_clone_to_array(this, dest, n_batches, pack, copy_data)
103 class(wfs_elec_t), intent(in) :: this
104 class(batch_t), allocatable, intent(out) :: dest(:)
105 integer, intent(in) :: n_batches
106 logical, optional, intent(in) :: pack
107 logical, optional, intent(in) :: copy_data
108
109 integer :: ib
110
112
113 assert(n_batches > 0)
115 safe_allocate_type_array(wfs_elec_t, dest, (1:n_batches))
116
117 select type (dest)
118 class is (wfs_elec_t)
119 do ib = 1, n_batches
120 call this%copy_to(dest(ib), pack, copy_data)
121 end do
122 class default
123 message(1) = "Internal error: imcompatible batches in wfs_elec_clone_to_array."
124 call messages_fatal(1)
125 end select
126
128 end subroutine wfs_elec_clone_to_array
129
130 !--------------------------------------------------------------
132 !
133 subroutine wfs_elec_copy_to(this, dest, pack, copy_data, new_np, special)
134 class(wfs_elec_t), intent(in) :: this
135 class(batch_t), intent(out) :: dest
136 logical, optional, intent(in) :: pack
137 logical, optional, intent(in) :: copy_data
138 integer, optional, intent(in) :: new_np
139 logical, optional, intent(in) :: special
141 push_sub(wfs_elec_copy_to)
142
143 select type (dest)
144 class is (wfs_elec_t)
145 dest%ik = this%ik
146 dest%has_phase = this%has_phase
147 call this%batch_t%copy_to(dest%batch_t, pack, copy_data, new_np, special=special)
148 class default
149 message(1) = "Internal error: imcompatible batches in wfs_elec_copy_to."
151 end select
154 end subroutine wfs_elec_copy_to
155
156 !--------------------------------------------------------------
161 subroutine wfs_elec_check_compatibility_with(this, target, only_check_dim)
162 class(wfs_elec_t), intent(in) :: this
163 class(batch_t), intent(in) :: target
164 logical, optional, intent(in) :: only_check_dim
165
167
168 select type (target)
169 class is (wfs_elec_t)
170 assert(this%ik == target%ik)
171 assert(this%has_phase .eqv. target%has_phase)
172 class default
173 message(1) = "Internal error: imcompatible batches in wfs_elec_check_compatibility_with."
174 call messages_fatal(1)
175 end select
176 call this%batch_t%check_compatibility_with(target, only_check_dim)
177
180
181 !--------------------------------------------------------------
183 !
184 subroutine wfs_elec_end(this, copy)
185 class(wfs_elec_t), intent(inout) :: this
186 logical, optional, intent(in) :: copy
187
188 push_sub(wfs_elec_end)
189
190 this%ik = -1
191 this%has_phase = .false.
192 call this%batch_t%end(copy)
193
194 pop_sub(wfs_elec_end)
195 end subroutine wfs_elec_end
196
197
198#include "real.F90"
199#include "wfs_elec_inc.F90"
200#include "undef.F90"
201
202#include "complex.F90"
203#include "wfs_elec_inc.F90"
204#include "undef.F90"
205
206end module wfs_elec_oct_m
207
208!! Local Variables:
209!! mode: f90
210!! coding: utf-8
211!! End:
This module implements batches of mesh functions.
Definition: batch.F90:133
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 zwfs_elec_init(this, dim, st_start, st_end, np, ik, special, packed)
initialize an empty wfs_elec_t object
Definition: wfs_elec.F90:529
subroutine dwfs_elec_init_with_memory_3(this, dim, st_start, st_end, psi, ik)
initialize a wfs_elec_t object with given memory
Definition: wfs_elec.F90:340
subroutine wfs_elec_check_compatibility_with(this, target, only_check_dim)
check whether the object is compatible with a target object
Definition: wfs_elec.F90:255
subroutine wfs_elec_copy_to(this, dest, pack, copy_data, new_np, special)
copy the data of teh contained batch to another (existing) wfs_elec_t object
Definition: wfs_elec.F90:227
subroutine zwfs_elec_init_with_memory_2(this, dim, st_start, st_end, psi, ik)
initialize a wfs_elec_t object with given memory
Definition: wfs_elec.F90:505
subroutine dwfs_elec_init_with_memory_2(this, dim, st_start, st_end, psi, ik)
initialize a wfs_elec_t object with given memory
Definition: wfs_elec.F90:362
subroutine wfs_elec_clone_to(this, dest, pack, copy_data, new_np)
clone to another wfs_elec_t object
Definition: wfs_elec.F90:171
subroutine wfs_elec_clone_to_array(this, dest, n_batches, pack, copy_data)
clone the data to multipe wfs_elec_t objects
Definition: wfs_elec.F90:196
subroutine, public dwfs_elec_init(this, dim, st_start, st_end, np, ik, special, packed)
initialize an empty wfs_elec_t object
Definition: wfs_elec.F90:386
subroutine wfs_elec_end(this, copy)
finalze the object and the contained batch
Definition: wfs_elec.F90:278
subroutine zwfs_elec_init_with_memory_3(this, dim, st_start, st_end, psi, ik)
initialize a wfs_elec_t object with given memory
Definition: wfs_elec.F90:483
Class defining batches of mesh functions.
Definition: batch.F90:159
batches of electronic states
Definition: wfs_elec.F90:138