SIRIUS 7.5.0
Electronic structure library and applications
generate_gvec_ylm.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_gvec_ylm.hpp
21 *
22 * \brief Generate complex spherical harmonics for the local set of G-vectors.
23 */
24
25#ifndef __GENERATE_GVEC_YLM_HPP__
26#define __GENERATE_GVEC_YLM_HPP__
27
28namespace sirius {
29
30/// Generate complex spherical harmonics for the local set of G-vectors.
31inline auto
32generate_gvec_ylm(Simulation_context const& ctx__, int lmax__)
33{
34 PROFILE("sirius::generate_gvec_ylm");
35
36 sddk::mdarray<std::complex<double>, 2> gvec_ylm(sf::lmmax(lmax__), ctx__.gvec().count(),
37 sddk::memory_t::host, "gvec_ylm");
38 #pragma omp parallel for schedule(static)
39 for (int igloc = 0; igloc < ctx__.gvec().count(); igloc++) {
40 auto rtp = r3::spherical_coordinates(ctx__.gvec().gvec_cart<index_domain_t::local>(igloc));
41 sf::spherical_harmonics(lmax__, rtp[1], rtp[2], &gvec_ylm(0, igloc));
42 }
43 return gvec_ylm;
44}
45
46}
47
48#endif
Simulation context is a set of parameters and objects describing a single simulation.
auto const & gvec() const
Return const reference to Gvec object.
Multidimensional array with the column-major (Fortran) order.
Definition: memory.hpp:660
auto spherical_coordinates(vector< double > vc)
Transform Cartesian coordinates [x,y,z] to spherical coordinates [r,theta,phi].
Definition: r3.hpp:590
int lmmax(int lmax)
Maximum number of combinations for a given .
Definition: specfunc.hpp:44
void spherical_harmonics(int lmax, double theta, double phi, std::complex< double > *ylm)
Optimized implementation of complex spherical harmonics.
Definition: specfunc.hpp:298
Namespace of the SIRIUS library.
Definition: sirius.f90:5
auto generate_gvec_ylm(Simulation_context const &ctx__, int lmax__)
Generate complex spherical harmonics for the local set of G-vectors.