52 logical,
public :: has_phase
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
90 if (.not.
allocated(dest))
then
93 message(1) =
"Internal error: destination batch in wfs_elec_clone_to has been previously allocated."
99 call this%copy_to(dest, pack, copy_data, new_np, special, dest_type)
101 message(1) =
"Internal error: imcompatible batches in wfs_elec_clone_to."
110 class(wfs_elec_t),
intent(in) :: this
111 class(batch_t),
allocatable,
intent(out) :: dest(:)
112 integer,
intent(in) :: n_batches
113 logical,
optional,
intent(in) :: pack
114 logical,
optional,
intent(in) :: copy_data
115 integer,
optional,
intent(in) :: new_np
116 logical,
optional,
intent(in) :: special
118 type(type_t),
optional,
intent(in) :: dest_type
124 if (.not.
allocated(dest))
then
125 safe_allocate_type_array(
wfs_elec_t, dest, (1:n_batches))
127 message(1) =
"Internal error: destination batch in wfs_elec_clone_to_array has been previously allocated."
134 call this%copy_to(dest(ib), pack, copy_data, new_np, special, dest_type)
137 message(1) =
"Internal error: incompatible batches in wfs_elec_clone_to_array."
146 subroutine wfs_elec_copy_to(this, dest, pack, copy_data, new_np, special, dest_type)
148 class(
batch_t),
intent(out) :: dest
149 logical,
optional,
intent(in) :: pack
150 logical,
optional,
intent(in) :: copy_data
151 integer,
optional,
intent(in) :: new_np
152 logical,
optional,
intent(in) :: special
153 type(
type_t),
optional,
intent(in) :: dest_type
160 dest%has_phase = this%has_phase
161 call this%batch_t%copy_to(dest%batch_t, pack, copy_data, new_np, special=special, dest_type=dest_type)
163 message(1) =
"Internal error: incompatible batches in wfs_elec_copy_to."
176 class(wfs_elec_t),
intent(in) :: this
177 class(batch_t),
intent(in) :: target
178 logical,
optional,
intent(in) :: only_check_dim
179 logical,
optional,
intent(in) :: type_check
185 assert(this%ik ==
target%ik)
186 assert(this%has_phase .eqv.
target%has_phase)
188 message(1) =
"Internal error: imcompatible batches in wfs_elec_check_compatibility_with."
191 call this%batch_t%check_compatibility_with(
target, only_check_dim, type_check)
200 class(wfs_elec_t),
intent(inout) :: this
201 logical,
optional,
intent(in) :: copy
206 this%has_phase = .false.
207 call this%batch_t%end(copy)
214#include "wfs_elec_inc.F90"
217#include "complex.F90"
218#include "wfs_elec_inc.F90"
This module implements batches of mesh functions.
character(len=256), dimension(max_lines), public message
to be output by fatal, warning
subroutine, public messages_fatal(no_lines, only_root_writes, namespace)
subroutine, public zwfs_elec_init(this, dim, st_start, st_end, np, ik, special, packed)
initialize an empty wfs_elec_t object
subroutine wfs_elec_clone_to(this, dest, pack, copy_data, new_np, special, dest_type)
clone to another wfs_elec_t object
subroutine dwfs_elec_init_with_memory_3(this, dim, st_start, st_end, psi, ik)
initialize a wfs_elec_t object with given memory
subroutine wfs_elec_clone_to_array(this, dest, n_batches, pack, copy_data, new_np, special, dest_type)
Clone the data to multiple wfs_elec_t objects.
subroutine zwfs_elec_init_with_memory_2(this, dim, st_start, st_end, psi, ik)
initialize a wfs_elec_t object with given memory
subroutine dwfs_elec_init_with_memory_2(this, dim, st_start, st_end, psi, ik)
initialize a wfs_elec_t object with given memory
subroutine wfs_elec_copy_to(this, dest, pack, copy_data, new_np, special, dest_type)
copy the data of the contained batch to another (existing) wfs_elec_t object
subroutine, public dwfs_elec_init(this, dim, st_start, st_end, np, ik, special, packed)
initialize an empty wfs_elec_t object
subroutine wfs_elec_end(this, copy)
finalze the object and the contained batch
subroutine wfs_elec_check_compatibility_with(this, target, only_check_dim, type_check)
check whether the object is compatible with a target object
subroutine zwfs_elec_init_with_memory_3(this, dim, st_start, st_end, psi, ik)
initialize a wfs_elec_t object with given memory
Class defining batches of mesh functions.
batches of electronic states