DFT-FE 1.1.0-pre
Density Functional Theory With Finite-Elements
Loading...
Searching...
No Matches
AtomCenteredSphericalFunctionContainer.h
Go to the documentation of this file.
1// ---------------------------------------------------------------------
2//
3// Copyright (c) 2017-2025 The Regents of the University of Michigan and DFT-FE
4// authors.
5//
6// This file is part of the DFT-FE code.
7//
8// The DFT-FE code is free software; you can use it, redistribute
9// it, and/or modify it under the terms of the GNU Lesser General
10// Public License as published by the Free Software Foundation; either
11// version 2.1 of the License, or (at your option) any later version.
12// The full text of the license can be found in the file LICENSE at
13// the top level of the DFT-FE distribution.
14//
15// ---------------------------------------------------------------------
16//
17// @author Vishal Subramanian, Kartick Ramakrishnan, Sambit Das
18//
19
20#ifndef DFTFE_ATOMCENTEREDSPHERICALFUNCTIONCONTAINERBASE_H
21#define DFTFE_ATOMCENTEREDSPHERICALFUNCTIONCONTAINERBASE_H
22
23#include "vector"
24#include "map"
26#include <memory>
27#include <MemorySpaceType.h>
28#include "FEBasisOperations.h"
29#include <headers.h>
30#include <TypeConfig.h>
31#include <dftUtils.h>
32
33
34namespace dftfe
35{
37 {
38 public:
39 /**
40 * @brief Initialises the class with the atomicNumbers of various atoms and the AtomCenteredSphericalFn of various spherical functions. This function is only called once per run.
41 * @param[in] atomicNumbers vector of size Natoms storing the Znumbers of
42 * various atoms present
43 * @param[in] listOfSphericalFunctions map of std::pain (Znum, l) to the
44 * sphericalFUnction class shared pointer.
45 */
46 void
47 init(const std::vector<dftfe::uInt> &atomicNumbers,
48 const std::map<std::pair<dftfe::uInt, dftfe::uInt>,
49 std::shared_ptr<AtomCenteredSphericalFunctionBase>>
50 &listOfSphericalFunctions);
51 /**
52 * @brief Initialises the position of atoms, the image posisiton and image ids after every update of atom positions.
53 * @param[in] atomCoords vector of size 3*Natoms storing the X,Y,Z
54 * coordiantes of atom in cell.
55 * @param[in] periodicCoords vector of vector storing the image coordinates
56 * @param[in] imageIds the image Id of image atoms present in periodicCoords
57 * input
58 */
59 void
60 initaliseCoordinates(const std::vector<double> &atomCoords,
61 const std::vector<std::vector<double>> &periodicCoords,
62 const std::vector<dftfe::Int> &imageIds);
63 /**
64 * @brief Returns the number of atoms present in domain
65 * @return Returns size of atomicNumbers vector
66 */
69
70
71 /**
72 * @brief Returns the cooridnates of atom present in domain
73 * @return Returns atomCoords vector
74 */
75 const std::vector<double> &
77 /**
78 * @brief Returns the map of atomId vs vector of image coordinates
79 * @return Returns d_periodicImageCoord
80 */
81 const std::map<dftfe::uInt, std::vector<double>> &
83
84 // This functions returns the number of spherical functions associated with
85 // an atomic number.
86 // If the atomic number does not exist, it returns a zero.
87 /**
88 * @brief Returns the he number of total spherical functions indexed by {ilm} associated with an atomic number. If the atomic number does not exist, it returns a zero.
89 * @return d_numSphericalFunctions.find(atomicNumber)->size()
90 */
93
94 /**
95 * @brief Returns the he number of radial spherical functions indexed by {i} associated with an atomic number. If the atomic number does not exist, it returns a zero.
96 * @return d_numRadialSphericalFunctions.find(atomicNumber)->size()
97 */
100 /**
101 * @brief Returns the total number of total spherical functions indexed by {ilm} present in the current processor. If the atomic number does not exist, it returns a zero.
102 */
105 /**
106 * @brief Returns the maximum number of total spherical functions indexed by {ilm} across all atom Types present in atomNumbers vector
107 */
110 /**
111 * @brief
112 * @param[out] totalAtomsInCurrentProcessor number of atoms in current
113 * processor based on compact support
114 * @param[out] totalNonLocalElements number of nonLocal elements in current
115 * processor
116 * @param[out] numberCellsForEachAtom number of cells associated which each
117 * atom in the current processor. vecot of size totalAtomsInCurrentProcessor
118 * @param[out] numberCellsAccumNonLocalAtoms number of cells accumulated
119 * till iatom in current processor. vector of size
120 * totalAtomsInCurrentProcessor
121 */
122 void
124 dftfe::uInt &totalAtomsInCurrentProcessor,
125 dftfe::uInt &totalNonLocalElements,
126 std::vector<dftfe::uInt> &numberCellsForEachAtom,
127 std::vector<dftfe::uInt> &numberCellsAccumNonLocalAtoms,
128 std::vector<dftfe::uInt> &iElemNonLocalToElemIndexMap);
129
130 /**
131 * @brief Returns the total number of total radial-spherical functions indexed by {i} present in atomicNumbers list.
132 */
135
136 /**
137 * @brief Returns the shared_ptr of AtomCenteredSphericalFunctionBase associated with std::pair(atomic Number and lQuantumNo)
138 */
139 const std::map<std::pair<dftfe::uInt, dftfe::uInt>,
140 std::shared_ptr<AtomCenteredSphericalFunctionBase>> &
142 /**
143 * @brief Returns the vector of size Natoms of all atoms in system
144 */
145 const std::vector<dftfe::uInt> &
147 /**
148 * @brief Returns the atomIds of atoms present in current processor
149 */
150 const std::vector<dftfe::uInt> &
152 /**
153 * @brief Returns the startIndex of spherical Function alpha associated with atomic number Znum
154 */
155 const dftfe::uInt
157 // COmputes the sparsity Pattern for the compact support Fn
158 // cutOffVal the max/min value to consider to be part of copact support
159 // cutOffType = 0 based on Fn Value, cutOffType = 1 based on Distance from
160 // atom
161 template <typename NumberType>
162 void
164 std::shared_ptr<
166 double,
168 &basisOperationsPtr,
169 const dftfe::uInt quadratureIndex,
170 const double cutOffVal = 1.0E-8,
171 const dftfe::uInt cutOffType = 0);
172
173
174 std::vector<std::vector<dftfe::uInt>> d_elementIndexesInAtomCompactSupport;
175 void
176 setImageCoordinates(const std::vector<dftfe::Int> &imageIds,
177 const std::vector<std::vector<double>> &periodicCoords);
178
179
180
181 const std::vector<dftfe::Int> &
183
184 const std::map<dftfe::uInt, std::vector<dftfe::Int>> &
186
187 bool
189
190 void
192 const std::map<dftfe::uInt, std::vector<dftfe::Int>> &sparsityPattern,
193 const std::vector<std::vector<dealii::CellId>>
194 &elementIdsInAtomCompactSupport,
195 const std::vector<std::vector<dftfe::uInt>>
196 &elementIndexesInAtomCompactSupport,
197 const std::vector<dftfe::uInt> &atomIdsInCurrentProcess,
198 dftfe::uInt numberElements);
199
200 const dftfe::uInt
202
203 private:
204 // A flattened vector that stores the coordinates of the atoms of interest
205 // in the unit cell
206 // Coord of atom I is stored at 3*I +0 ( x-coord),3*I+1 ( y-coord),3*I+2 (
207 // z-coord)
208 std::vector<double> d_atomCoords;
209
210 // A vector of size = number of atoms of interest
211 // the Ith atom in d_atomicNumbers has its coordinates
212 // in d_atomCoords[3*I+0], d_atomCoords[3*I+1], d_atomCoords[3*I+2]
213 std::vector<dftfe::uInt> d_atomicNumbers;
214
215 // This maps the atom I in the unit cell to all its image atoms.
216 // number of image atoms of Ith atom = d_periodicImageCoord[I].size()/ dim
217 // with dim=3 The coordinates are stored as a flattened vector
218 std::map<dftfe::uInt, std::vector<double>> d_periodicImageCoord;
219
220
221 // This maps, from std::pair<atomic number, \alpha> to S_{z,\alpha},
222 // where \alpha is the index for unique radial function
223 std::map<std::pair<dftfe::uInt, dftfe::uInt>,
224 std::shared_ptr<AtomCenteredSphericalFunctionBase>>
226 // Stores the number of distinct Radial Functions for a particular AtomType
227 std::map<dftfe::uInt, dftfe::uInt> d_numRadialSphericalFunctions;
228 // Stores the number of distinct Functions include m for a particular
229 // AtomType
230 std::map<dftfe::uInt, dftfe::uInt> d_numSphericalFunctions;
231 // This maps is between atomId in unit cell and the sparsity pattern of the
232 // atom and its images in the unitcell domain.
233 std::map<dftfe::uInt, std::vector<dftfe::Int>> d_sparsityPattern;
234 //
235 std::vector<std::vector<dealii::CellId>> d_elementIdsInAtomCompactSupport;
236 // std::vector<std::vector<dftfe::uInt>>
237 // d_elementIndexesInAtomCompactSupport;
238 std::vector<std::vector<dealii::DoFHandler<3>::active_cell_iterator>>
240 std::vector<dftfe::uInt> d_AtomIdsInCurrentProcess;
241 std::vector<dftfe::uInt> d_offsetLocation;
242 std::vector<std::vector<dftfe::Int>> d_AtomIdsInElement;
243 std::map<dftfe::uInt, std::vector<dftfe::uInt>>
245
246 }; // end of class AtomCenteredSphericalFunctionContainerBase
247} // end of namespace dftfe
248
249#endif // DFTFE_ATOMCENTEREDSPHERICALFUNCTIONCONTAINERBASE_H
Definition AtomCenteredSphericalFunctionContainer.h:37
void init(const std::vector< dftfe::uInt > &atomicNumbers, const std::map< std::pair< dftfe::uInt, dftfe::uInt >, std::shared_ptr< AtomCenteredSphericalFunctionBase > > &listOfSphericalFunctions)
Initialises the class with the atomicNumbers of various atoms and the AtomCenteredSphericalFn of vari...
void computeSparseStructure(std::shared_ptr< dftfe::basis::FEBasisOperations< NumberType, double, dftfe::utils::MemorySpace::HOST > > &basisOperationsPtr, const dftfe::uInt quadratureIndex, const double cutOffVal=1.0E-8, const dftfe::uInt cutOffType=0)
std::map< dftfe::uInt, dftfe::uInt > d_numRadialSphericalFunctions
Definition AtomCenteredSphericalFunctionContainer.h:227
std::vector< std::vector< dealii::CellId > > d_elementIdsInAtomCompactSupport
Definition AtomCenteredSphericalFunctionContainer.h:235
const std::vector< dftfe::Int > & getAtomIdsInElement(dftfe::uInt iElem)
const std::vector< dftfe::uInt > & getAtomIdsInCurrentProcess() const
Returns the atomIds of atoms present in current processor.
void setImageCoordinates(const std::vector< dftfe::Int > &imageIds, const std::vector< std::vector< double > > &periodicCoords)
void initaliseCoordinates(const std::vector< double > &atomCoords, const std::vector< std::vector< double > > &periodicCoords, const std::vector< dftfe::Int > &imageIds)
Initialises the position of atoms, the image posisiton and image ids after every update of atom posit...
std::vector< dftfe::uInt > d_atomicNumbers
Definition AtomCenteredSphericalFunctionContainer.h:213
std::vector< dftfe::uInt > d_offsetLocation
Definition AtomCenteredSphericalFunctionContainer.h:241
const std::map< dftfe::uInt, std::vector< dftfe::Int > > & getSparsityPattern()
void getTotalAtomsAndNonLocalElementsInCurrentProcessor(dftfe::uInt &totalAtomsInCurrentProcessor, dftfe::uInt &totalNonLocalElements, std::vector< dftfe::uInt > &numberCellsForEachAtom, std::vector< dftfe::uInt > &numberCellsAccumNonLocalAtoms, std::vector< dftfe::uInt > &iElemNonLocalToElemIndexMap)
const std::vector< dftfe::uInt > & getAtomicNumbers() const
Returns the vector of size Natoms of all atoms in system.
dftfe::uInt getTotalNumberOfSphericalFunctionsPerAtom(dftfe::uInt atomicNumber)
Returns the he number of total spherical functions indexed by {ilm} associated with an atomic number....
std::map< dftfe::uInt, std::vector< dftfe::Int > > d_sparsityPattern
Definition AtomCenteredSphericalFunctionContainer.h:233
std::vector< dftfe::uInt > d_AtomIdsInCurrentProcess
Definition AtomCenteredSphericalFunctionContainer.h:240
dftfe::uInt getTotalNumberOfSphericalFunctionsInCurrentProcessor()
Returns the total number of total spherical functions indexed by {ilm} present in the current process...
void getDataForSparseStructure(const std::map< dftfe::uInt, std::vector< dftfe::Int > > &sparsityPattern, const std::vector< std::vector< dealii::CellId > > &elementIdsInAtomCompactSupport, const std::vector< std::vector< dftfe::uInt > > &elementIndexesInAtomCompactSupport, const std::vector< dftfe::uInt > &atomIdsInCurrentProcess, dftfe::uInt numberElements)
const dftfe::uInt getOffsetLocation(const dftfe::uInt iAtom)
std::map< dftfe::uInt, std::vector< dftfe::uInt > > d_totalSphericalFunctionIndexStart
Definition AtomCenteredSphericalFunctionContainer.h:244
const std::map< std::pair< dftfe::uInt, dftfe::uInt >, std::shared_ptr< AtomCenteredSphericalFunctionBase > > & getSphericalFunctions() const
Returns the shared_ptr of AtomCenteredSphericalFunctionBase associated with std::pair(atomic Number a...
const dftfe::uInt getTotalSphericalFunctionIndexStart(dftfe::uInt Znum, dftfe::uInt alpha)
Returns the startIndex of spherical Function alpha associated with atomic number Znum.
std::vector< double > d_atomCoords
Definition AtomCenteredSphericalFunctionContainer.h:208
const std::vector< double > & getAtomCoordinates() const
Returns the cooridnates of atom present in domain.
dftfe::uInt getTotalNumberOfRadialSphericalFunctions()
Returns the total number of total radial-spherical functions indexed by {i} present in atomicNumbers ...
std::map< std::pair< dftfe::uInt, dftfe::uInt >, std::shared_ptr< AtomCenteredSphericalFunctionBase > > d_sphericalFunctionsContainer
Definition AtomCenteredSphericalFunctionContainer.h:225
dftfe::uInt getTotalNumberOfRadialSphericalFunctionsPerAtom(dftfe::uInt atomicNumber)
Returns the he number of radial spherical functions indexed by {i} associated with an atomic number....
std::map< dftfe::uInt, std::vector< double > > d_periodicImageCoord
Definition AtomCenteredSphericalFunctionContainer.h:218
std::map< dftfe::uInt, dftfe::uInt > d_numSphericalFunctions
Definition AtomCenteredSphericalFunctionContainer.h:230
std::vector< std::vector< dftfe::uInt > > d_elementIndexesInAtomCompactSupport
Definition AtomCenteredSphericalFunctionContainer.h:174
dftfe::uInt getMaximumNumberOfSphericalFunctions()
Returns the maximum number of total spherical functions indexed by {ilm} across all atom Types presen...
dftfe::uInt getNumAtomCentersSize()
Returns the number of atoms present in domain.
std::vector< std::vector< dftfe::Int > > d_AtomIdsInElement
Definition AtomCenteredSphericalFunctionContainer.h:242
const std::map< dftfe::uInt, std::vector< double > > & getPeriodicImageCoordinatesList() const
Returns the map of atomId vs vector of image coordinates.
std::vector< std::vector< dealii::DoFHandler< 3 >::active_cell_iterator > > d_elementOneFieldIteratorsInAtomCompactSupport
Definition AtomCenteredSphericalFunctionContainer.h:239
Definition FEBasisOperations.h:85
@ HOST
Definition MemorySpaceType.h:34
Definition pseudoPotentialToDftfeConverter.cc:34
std::uint32_t uInt
Definition TypeConfig.h:10