SIRIUS
7.5.0
Electronic structure library and applications
src
lapw
generate_sbessel_mt.hpp
Go to the documentation of this file.
1
// Copyright (c) 2013-2023 Anton Kozhevnikov, Thomas Schulthess
2
// All rights reserved.
3
//
4
// Redistribution and use in source and binary forms, with or without modification, are permitted provided that
5
// the following conditions are met:
6
//
7
// 1. Redistributions of source code must retain the above copyright notice, this list of conditions and the
8
// following disclaimer.
9
// 2. Redistributions in binary form must reproduce the above copyright notice, this list of conditions
10
// and the following disclaimer in the documentation and/or other materials provided with the distribution.
11
//
12
// THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS "AS IS" AND ANY EXPRESS OR IMPLIED
13
// WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A
14
// PARTICULAR PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL THE COPYRIGHT HOLDER OR CONTRIBUTORS BE LIABLE FOR
15
// ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL, EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
16
// PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR PROFITS; OR BUSINESS INTERRUPTION) HOWEVER
17
// CAUSED AND ON ANY THEORY OF LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING NEGLIGENCE OR
18
// OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
19
20
/** \file generate_sbessel_mt.hpp
21
*
22
* \brief Generate spherical Bessel functions at the muffin-tin boundary for the local set of G-vectors.
23
*/
24
25
#ifndef __GENERATE_SBESSEL_MT_HPP__
26
#define __GENERATE_SBESSEL_MT_HPP__
27
28
namespace
sirius
{
29
30
/// Compute values of spherical Bessel functions at MT boundary.
31
inline
auto
32
generate_sbessel_mt
(
Simulation_context
const
& ctx__,
int
lmax__)
33
{
34
PROFILE(
"sirius::generate_sbessel_mt"
);
35
36
sddk::mdarray<double, 3>
sbessel_mt(lmax__ + 1, ctx__.
gvec
().count(), ctx__.unit_cell().num_atom_types());
37
for
(
int
iat = 0; iat < ctx__.unit_cell().num_atom_types(); iat++) {
38
#pragma omp parallel for schedule(static)
39
for
(
int
igloc = 0; igloc < ctx__.
gvec
().count(); igloc++) {
40
auto
gv = ctx__.
gvec
().gvec_cart<
index_domain_t::local
>(igloc);
41
gsl_sf_bessel_jl_array(lmax__, gv.length() * ctx__.unit_cell().atom_type(iat).mt_radius(),
42
&sbessel_mt(0, igloc, iat));
43
}
44
}
45
return
sbessel_mt;
46
}
47
48
}
49
50
#endif
sirius::Simulation_context
Simulation context is a set of parameters and objects describing a single simulation.
Definition:
simulation_context.hpp:183
sirius::Simulation_context::gvec
auto const & gvec() const
Return const reference to Gvec object.
Definition:
simulation_context.hpp:404
sirius::sddk::mdarray< double, 3 >
sirius
Namespace of the SIRIUS library.
Definition:
sirius.f90:5
sirius::index_domain_t::local
@ local
Local index.
sirius::generate_sbessel_mt
auto generate_sbessel_mt(Simulation_context const &ctx__, int lmax__)
Compute values of spherical Bessel functions at MT boundary.
Definition:
generate_sbessel_mt.hpp:32
Generated on Wed Nov 22 2023 17:00:16 for SIRIUS by
1.9.3