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 if (.not. allocated(dest)) then
87 safe_allocate_type(wfs_elec_t, dest)
88 else
89 message(1) = "Internal error: destination batch in wfs_elec_clone_to has been previously allocated."
90 call messages_fatal(1)
91 end if
92
93 select type (dest)
94 class is (wfs_elec_t)
95 call this%copy_to(dest, pack, copy_data, new_np)
96 class default
97 message(1) = "Internal error: imcompatible batches in wfs_elec_clone_to."
98 call messages_fatal(1)
99 end select
100
101 pop_sub(wfs_elec_clone_to)
102 end subroutine wfs_elec_clone_to
103
104 !--------------------------------------------------------------
106 !
107 subroutine wfs_elec_clone_to_array(this, dest, n_batches, pack, copy_data)
108 class(wfs_elec_t), intent(in) :: this
109 class(batch_t), allocatable, intent(out) :: dest(:)
110 integer, intent(in) :: n_batches
111 logical, optional, intent(in) :: pack
112 logical, optional, intent(in) :: copy_data
113
114 integer :: ib
115
117
118 if (.not. allocated(dest)) then
119 safe_allocate_type_array(wfs_elec_t, dest, (1:n_batches))
120 else
121 message(1) = "Internal error: destination batch in wfs_elec_clone_to_array has been previously allocated."
122 call messages_fatal(1)
123 end if
124
125 select type (dest)
126 class is (wfs_elec_t)
127 do ib = 1, n_batches
128 call this%copy_to(dest(ib), pack, copy_data)
129 end do
130 class default
131 message(1) = "Internal error: imcompatible batches in wfs_elec_clone_to_array."
132 call messages_fatal(1)
133 end select
134
136 end subroutine wfs_elec_clone_to_array
137
138 !--------------------------------------------------------------
141 subroutine wfs_elec_copy_to(this, dest, pack, copy_data, new_np, special)
142 class(wfs_elec_t), intent(in) :: this
143 class(batch_t), intent(out) :: dest
144 logical, optional, intent(in) :: pack
145 logical, optional, intent(in) :: copy_data
146 integer, optional, intent(in) :: new_np
147 logical, optional, intent(in) :: special
148
151 select type (dest)
152 class is (wfs_elec_t)
153 dest%ik = this%ik
154 dest%has_phase = this%has_phase
155 call this%batch_t%copy_to(dest%batch_t, pack, copy_data, new_np, special=special)
156 class default
157 message(1) = "Internal error: imcompatible batches in wfs_elec_copy_to."
158 call messages_fatal(1)
159 end select
160
161 pop_sub(wfs_elec_copy_to)
162 end subroutine wfs_elec_copy_to
163
164 !--------------------------------------------------------------
169 subroutine wfs_elec_check_compatibility_with(this, target, only_check_dim)
170 class(wfs_elec_t), intent(in) :: this
171 class(batch_t), intent(in) :: target
172 logical, optional, intent(in) :: only_check_dim
173
175
176 select type (target)
177 class is (wfs_elec_t)
178 assert(this%ik == target%ik)
179 assert(this%has_phase .eqv. target%has_phase)
180 class default
181 message(1) = "Internal error: imcompatible batches in wfs_elec_check_compatibility_with."
182 call messages_fatal(1)
183 end select
184 call this%batch_t%check_compatibility_with(target, only_check_dim)
185
188
189 !--------------------------------------------------------------
191 !
192 subroutine wfs_elec_end(this, copy)
193 class(wfs_elec_t), intent(inout) :: this
194 logical, optional, intent(in) :: copy
195
196 push_sub(wfs_elec_end)
197
198 this%ik = -1
199 this%has_phase = .false.
200 call this%batch_t%end(copy)
201
202 pop_sub(wfs_elec_end)
203 end subroutine wfs_elec_end
204
205
206#include "real.F90"
207#include "wfs_elec_inc.F90"
208#include "undef.F90"
209
210#include "complex.F90"
211#include "wfs_elec_inc.F90"
212#include "undef.F90"
213
214end module wfs_elec_oct_m
215
216!! Local Variables:
217!! mode: f90
218!! coding: utf-8
219!! 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:160
subroutine, public messages_fatal(no_lines, only_root_writes, namespace)
Definition: messages.F90:420
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:537
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:348
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:263
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:235
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:513
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:370
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:201
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:394
subroutine wfs_elec_end(this, copy)
finalze the object and the contained batch
Definition: wfs_elec.F90:286
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:491
Class defining batches of mesh functions.
Definition: batch.F90:159
batches of electronic states
Definition: wfs_elec.F90:138