SIRIUS 7.5.0
Electronic structure library and applications
beta_projectors_gradient.hpp
Go to the documentation of this file.
1// Copyright (c) 2013-2017 Anton Kozhevnikov, Ilia Sivkov, 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 beta_projectors_gradient.hpp
21 *
22 * \brief Contains declaration and implementation of sirius::Beta_projectors_gradient class.
23 */
24
25#ifndef __BETA_PROJECTORS_GRADIENT_HPP__
26#define __BETA_PROJECTORS_GRADIENT_HPP__
27
29#include "beta_projectors.hpp"
30
31namespace sirius {
32
33/// Compute gradient of beta-projectors over atomic positions \f$ d \langle {\bf G+k} | \beta \rangle / d \tau_{\alpha} \f$.
34template <typename T>
36{
37 private:
38 void generate_pw_coefs_t(Beta_projectors<T>& beta__)
39 {
40 if (!this->num_beta_t()) {
41 return;
42 }
43
44 for (int x = 0; x < 3; x++) {
45 #pragma omp parallel for
46 for (int i = 0; i < this->num_beta_t(); i++) {
47 for (int igkloc = 0; igkloc < this->num_gkvec_loc(); igkloc++) {
48 auto vgc = this->gkvec_.template gkvec_cart<index_domain_t::local>(igkloc);
49 this->pw_coeffs_t_(igkloc, i, x) = std::complex<T>(0, -vgc[x]) * beta__.pw_coeffs_t(igkloc, i, 0);
50 }
51 }
52 }
53 }
54
55 public:
57 : Beta_projectors_base<T>(ctx__, gkvec__, 3)
58 {
59 generate_pw_coefs_t(beta__);
60 }
61};
62
63} // namespace sirius
64
65#endif // __BETA_PROJECTORS_GRADIENT_HPP__
Contains declaration and implementation of sirius::Beta_projectors class.
Contains declaration and implementation of sirius::Beta_projectors_base class.
Base class for beta-projectors, gradient of beta-projectors and strain derivatives of beta-projectors...
fft::Gvec const & gkvec_
List of G+k vectors.
sddk::mdarray< std::complex< T >, 3 > pw_coeffs_t_
Phase-factor independent coefficients of |beta> functions for atom types.
Compute gradient of beta-projectors over atomic positions .
Stores <G+k | beta> expansion.
Simulation context is a set of parameters and objects describing a single simulation.
A set of G-vectors for FFTs and G+k basis functions.
Definition: gvec.hpp:130
Namespace of the SIRIUS library.
Definition: sirius.f90:5