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 use types_oct_m
28
29 implicit none
30
31 private
32 public :: &
33 wfs_elec_t, &
37
45 !
46 type, extends(batch_t) :: wfs_elec_t
47 private
48 integer, public :: ik
52 logical, public :: has_phase
56 contains
57 procedure :: clone_to => wfs_elec_clone_to
58 procedure :: clone_to_array => wfs_elec_clone_to_array
59 procedure :: copy_to => wfs_elec_copy_to
60 procedure :: check_compatibility_with => wfs_elec_check_compatibility_with
61 procedure :: end => wfs_elec_end
62 end type wfs_elec_t
63
64 !--------------------------------------------------------------
65 interface wfs_elec_init
66 module procedure dwfs_elec_init_with_memory_3
67 module procedure zwfs_elec_init_with_memory_3
68 module procedure dwfs_elec_init_with_memory_2
69 module procedure zwfs_elec_init_with_memory_2
70 end interface wfs_elec_init
71
72contains
73
74 !--------------------------------------------------------------
78 subroutine wfs_elec_clone_to(this, dest, pack, copy_data, new_np, special, dest_type)
79 class(wfs_elec_t), intent(in) :: this
80 class(batch_t), allocatable, intent(out) :: dest
81 logical, optional, intent(in) :: pack
82 logical, optional, intent(in) :: copy_data
83 integer, optional, intent(in) :: new_np
84 logical, optional, intent(in) :: special
86 type(type_t), optional, intent(in) :: dest_type
87
88 push_sub(wfs_elec_clone_to)
89
90 if (.not. allocated(dest)) then
91 safe_allocate_type(wfs_elec_t, dest)
92 else
93 message(1) = "Internal error: destination batch in wfs_elec_clone_to has been previously allocated."
94 call messages_fatal(1)
95 end if
96
97 select type (dest)
98 class is (wfs_elec_t)
99 call this%copy_to(dest, pack, copy_data, new_np, special, dest_type)
100 class default
101 message(1) = "Internal error: imcompatible batches in wfs_elec_clone_to."
102 call messages_fatal(1)
103 end select
104
105 pop_sub(wfs_elec_clone_to)
106 end subroutine wfs_elec_clone_to
107
108 !--------------------------------------------------------------
110 !
111 subroutine wfs_elec_clone_to_array(this, dest, n_batches, pack, copy_data, new_np, special, dest_type)
112 class(wfs_elec_t), intent(in) :: this
113 class(batch_t), allocatable, intent(out) :: dest(:)
114 integer, intent(in) :: n_batches
115 logical, optional, intent(in) :: pack
116 logical, optional, intent(in) :: copy_data
117 integer, optional, intent(in) :: new_np
118 logical, optional, intent(in) :: special
120 type(type_t), optional, intent(in) :: dest_type
121
122 integer :: ib
123
125
126 if (.not. allocated(dest)) then
127 safe_allocate_type_array(wfs_elec_t, dest, (1:n_batches))
128 else
129 message(1) = "Internal error: destination batch in wfs_elec_clone_to_array has been previously allocated."
130 call messages_fatal(1)
131 end if
132
133 select type (dest)
134 class is (wfs_elec_t)
135 do ib = 1, n_batches
136 call this%copy_to(dest(ib), pack, copy_data, new_np, special, dest_type)
137 end do
138 class default
139 message(1) = "Internal error: imcompatible batches in wfs_elec_clone_to_array."
140 call messages_fatal(1)
141 end select
142
144 end subroutine wfs_elec_clone_to_array
146 !--------------------------------------------------------------
148 !
149 subroutine wfs_elec_copy_to(this, dest, pack, copy_data, new_np, special, dest_type)
150 class(wfs_elec_t), intent(in) :: this
151 class(batch_t), intent(out) :: dest
152 logical, optional, intent(in) :: pack
153 logical, optional, intent(in) :: copy_data
154 integer, optional, intent(in) :: new_np
155 logical, optional, intent(in) :: special
156 type(type_t), optional, intent(in) :: dest_type
157
159
160 select type (dest)
161 class is (wfs_elec_t)
162 dest%ik = this%ik
163 dest%has_phase = this%has_phase
164 call this%batch_t%copy_to(dest%batch_t, pack, copy_data, new_np, special=special, dest_type=dest_type)
165 class default
166 message(1) = "Internal error: imcompatible batches in wfs_elec_copy_to."
167 call messages_fatal(1)
168 end select
169
170 pop_sub(wfs_elec_copy_to)
171 end subroutine wfs_elec_copy_to
172
173 !--------------------------------------------------------------
178 subroutine wfs_elec_check_compatibility_with(this, target, only_check_dim, type_check)
179 class(wfs_elec_t), intent(in) :: this
180 class(batch_t), intent(in) :: target
181 logical, optional, intent(in) :: only_check_dim
182 logical, optional, intent(in) :: type_check
183
185
186 select type (target)
187 class is (wfs_elec_t)
188 assert(this%ik == target%ik)
189 assert(this%has_phase .eqv. target%has_phase)
190 class default
191 message(1) = "Internal error: imcompatible batches in wfs_elec_check_compatibility_with."
192 call messages_fatal(1)
193 end select
194 call this%batch_t%check_compatibility_with(target, only_check_dim, type_check)
195
198
199 !--------------------------------------------------------------
201 !
202 subroutine wfs_elec_end(this, copy)
203 class(wfs_elec_t), intent(inout) :: this
204 logical, optional, intent(in) :: copy
205
206 push_sub(wfs_elec_end)
207
208 this%ik = -1
209 this%has_phase = .false.
210 call this%batch_t%end(copy)
211
212 pop_sub(wfs_elec_end)
213 end subroutine wfs_elec_end
214
215
216#include "real.F90"
217#include "wfs_elec_inc.F90"
218#include "undef.F90"
219
220#include "complex.F90"
221#include "wfs_elec_inc.F90"
222#include "undef.F90"
223
224end module wfs_elec_oct_m
225
226!! Local Variables:
227!! mode: f90
228!! coding: utf-8
229!! 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:414
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:548
subroutine wfs_elec_clone_to(this, dest, pack, copy_data, new_np, special, dest_type)
clone to another wfs_elec_t object
Definition: wfs_elec.F90:172
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:358
subroutine wfs_elec_clone_to_array(this, dest, n_batches, pack, copy_data, new_np, special, dest_type)
clone the data to multipe wfs_elec_t objects
Definition: wfs_elec.F90:205
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:524
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:380
subroutine wfs_elec_copy_to(this, dest, pack, copy_data, new_np, special, dest_type)
copy the data of teh contained batch to another (existing) wfs_elec_t object
Definition: wfs_elec.F90:243
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:404
subroutine wfs_elec_end(this, copy)
finalze the object and the contained batch
Definition: wfs_elec.F90:296
subroutine wfs_elec_check_compatibility_with(this, target, only_check_dim, type_check)
check whether the object is compatible with a target object
Definition: wfs_elec.F90:272
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:502
Class defining batches of mesh functions.
Definition: batch.F90:159
batches of electronic states
Definition: wfs_elec.F90:139