DFT-FE 1.3.0-pre
Density Functional Theory With Finite-Elements
Loading...
Searching...
No Matches
dft.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
18#ifndef dft_H_
19#define dft_H_
20#include <constants.h>
22#include <elpaScalaManager.h>
23#include <headers.h>
24#include <MemorySpaceType.h>
25#include <MemoryStorage.h>
26#include <FEBasisOperations.h>
27#include <BLASWrapper.h>
28#include <AuxDensityMatrix.h>
30
31#include <complex>
32#include <deque>
33#include <iomanip>
34#include <iostream>
35#include <numeric>
36#include <sstream>
37
38#ifdef DFTFE_WITH_DEVICE
40# include "deviceKernelsGeneric.h"
41# include <linearSolverCGDevice.h>
43#endif
44
45// ************* For debugging purposes only. Remove afterwards
46#include "hubbardClass.h"
47
49#include <dealiiLinearSolver.h>
50#include <dftParameters.h>
51#include <eigenSolver.h>
52#include <interpolation.h>
58#include <vselfBinsManager.h>
59#include <excManager.h>
60#include <dftd.h>
61#include "dftBase.h"
62#ifdef USE_PETSC
63# include <petsc.h>
64
65# include <slepceps.h>
66#endif
67
68#include <mixingClass.h>
69#include <oncvClass.h>
70#include <AuxDensityMatrix.h>
75#include <groupSymmetry.h>
76
77namespace dftfe
78{
79 //
80 // Initialize Namespace
81 //
82
83
84
85#ifndef DOXYGEN_SHOULD_SKIP_THIS
86
87 struct orbital
88 {
89 dftfe::uInt atomID;
90 dftfe::uInt waveID;
91 dftfe::uInt Z, n, l;
92 dftfe::Int m;
93 alglib::spline1dinterpolant psi;
94 };
95
96#endif /* DOXYGEN_SHOULD_SKIP_THIS */
97
98 /**
99 * @brief This class is the primary interface location of all other parts of the DFT-FE code
100 * for all steps involved in obtaining the Kohn-Sham DFT ground-state
101 * solution.
102 *
103 * @author Shiva Rudraraju, Phani Motamarri, Sambit Das
104 */
105 template <dftfe::utils::MemorySpace memorySpace>
106 class dftClass : public dftBase
107 {
108 public:
109 /**
110 * @brief dftClass constructor
111 *
112 * @param[in] mpi_comm_parent parent communicator
113 * @param[in] mpi_comm_domain mpi_communicator for domain decomposition
114 * parallelization
115 * @param[in] interpoolcomm mpi_communicator for parallelization over k
116 * points
117 * @param[in] interBandGroupComm mpi_communicator for parallelization over
118 * bands
119 * @param[in] scratchFolderName scratch folder name
120 * @param[in] dftParams dftParameters object containg parameter values
121 * parsed from an input parameter file in dftfeWrapper class
122 */
123 dftClass(const MPI_Comm &mpiCommParent,
124 const MPI_Comm &mpi_comm_domain,
125 const MPI_Comm &interpoolcomm,
126 const MPI_Comm &interBandGroupComm,
127 const std::string &scratchFolderName,
128 dftParameters &dftParams);
129
130 /**
131 * @brief dftClass destructor
132 */
134
135 /**
136 * @brief atomic system pre-processing steps.
137 *
138 * Reads the coordinates of the atoms.
139 * If periodic calculation, reads fractional coordinates of atoms in the
140 * unit-cell, lattice vectors, kPoint quadrature rules to be used and also
141 * generates image atoms. Also determines orbital-ordering
142 */
143 void
145
146 /**
147 * @brief Does KSDFT problem pre-processing steps including mesh generation calls.
148 */
149 void
151
152 /**
153 * @brief Does KSDFT problem pre-processing steps but without remeshing.
154 */
155 void
156 initNoRemesh(const bool updateImagesAndKPointsAndVselfBins = true,
157 const bool checkSmearedChargeWidthsForOverlap = true,
158 const bool useSingleAtomSolutionOverride = false,
159 const bool isMeshDeformed = false);
160
161
162
163 /**
164 * @brief FIXME: legacy call, move to main.cc
165 */
166 void
168
169
170 void
172
173 /**
174 * @brief Writes inital density and mesh to file.
175 */
176 void
178
179 /**
180 * @brief compute approximation to ground-state without solving the SCF iteration
181 */
182 void
184 /**
185 * @brief compute bands without solving the SCF iteration
186 */
187 void
189 /**
190 * @brief Kohn-Sham ground-state solve using SCF iteration
191 *
192 * @return tuple of boolean flag on whether scf converged,
193 * and L2 norm of residual electron-density of the last SCF iteration step
194 *
195 */
196 std::tuple<bool, double>
197 solve(const bool computeForces = true,
198 const bool computestress = true,
199 const bool restartGroundStateCalcFromChk = false);
200
201 void
203
204
205 void
213
214 /**
215 * @brief Copies the residual residualValues=outValues-inValues
216 */
217 double
220 &outValues,
222 &inValues,
224 &residualValues,
226 &JxW,
227 const bool computeNorm);
228
229
230 double
232 const distributedCPUVec<double> &inValues,
233 distributedCPUVec<double> &residualValues);
234
235
236 /**
237 * @brief Computes the diagonal mass matrix for rho nodal grid, used for nodal mixing
238 */
239 void
242 &massVec);
243
244 void
245 initializeKohnShamDFTOperator(const bool initializeCublas = true);
246
247
248 void
250
251
252 void
254
255
256 double
258
259 double
261
262 double
264
267
270
271 double
273
274 void
276
277 virtual void
279
280 /**
281 * @brief Number of Kohn-Sham eigen values to be computed
282 */
284
285
286 /**
287 * @brief Number of random wavefunctions
288 */
290
291 void
293
294 /**
295 *@brief Get local dofs global indices real
296 */
297 const std::vector<dealii::types::global_dof_index> &
299
300 /**
301 *@brief Get local dofs global indices imag
302 */
303 const std::vector<dealii::types::global_dof_index> &
305
306 /**
307 *@brief Get local dofs local proc indices real
308 */
309 const std::vector<dealii::types::global_dof_index> &
311
312 /**
313 *@brief Get local dofs local proc indices imag
314 */
315 const std::vector<dealii::types::global_dof_index> &
317
318 /**
319 *@brief Get matrix free data object
320 */
321 const dealii::MatrixFree<3, double> &
323
324
325 /** @brief Updates atom positions, remeshes/moves mesh and calls appropriate reinits.
326 *
327 * Function to update the atom positions and mesh based on the provided
328 * displacement input. Depending on the maximum displacement magnitude this
329 * function decides wether to do auto remeshing or move mesh using Gaussian
330 * functions. Additionaly this function also wraps the atom position across
331 * the periodic boundary if the atom moves across it beyond a certain
332 * magnitude. In case of floating atoms, only the atomic positions are
333 * updated keeping the mesh fixed. This function also calls initNoRemesh to
334 * reinitialize all the required FEM and KSDFT objects.
335 *
336 * @param[in] globalAtomsDisplacements vector containing the displacements
337 * (from current position) of all atoms (global).
338 * @return void.
339 */
340 void
342 const std::vector<dealii::Tensor<1, 3, double>> &globalAtomsDisplacements,
343 const double maxJacobianRatioFactor = 1.25,
344 const bool useSingleAtomSolutionsOverride = false);
345
346
347 /**
348 * @brief writes the current domain bounding vectors and atom coordinates to files, which are required for
349 * geometry relaxation restart
350
351 */
352 void
354
355 /**
356 * @brief writes the current domain bounding vectors and atom coordinates to files for
357 * structural optimization and dynamics restarts. The coordinates are stored
358 * as: 1. fractional for semi-periodic/periodic 2. Cartesian for
359 * non-periodic.
360 * @param[in] Path The folder path to store the atom coordinates required
361 * during restart.
362 */
363 void
364 writeDomainAndAtomCoordinates(const std::string Path) const;
365
366 /**
367 * @brief writes atomistics data for subsequent post-processing. Related to
368 * WRITE STRUCTURE ENERGY FORCES DATA POST PROCESS input parameter.
369 * @param[in] Path The folder path to store the atomistics data.
370 */
371 void
372 writeStructureEnergyForcesDataPostProcess(const std::string Path) const;
373
374 /**
375 * @brief writes quadrature grid information and associated spin-up
376 * and spin-down electron-density for post-processing
377 * @param[in] Path The folder path to store the atomistics data.
378 */
379 virtual void
380 writeGSElectronDensity(const std::string Path) const;
381
382
383 /**
384 * @brief Gets the current atom Locations in cartesian form
385 * (origin at center of domain) from dftClass
386 */
387 const std::vector<std::vector<double>> &
389
390 /**
391 * @brief Gets the nearest atom distance for each atom
392 */
393 const std::vector<double> &
395
396 /**
397 * @brief Gets the current image atom Locations in cartesian form
398 * (origin at center of domain) from dftClass
399 */
400 const std::vector<std::vector<double>> &
402
403 /**
404 * @brief Gets the current image atom ids from dftClass
405 */
406 const std::vector<dftfe::Int> &
408
409 /**
410 * @brief Gets the current atom Locations in fractional form
411 * from dftClass (only applicable for periodic and semi-periodic BCs)
412 */
413 const std::vector<std::vector<double>> &
415
416
417
418 /**
419 * @brief Gets the current cell lattice vectors
420 *
421 * @return std::vector<std::vector<double>> 3 \times 3 matrix,lattice[i][j]
422 * corresponds to jth component of ith lattice vector
423 */
424 const std::vector<std::vector<double>> &
425 getCell() const;
426
427 /**
428 * @brief Gets the current cell volume
429 *
430 */
431 double
433
434 /**
435 * @brief Gets the current atom types from dftClass
436 */
437 const std::set<dftfe::uInt> &
439
440 /**
441 * @brief Gets the current atomic forces from dftClass
442 */
443 const std::vector<double> &
445
446 /**
447 * @brief Gets the current cell stress from dftClass
448 */
449 const dealii::Tensor<2, 3, double> &
451
452 /**
453 * @brief Get reference to dftParameters object
454 */
457
458
459 /**
460 * @brief Get reference the memorySpace templated eigen vectors
461 */
464
465 /**
466 * @brief Get reference the host eigen vectors
467 */
471
472 /**
473 * @brief Get reference to the eigen values
474 */
475 const std::vector<std::vector<double>> &
477
478 /**
479 * @brief Get the value of fermi energy
480 */
481 double
483
484 /**
485 * @brief Get the number of electrons
486 */
487 double
489
490 void
491 setNumElectrons(dftfe::uInt inputNumElectrons);
492
495#ifdef DFTFE_WITH_DEVICE
496
497 /**
498 * @brief Get the Ptr to Chebyshev solver in device
499 */
500 chebyshevOrthogonalizedSubspaceIterationSolverDevice *
501 getSubspaceIterationSolverDevice();
502#endif
503 /**
504 * @brief Get the Ptr to Chebyshev solver in host
505 */
508
509 /**
510 * @brief Function that computes the Eigen space for the Kohn Sham operator.
511 */
512 void
514 const dftfe::uInt s,
515 const dftfe::uInt kPointIndex,
517 &kohnShamDFTEigenOperator,
518 elpaScalaManager &elpaScala,
519 chebyshevOrthogonalizedSubspaceIterationSolver &subspaceIterationSolver,
520 std::vector<double> &residualNormWaveFunctions,
521 const bool computeResidual,
522 const bool useMixedPrec = false,
523 const bool isFirstScf = false);
524
525
526#ifdef DFTFE_WITH_DEVICE
527 /**
528 * @brief Function that computes the Eigen space for the Kohn Sham operator in device.
529 */
530 void
532 const dftfe::uInt s,
533 const dftfe::uInt kPointIndex,
535 &kohnShamDFTEigenOperator,
536 elpaScalaManager &elpaScala,
537 chebyshevOrthogonalizedSubspaceIterationSolverDevice
538 &subspaceIterationSolverDevice,
539 std::vector<double> &residualNormWaveFunctions,
540 const bool computeResidual,
541 const dftfe::uInt numberRayleighRitzAvoidancePasses = 0,
542 const bool useMixedPrec = false,
543 const bool isFirstScf = false);
544#endif
545
546 /**
547 *@brief Computes Fermi-energy obtained by imposing constraint on the number of electrons
548 */
549 void
551 const std::vector<std::vector<double>> &eigenValuesInput,
552 const double numElectronsInput);
553
554 /**
555 *@brief find HOMO eigenvalue for pure state
556 */
557 void
559 const std::vector<std::vector<double>> &eigenValuesInput,
560 const double numElectronsInput);
561
562 /**
563 *@brief Computes the kinetic energy
564 */
565 double
568 &kineticEnergyDensityValues);
569
570 /**
571 *@brief get the Ptr to the operator class ( Kohn Sham Base Operator)
572 */
575
576 /**
577 *@brief get the index of the DoF Handler corresponding to
578 *
579 */
582
585
586 const std::vector<double> &
588
591
594
595 const dealii::MatrixFree<3, double> &
597
598 dealii::AffineConstraints<double> *
600
603
606
609
610
611 std::shared_ptr<
613 double,
616
617 std::shared_ptr<
620
621 std::shared_ptr<
622 dftfe::basis::
623 FEBasisOperations<double, double, dftfe::utils::MemorySpace::HOST>>
625
626 std::shared_ptr<
629
630 std::shared_ptr<dftfe::linearAlgebra::BLASWrapper<memorySpace>>
632
633
634 std::shared_ptr<
637
638 std::vector<
641
642 std::vector<
645
646 /**
647 *@brief l2 projection
648 */
649 void
651 const std::shared_ptr<
653 FEBasisOperations<double, double, dftfe::utils::MemorySpace::HOST>>
654 &basisOperationsPtr,
655 const dealii::AffineConstraints<double> &constraintMatrix,
656 const dftfe::uInt dofHandlerId,
657 const dftfe::uInt quadratureId,
659 &quadratureValueData,
660 distributedCPUVec<double> &nodalField);
661
662 /// map of atom node number and atomic weight
663 std::map<dealii::types::global_dof_index, double> &
665
666 /// non-intersecting smeared charges of all atoms at quad points
667 std::map<dealii::CellId, std::vector<double>> &
669
672
673 const dealii::AffineConstraints<double> *
675
676 const std::vector<std::vector<double>> &
678
679 const MPI_Comm &
680 getMPIDomain() const override;
681
682 const MPI_Comm &
683 getMPIParent() const override;
684
685 const MPI_Comm &
686 getMPIInterPool() const override;
687
688 const MPI_Comm &
689 getMPIInterBand() const override;
690
691 const std::map<dealii::CellId, std::vector<dftfe::uInt>> &
693
694 void
696
697 void
699 const std::shared_ptr<
701 FEBasisOperations<double, double, dftfe::utils::MemorySpace::HOST>>
702 &basisOperationsPtr,
703 const dftfe::uInt densityQuadratureId,
705 &rhoQuadValues,
706 const std::map<dealii::CellId, std::vector<double>> *bQuadValues);
707
710
711 /**
712 *@brief Returns the shared ptr to hubbard class
713 */
714 std::shared_ptr<hubbard<dataTypes::number, memorySpace>>
716
717 /**
718 *@brief Function to check if hubbard corrections is being used
719 */
720 bool
722
723 /**
724 *@brief write wavefunction solution fields
725 */
726 void
727 outputWfc(const std::string outputFileName = "wfcOutput");
728
729 /**
730 *@brief return the pseudo potential field
731 */
732 const std::map<dealii::CellId, std::vector<double>> &
734
735 private:
736 /**
737 * @brief generate image charges and update k point cartesian coordinates based
738 * on current lattice vectors
739 */
740 void
742
743 /**
744 * @brief Checks if the Exc functional requires Hubbard correction
745 * and sets up the Hubbard class if required.
746 */
747 void
749
750 /**
751 *@brief project ground state electron density from previous mesh into
752 * the new mesh to be used as initial guess for the new ground state solve
753 */
754 void
756
757 /**
758 *@brief save triangulation information and rho quadrature data to checkpoint file for restarts
759 */
760 void
762
763 /**
764 *@brief load triangulation information rho quadrature data from checkpoint file for restarted run
765 */
766 void
768 /**
769 * @brief save data of quad points to checkpoint file. Used for restart calculations, nscf and bands.
770 *
771 * @param[in] basisOperationsPtr basisoperationsPtr object
772 * @param[in] quadratureId quadrature Id of quad point used in checkpoint
773 * file
774 * @param[out] quadratureValueData quadrature data of field that is to be
775 * saved
776 * @param[in] fieldDimension dimension of field.
777 * @param[in] fieldName file name of checkpoint data to be saved
778 * @param[in] folderPath restart folder name
779 * @param[in] mpi_comm_parent parent communicator
780 * @param[in] mpi_comm_domain mpi_communicator for domain decomposition
781 * parallelization
782 * @param[in] interpoolcomm mpi_communicator for parallelization over k
783 * points
784 * @param[in] interBandGroupComm mpi_communicator for parallelization over
785 * bands
786 */
787 void
789 const std::shared_ptr<
791 double,
793 &basisOperationsPtr,
794 const dftfe::uInt quadratureId,
796 &quadratureValueData,
797 const dftfe::uInt fieldDimension,
798 const std::string &fieldName,
799 const std::string &folderPath,
800 const MPI_Comm &mpi_comm_parent,
801 const MPI_Comm &mpi_comm_domain,
802 const MPI_Comm &interpoolcomm,
803 const MPI_Comm &interBandGroupComm);
804 /**
805 * @brief loads data from quad points of checkpoint file. Used for restart calculations, nscf and bands.
806 *
807 * @param[in] basisOperationsPtr basisoperationsPtr object
808 * @param[in] quadratureId quadrature Id of quad point used in checkpoint
809 * file
810 * @param[out] quadratureValueData quadrature data of field that is to be
811 * loaded
812 * @param[in] fieldDimension dimension of field.
813 * @param[in] fieldName file name containing checkpoint data
814 * @param[in] folderPath restart folder name
815 * @param[in] mpi_comm_parent parent communicator
816 * @param[in] mpi_comm_domain mpi_communicator for domain decomposition
817 * parallelization
818 * @param[in] interpoolcomm mpi_communicator for parallelization over k
819 * points
820 * @param[in] interBandGroupComm mpi_communicator for parallelization over
821 * bands
822 */
823 void
825 const std::shared_ptr<
827 double,
829 &basisOperationsPtr,
830 const dftfe::uInt quadratureId,
832 &quadratureValueData,
833 const dftfe::uInt fieldDimension,
834 const std::string &fieldName,
835 const std::string &folderPath,
836 const MPI_Comm &mpi_comm_parent,
837 const MPI_Comm &mpi_comm_domain,
838 const MPI_Comm &interpoolcomm,
839 const MPI_Comm &interBandGroupComm);
840
841 void
843 void
844 writeMesh(std::string meshFileName);
845
846 /// creates datastructures related to periodic image charges
847 void
848 generateImageCharges(const double pspCutOff,
849 std::vector<dftfe::Int> &imageIds,
850 std::vector<double> &imageCharges,
851 std::vector<std::vector<double>> &imagePositions);
852
853 void
855 const double pspCutOff,
856 const std::vector<dftfe::Int> &imageIds,
857 const std::vector<std::vector<double>> &imagePositions,
858 std::vector<std::vector<dftfe::Int>> &globalChargeIdToImageIdMap);
859
860 void
862
863 //
864 // generate mesh using a-posteriori error estimates
865 //
867 computeTraceXtHX(dftfe::uInt numberWaveFunctionsEstimate);
868 double
869 computeTraceXtKX(dftfe::uInt numberWaveFunctionsEstimate);
870
871
872 /**
873 *@brief moves the triangulation vertices using Gaussians such that the all atoms are on triangulation vertices
874 */
875 void
876 moveMeshToAtoms(dealii::Triangulation<3, 3> &triangulationMove,
877 dealii::Triangulation<3, 3> &triangulationSerial,
878 bool reuseFlag = false,
879 bool moveSubdivided = false);
880
881 /**
882 *@brief a
883 */
884 void
886
887 /**
888 *@brief a
889 */
890 void
892
893 /**
894 * Initializes the guess of electron-density and single-atom wavefunctions
895 * on the mesh, maps finite-element nodes to given atomic positions,
896 * initializes pseudopotential files and exchange-correlation functionals to
897 * be used based on user-choice. In periodic problems, periodic faces are
898 * mapped here. Further finite-element nodes to be pinned for solving the
899 * Poisson problem electro-static potential is set here
900 */
901 void
903 dealii::parallel::distributed::Triangulation<3> &triangulation);
904 void
905 initBoundaryConditions(const bool recomputeBasisData = true,
906 const bool meshOnlyDeformed = false,
907 const bool vselfPerturbationUpdateForStress = false);
908 void
910 void
911 initPseudoPotentialAll(const bool updateNonlocalSparsity = true);
912
913 /**
914 * create a dofHandler containing finite-element interpolating polynomial
915 * twice of the original polynomial required for Kerker mixing and
916 * initialize various objects related to this refined dofHandler
917 */
918 void
920 dealii::parallel::distributed::Triangulation<3> &triangulation);
921 void
922 initpRefinedObjects(const bool recomputeBasisData,
923 const bool meshOnlyDeformed,
924 const bool vselfPerturbationUpdateForStress = false);
925
926
927 /**
928 *@brief Sets inhomegeneous dirichlet boundary conditions upto quadrupole for total potential constraints on
929 * non-periodic boundary (boundary id==0).
930 *
931 * @param[in] dofHandler
932 * @param[out] constraintMatrix dealii::AffineConstraints<double> object
933 *with inhomogeneous Dirichlet boundary condition entries added
934 */
935 void
937 const dealii::DoFHandler<3> &_dofHandler,
938 const dealii::AffineConstraints<double> &onlyHangingNodeConstraints,
939 dealii::AffineConstraints<double> &constraintMatrix);
940
941
942 /**
943 *@brief add atomic densities at quadrature points
944 *
945 */
946 void
949 &quadratureValueData,
951 &quadratureGradValueData,
952 const bool isConsiderGradData = false);
953
954
955 /**
956 *@brief Finds the global dof ids of the nodes containing atoms.
957 *
958 * @param[in] dofHandler
959 * @param[out] atomNodeIdToChargeValueMap local map of global dof id to atom
960 *charge id
961 */
962 void
963 locateAtomCoreNodes(const dealii::DoFHandler<3> &_dofHandler,
964 std::map<dealii::types::global_dof_index, double>
965 &atomNodeIdToChargeValueMap);
966
967 /**
968 *@brief Sets homogeneous dirichlet boundary conditions on a node farthest from
969 * all atoms (pinned node). This is only done in case of periodic boundary
970 *conditions to get an unique solution to the total electrostatic potential
971 *problem.
972 *
973 * @param[in] dofHandler
974 * @param[in] constraintMatrixBase base dealii::AffineConstraints<double>
975 *object
976 * @param[out] constraintMatrix dealii::AffineConstraints<double> object
977 *with homogeneous Dirichlet boundary condition entries added
978 */
979 void
981 const dealii::DoFHandler<3> &_dofHandler,
982 const dealii::AffineConstraints<double> &constraintMatrixBase,
983 dealii::AffineConstraints<double> &constraintMatrix);
984
985 void
987
989
990 void
992
993 void
995 void
997 void
999 std::vector<std::vector<distributedCPUVec<double>>> eigenVectors);
1000 void
1002
1003
1004 /**
1005 *@brief computes density nodal data from wavefunctions
1006 */
1007 void
1009
1010
1011 void
1015 distributedCPUVec<double> &fvSpin1);
1016
1017 void
1019 void
1021 void
1023 void
1025 dftfe::uInt n,
1026 dftfe::uInt l,
1027 dftfe::uInt &flag);
1028 void
1030 const dealii::DoFHandler<3> &_dofHandler,
1031 const dftfe::uInt lpspQuadratureId,
1032 const dealii::MatrixFree<3, double> &_matrix_free_data,
1033 const dftfe::uInt _phiExtDofHandlerIndex,
1034 const dealii::AffineConstraints<double> &phiExtConstraintMatrix,
1035 const std::map<dealii::types::global_dof_index, dealii::Point<3>>
1036 &supportPoints,
1037 const vselfBinsManager &vselfBinManager,
1039 std::map<dealii::CellId, std::vector<double>> &_pseudoValues,
1040 std::map<dftfe::uInt, std::map<dealii::CellId, std::vector<double>>>
1041 &_pseudoValuesAtoms);
1042
1043
1044
1045 /**
1046 *@brief Sets homegeneous dirichlet boundary conditions for total potential constraints on
1047 * non-periodic boundary (boundary id==0).
1048 *
1049 * @param[in] dofHandler
1050 * @param[out] constraintMatrix dealii::AffineConstraints<double> object
1051 *with homogeneous Dirichlet boundary condition entries added
1052 */
1053 void
1055 const dealii::DoFHandler<3> &_dofHandler,
1056 const dealii::AffineConstraints<double> &onlyHangingNodeConstraints,
1057 dealii::AffineConstraints<double> &constraintMatrix);
1058
1059
1060
1061 /**
1062 *@brief Computes total charge by integrating the electron-density
1063 */
1064 double
1065 totalCharge(const dealii::DoFHandler<3> &dofHandlerOfField,
1066 const distributedCPUVec<double> &rhoNodalField);
1067
1068
1069 double
1071 const dealii::DoFHandler<3> &dofHandlerOfField,
1072 const std::map<dealii::CellId, std::vector<double>> *rhoQuadValues);
1073
1074 double
1076 const dealii::DoFHandler<3> &dofHandlerOfField,
1078 &rhoQuadValues);
1079
1080
1081 double
1082 totalCharge(const dealii::MatrixFree<3, double> &matrixFreeDataObject,
1083 const distributedCPUVec<double> &rhoNodalField);
1084
1085
1086
1087 double
1088 rhofieldl2Norm(const dealii::MatrixFree<3, double> &matrixFreeDataObject,
1089 const distributedCPUVec<double> &rhoNodalField,
1090 const dftfe::uInt dofHandlerId,
1091 const dftfe::uInt quadratureId);
1092
1093 double
1095 const dealii::MatrixFree<3, double> &matrixFreeDataObject,
1096 const distributedCPUVec<double> &rhoNodalField1,
1097 const distributedCPUVec<double> &rhoNodalField2,
1098 const dftfe::uInt dofHandlerId,
1099 const dftfe::uInt quadratureId);
1100
1101
1102 double
1103 fieldGradl2Norm(const dealii::MatrixFree<3, double> &matrixFreeDataObject,
1104 const distributedCPUVec<double> &field);
1105
1106
1107 /**
1108 *@brief l2 projection
1109 */
1110 void
1112 const std::shared_ptr<
1113 dftfe::basis::
1114 FEBasisOperations<double, double, dftfe::utils::MemorySpace::HOST>>
1115 &basisOperationsPtr,
1116 const dealii::AffineConstraints<double> &constraintMatrix,
1117 const dftfe::uInt dofHandlerId,
1118 const dftfe::uInt quadratureId,
1120 &quadratureValueData,
1121 distributedCPUVec<double> &nodalField);
1122
1123 /**
1124 *@brief Computes net magnetization from the difference of local spin densities
1125 */
1126 void
1129 &magQuadValues);
1130
1131 void
1133 const std::vector<
1135 &densityQuadValues);
1136
1137 void
1139 const std::vector<
1141 &densityQuadValues);
1142
1143 void
1145 const std::vector<
1147 &densityQuadValues);
1148
1149 /**
1150 *@brief normalize the input electron density
1151 */
1152 void
1154
1155 /**
1156 *@brief normalize input mag electron density to total magnetization
1157 *for use in constraint magnetization case (only for initial guess)
1158 */
1159 void
1161
1162
1163 /**
1164 *@brief normalize the output total electron density in each scf
1165 */
1166 void
1168
1169 /**
1170 *@brief normalize the electron density
1171 */
1172 void
1174
1175 /**
1176 *@brief Computes output electron-density from wavefunctions
1177 */
1178 void
1179 compute_rhoOut(const bool isGroundState = false);
1180
1181 /**
1182 *@brief Mixing schemes for mixing electron-density
1183 */
1184
1185 void
1187#ifdef DFTFE_WITH_DEVICE
1188 kerkerSolverProblemDeviceWrapperClass
1189 &kerkerPreconditionedResidualSolverProblemDevice,
1190 linearSolverCGDevice &CGSolverDevice,
1191#endif
1193 &kerkerPreconditionedResidualSolverProblem,
1194 dealiiLinearSolver &CGSolver,
1195 distributedCPUVec<double> &residualRho,
1196 distributedCPUVec<double> &preCondTotalDensityResidualVector);
1197
1198 double
1200
1201 double
1203 /**
1204 *@brief Computes Fermi-energy obtained by imposing separate constraints on the number of spin-up and spin-down electrons
1205 */
1206 void
1208 const std::vector<std::vector<double>> &eigenValuesInput);
1209
1210
1211 /**
1212 *@brief Find spin-up and spin-down channel HOMO eigenvalues
1213 */
1214 void
1216 const std::vector<std::vector<double>> &eigenValuesInput);
1217
1218
1219 /**
1220 *@brief compute density of states and local density of states
1221 */
1222 void
1223 compute_tdos(const std::vector<std::vector<double>> &eigenValuesInput,
1224 const std::string &fileName);
1225
1226 void
1227 compute_ldos(const std::vector<std::vector<double>> &eigenValuesInput,
1228 const std::string &fileName);
1229
1230 /**
1231 *@brief compute localization length
1232 */
1233 void
1234 compute_localizationLength(const std::string &locLengthFileName);
1235
1236 /**
1237 *@brief write electron density solution fields
1238 */
1239 void
1241
1242 /**
1243 *@brief write the KS eigen values for given BZ sampling/path
1244 */
1245 void
1247
1248 /**
1249 *@brief Computes the volume of the domain
1250 */
1251 double
1252 computeVolume(const dealii::DoFHandler<3> &_dofHandler);
1253
1254 /**
1255 *@brief Deforms the domain by the given deformation gradient and reinitializes the
1256 * dftClass datastructures.
1257 */
1258 void
1259 deformDomain(const dealii::Tensor<2, 3, double> &deformationGradient,
1260 const bool vselfPerturbationUpdateForStress = false,
1261 const bool useSingleAtomSolutionsOverride = false,
1262 const bool print = true);
1263
1264 /**
1265 *@brief Computes inner Product and Y = alpha*X + Y for complex vectors used during
1266 * periodic boundary conditions
1267 */
1268
1269#ifdef USE_COMPLEX
1270 std::complex<double>
1272
1273 void
1274 alphaTimesXPlusY(std::complex<double> alpha,
1277
1278#endif
1279 /**
1280 *@brief Sets dirichlet boundary conditions for total potential constraints on
1281 * non-periodic boundary (boundary id==0). Currently setting homogeneous bc
1282 *
1283 */
1284 void
1286
1287
1288 void
1290 const std::vector<
1292 &densityQuadValues,
1293 const std::vector<
1295 &gradDensityQuadValues,
1296 const std::vector<
1298 &tauQuadValues,
1299 const std::map<dealii::CellId, std::vector<double>> &rhoCore,
1300 const std::map<dealii::CellId, std::vector<double>> &gradRhoCore,
1302 &eigenVectorsFlattenedMemSpace,
1303 const std::vector<std::vector<double>> &eigenValues,
1304 const double fermiEnergy_,
1305 const double fermiEnergyUp_,
1306 const double fermiEnergyDown_,
1307 std::shared_ptr<AuxDensityMatrix<memorySpace>> auxDensityMatrixXCPtr);
1308
1309 std::shared_ptr<excManager<memorySpace>> d_excManagerPtr;
1311
1312 /**
1313 * stores required data for Kohn-Sham problem
1314 */
1317 std::set<dftfe::uInt> atomTypes;
1318
1319 /// FIXME: eventually it should be a map of atomic number to struct-
1320 /// {valence number, mesh input etc}
1321 std::map<dftfe::uInt, dftfe::uInt> d_atomTypeAtributes;
1322
1323 /// FIXME: remove atom type atributes from atomLocations
1324 std::vector<std::vector<double>> atomLocations, atomLocationsFractional,
1326 std::vector<std::vector<double>> d_atomLocationsAutoMesh;
1327 std::vector<std::vector<double>> d_imagePositionsAutoMesh;
1328
1329 /// Gaussian displacements of atoms read from file
1330 std::vector<dealii::Tensor<1, 3, double>> d_atomsDisplacementsGaussianRead;
1331
1332 ///
1334
1335 ///
1337
1339
1340 /// Gaussian generator parameter for force computation and Gaussian
1341 /// deformation of atoms and FEM mesh Gaussian generator: Gamma(r)=
1342 /// exp(-(r/d_gaussianConstant)^2) Stored for all domain atoms
1343 std::vector<double> d_gaussianConstantsForce;
1344
1345 /// Gaussian constants for automesh mesh movement stored for all domain
1346 /// atoms
1347 std::vector<double> d_gaussianConstantsAutoMesh;
1348
1349 /// composite generator flat top widths for all domain atoms
1350 std::vector<double> d_generatorFlatTopWidths;
1351
1352 /// flat top widths for all domain atoms in case of automesh mesh movement
1353 /// composite gaussian
1354 std::vector<double> d_flatTopWidthsAutoMeshMove;
1355
1356 /// smeared charge widths for all domain atoms
1357 std::vector<double> d_smearedChargeWidths;
1358
1359 /// smeared charge normalization scaling for all domain atoms
1360 std::vector<double> d_smearedChargeScaling;
1361
1362 /// nearest atom ids for all domain atoms
1363 std::vector<dftfe::uInt> d_nearestAtomIds;
1364
1365 /// nearest atom distances for all domain atoms
1366 std::vector<double> d_nearestAtomDistances;
1367
1368 ///
1370
1371 /// vector of lendth number of periodic image charges with corresponding
1372 /// master chargeIds
1373 std::vector<dftfe::Int> d_imageIds;
1374 // std::vector<dftfe::Int> d_imageIdsAutoMesh;
1375
1376
1377 /// vector of length number of periodic image charges with corresponding
1378 /// charge values
1379 std::vector<double> d_imageCharges;
1380
1381 /// vector of length number of periodic image charges with corresponding
1382 /// positions in cartesian coordinates
1383 std::vector<std::vector<double>> d_imagePositions;
1384
1385 /// globalChargeId to ImageChargeId Map
1386 std::vector<std::vector<dftfe::Int>> d_globalChargeIdToImageIdMap;
1387
1388 /// vector of lendth number of periodic image charges with corresponding
1389 /// master chargeIds , generated with a truncated pspCutoff
1390 std::vector<dftfe::Int> d_imageIdsTrunc;
1391
1392 /// vector of length number of periodic image charges with corresponding
1393 /// charge values , generated with a truncated pspCutoff
1394 std::vector<double> d_imageChargesTrunc;
1395
1396 /// vector of length number of periodic image charges with corresponding
1397 /// positions in cartesian coordinates, generated with a truncated pspCutOff
1398 std::vector<std::vector<double>> d_imagePositionsTrunc;
1399
1400 /// globalChargeId to ImageChargeId Map generated with a truncated pspCutOff
1401 std::vector<std::vector<dftfe::Int>> d_globalChargeIdToImageIdMapTrunc;
1402
1403 /// distance from the domain till which periodic images will be considered
1404 double d_pspCutOff = 15.0;
1405
1406 /// distance from the domain till which periodic images will be considered
1407 const double d_pspCutOffTrunc = 15.0;
1408
1409 /// cut-off distance from atom till which non-local projectors are
1410 /// non-trivial
1411 double d_nlPSPCutOff = 8.0;
1412
1413 /// non-intersecting smeared charges of all atoms at quad points
1414 std::map<dealii::CellId, std::vector<double>> d_bQuadValuesAllAtoms;
1415
1416 /// non-intersecting smeared charge gradients of all atoms at quad points
1417 std::map<dealii::CellId, std::vector<double>> d_gradbQuadValuesAllAtoms;
1418
1419 /// non-intersecting smeared charges atom ids of all atoms at quad points
1420 std::map<dealii::CellId, std::vector<dftfe::Int>> d_bQuadAtomIdsAllAtoms;
1421
1422 /// non-intersecting smeared charges atom ids of all atoms (with image atom
1423 /// ids separately accounted) at quad points
1424 std::map<dealii::CellId, std::vector<dftfe::Int>>
1426
1427 /// map of cell and non-trivial global atom ids (no images) for smeared
1428 /// charges for each bin
1429 std::map<dealii::CellId, std::vector<dftfe::uInt>> d_bCellNonTrivialAtomIds;
1430
1431 /// map of cell and non-trivial global atom ids (no images) for smeared
1432 /// charge for each bin
1433 std::vector<std::map<dealii::CellId, std::vector<dftfe::uInt>>>
1435
1436 /// map of cell and non-trivial global atom and image ids for smeared
1437 /// charges for each bin
1438 std::map<dealii::CellId, std::vector<dftfe::uInt>>
1440
1441 /// map of cell and non-trivial global atom and image ids for smeared charge
1442 /// for each bin
1443 std::vector<std::map<dealii::CellId, std::vector<dftfe::uInt>>>
1445
1446 /// minimum smeared charge width
1447 const double d_smearedChargeWidthMin = 0.4;
1448
1449 std::vector<orbital> waveFunctionsVector;
1450 std::map<
1452 std::map<dftfe::uInt, std::map<dftfe::uInt, alglib::spline1dinterpolant>>>
1454 std::map<dftfe::uInt, std::map<dftfe::uInt, std::map<dftfe::uInt, double>>>
1456
1457 /**
1458 * meshGenerator based object
1459 */
1461
1464
1465
1466 /// affine transformation object
1468
1469 /// meshMovementGaussianClass object
1471
1472 std::vector<dealii::Tensor<1, 3, double>>
1474 std::vector<dealii::Point<3>> d_controlPointLocationsCurrentMove;
1475
1476 /// volume of the domain
1478
1479 /// init wfc trunctation radius
1481
1482 /**
1483 * dealii based FE data structres
1484 */
1485 dealii::FESystem<3> FE, FEEigen;
1510 dealii::MatrixFree<3, double> matrix_free_data, d_matrixFreeDataPRefined;
1511 std::shared_ptr<
1513 double,
1516 std::shared_ptr<
1517 dftfe::basis::
1518 FEBasisOperations<double, double, dftfe::utils::MemorySpace::HOST>>
1520#if defined(DFTFE_WITH_DEVICE)
1521 std::shared_ptr<
1523 double,
1525 d_basisOperationsPtrDevice;
1526 std::shared_ptr<
1527 dftfe::basis::
1528 FEBasisOperations<double, double, dftfe::utils::MemorySpace::DEVICE>>
1529 d_basisOperationsPtrElectroDevice;
1530#endif
1531
1532 std::shared_ptr<
1535
1536 std::shared_ptr<dftfe::oncvClass<dataTypes::number, memorySpace>>
1538
1539 std::shared_ptr<
1542
1543 std::shared_ptr<
1544#if defined(DFTFE_WITH_DEVICE)
1546#else
1548#endif
1550
1551 std::map<dealii::types::global_dof_index, dealii::Point<3>> d_supportPoints,
1553 std::vector<const dealii::AffineConstraints<double> *> d_constraintsVector;
1554 std::vector<const dealii::AffineConstraints<double> *>
1556
1557 /**
1558 * parallel objects
1559 */
1560 const MPI_Comm mpi_communicator;
1561#if defined(DFTFE_WITH_DEVICE)
1562 utils::DeviceCCLWrapper *d_devicecclMpiCommDomainPtr;
1563#endif
1564 const MPI_Comm d_mpiCommParent;
1565 const MPI_Comm interpoolcomm;
1566 const MPI_Comm interBandGroupComm;
1573 std::vector<dealii::types::global_dof_index> local_dof_indicesReal,
1575 std::vector<dealii::types::global_dof_index> localProc_dof_indicesReal,
1577 std::vector<bool> selectedDofsHanging;
1578
1579 std::shared_ptr<dftfe::groupSymmetryClass> groupSymmetryPtr;
1580 std::shared_ptr<configurationalForceClass<memorySpace>> d_configForcePtr;
1581
1583
1585
1587#ifdef DFTFE_WITH_DEVICE
1588 poissonSolverProblemDeviceWrapperClass d_phiTotalSolverProblemDevice;
1589
1590 poissonSolverProblemDeviceWrapperClass d_phiPrimeSolverProblemDevice;
1591#endif
1592
1594
1596
1597 const std::string d_dftfeScratchFolderName;
1598
1599 /**
1600 * chebyshev subspace iteration solver objects
1601 *
1602 */
1604#ifdef DFTFE_WITH_DEVICE
1605 chebyshevOrthogonalizedSubspaceIterationSolverDevice
1606 d_subspaceIterationSolverDevice;
1607#endif
1608
1609 /**
1610 * constraint Matrices
1611 */
1612
1613 /**
1614 *object which is used to store dealii constraint matrix information
1615 *using STL vectors. The relevant dealii constraint matrix
1616 *has hanging node constraints and periodic constraints(for periodic
1617 *problems) used in eigen solve
1618 */
1621
1622
1623 /**
1624 *object which is used to store dealii constraint matrix information
1625 *using STL vectors. The relevant dealii constraint matrix
1626 *has hanging node constraints used in Poisson problem solution
1627 *
1628 */
1631
1632
1633#ifdef DFTFE_WITH_DEVICE
1635 d_constraintsNoneDataInfoDevice;
1636#endif
1637
1638
1639 dealii::AffineConstraints<double> constraintsNone, constraintsNoneEigen,
1641
1642 dealii::AffineConstraints<double> d_constraintsForTotalPotentialElectro;
1643
1644 dealii::AffineConstraints<double> d_constraintsForPhiPrimeElectro;
1645
1646 dealii::AffineConstraints<double> d_constraintsForHelmholtzRhoNodal;
1647
1648 dealii::AffineConstraints<double> d_constraintsPRefined;
1649
1650 dealii::AffineConstraints<double> d_constraintsPRefinedOnlyHanging;
1651
1652 dealii::AffineConstraints<double> d_constraintsRhoNodal;
1653
1654 dealii::AffineConstraints<double> d_constraintsRhoNodalOnlyHanging;
1655
1658
1659 /**
1660 * data storage for Kohn-Sham eigenvalues and partial occupancies
1661 */
1662 std::vector<std::vector<double>> eigenValues;
1663 std::vector<std::vector<double>> d_partialOccupancies;
1664
1665 /**
1666 * data storage for the occupancy of Kohn-Sham wavefunctions
1667 */
1668 std::vector<std::vector<double>> d_fracOccupancy;
1669
1670 std::vector<std::vector<double>> d_densityMatDerFermiEnergy;
1671
1672 /**
1673 * The indexing of d_eigenVectorsFlattenedHost and
1674 * d_eigenVectorsFlattenedDevice [kPoint * numSpinComponents *
1675 * numLocallyOwnedNodes * numWaveFunctions + iSpin * numLocallyOwnedNodes *
1676 * numWaveFunctions + iNode * numWaveFunctions + iWaveFunction]
1677 */
1681
1685
1686 /// device eigenvectors
1687#ifdef DFTFE_WITH_DEVICE
1690 d_eigenVectorsFlattenedDevice;
1693 d_eigenVectorsDensityMatrixPrimeFlattenedDevice;
1694#endif
1695
1696 /// parallel message stream
1697 dealii::ConditionalOStream pcout;
1698
1699 /// compute-time logger
1700 dealii::TimerOutput computing_timer;
1701 dealii::TimerOutput computingTimerStandard;
1702
1703 /// A plain global timer to track only the total elapsed time after every
1704 /// ground-state solve
1705 dealii::Timer d_globalTimer;
1706
1707 // dft related objects
1708 std::vector<
1712 std::vector<distributedCPUVec<double>> d_densityInNodalValues,
1714 std::vector<distributedCPUVec<double>> d_tauOutNodalValues;
1715 std::vector<
1718
1719 // std::map<dealii::CellId, std::vector<double>> d_phiInValues,
1720 // d_phiOutValues;
1726
1730
1731
1734
1735
1739
1740 std::shared_ptr<AuxDensityMatrix<memorySpace>> d_auxDensityMatrixXCInPtr;
1741 std::shared_ptr<AuxDensityMatrix<memorySpace>> d_auxDensityMatrixXCOutPtr;
1742
1743 // For multipole boundary conditions
1745 std::vector<double> d_dipole;
1746 std::vector<double> d_quadrupole;
1747 std::vector<double> d_smearedChargeMoments;
1749
1750
1751 /// for low rank jacobian inverse approximation
1752 std::deque<distributedCPUVec<double>> d_vcontainerVals;
1753 std::deque<distributedCPUVec<double>> d_fvcontainerVals;
1754 std::deque<distributedCPUVec<double>> d_vSpin0containerVals;
1755 std::deque<distributedCPUVec<double>> d_fvSpin0containerVals;
1756 std::deque<distributedCPUVec<double>> d_vSpin1containerVals;
1757 std::deque<distributedCPUVec<double>> d_fvSpin1containerVals;
1763 static constexpr double d_tikhonovRegularizationConstantLRD = 1.0e-6;
1764
1765 /// for xl-bomd
1766 std::map<dealii::CellId, std::vector<double>> d_rhoAtomsValues,
1768 std::map<dftfe::uInt, std::map<dealii::CellId, std::vector<double>>>
1771
1772 std::vector<
1776
1777 // storage for total electrostatic potential solution vector corresponding
1778 // to input scf electron density
1780
1781 // storage for total electrostatic potential solution vector corresponding
1782 // to output scf electron density
1784
1785 // storage for electrostatic potential Gateaux derivate corresponding
1786 // to electron number preserving electron-density peturbation (required for
1787 // LRDM)
1789
1790 // storage for sum of nuclear electrostatic potential
1792
1793 // storage of densities for xl-bomd
1794 std::deque<distributedCPUVec<double>> d_groundStateDensityHistory;
1795
1796 std::map<dealii::CellId, std::vector<double>> d_pseudoVLoc;
1797
1798 /// Internal data:: map for cell id to Vpseudo local of individual atoms.
1799 /// Only for atoms whose psp tail intersects the local domain.
1800 std::map<dftfe::uInt, std::map<dealii::CellId, std::vector<double>>>
1802
1803
1804 std::vector<std::vector<double>> d_localVselfs;
1805
1806 // nonlocal pseudopotential related objects used only for pseudopotential
1807 // calculation
1808 std::map<dealii::CellId, std::vector<double>> d_rhoCore;
1809
1810 std::map<dealii::CellId, std::vector<double>> d_gradRhoCore;
1811
1812 std::map<dftfe::uInt, std::map<dealii::CellId, std::vector<double>>>
1814
1815 std::map<dealii::CellId, std::vector<double>> d_hessianRhoCore;
1816
1817 std::map<dftfe::uInt, std::map<dealii::CellId, std::vector<double>>>
1819
1820
1821
1822 /// map of atom node number and atomic weight
1823 std::map<dealii::types::global_dof_index, double> d_atomNodeIdToChargeMap;
1824
1825 /// vselfBinsManager object
1827
1828 /// Gateaux derivative of vself field with respect to affine strain tensor
1829 /// components using central finite difference. This is used for cell stress
1830 /// computation
1831 std::vector<distributedCPUVec<double>> d_vselfFieldGateauxDerStrainFDBins;
1832
1833 /// Compute Gateaux derivative of vself field in bins with respect to affine
1834 /// strain tensor components
1835 void
1837
1838 /// dftParameters object
1840
1841 /// kPoint cartesian coordinates
1842 std::vector<double> d_kPointCoordinates;
1843
1844 /// k point crystal coordinates
1845 std::vector<double> d_kPointCoordinatesFrac;
1846
1847 /// k point weights
1848 std::vector<double> d_kPointWeights;
1849
1850 /// closest tria vertex
1851 std::vector<dealii::Point<3>> d_closestTriaVertexToAtomsLocation;
1852 std::vector<dealii::Tensor<1, 3, double>> d_dispClosestTriaVerticesToAtoms;
1853
1854 /// global k index of lower bound of the local k point set
1856 /**
1857 * Recomputes the k point cartesian coordinates from the crystal k point
1858 * coordinates and the current lattice vectors, which can change in each
1859 * ground state solve dutring cell optimization.
1860 */
1861 void
1863
1864 /// fermi energy
1866
1868
1870
1871 /// entropic energy
1873
1874 // chebyshev filter variables and functions
1875 // dftfe::Int numPass ; // number of filter passes
1876
1877 std::vector<double> a0;
1878 std::vector<double> bLow;
1879
1880 /// stores flag for first ever call to chebyshev filtering for a given FEM
1881 /// mesh vector for each k point and spin
1882 std::vector<bool> d_isFirstFilteringCall;
1883
1885
1887
1889
1890
1891 /**
1892 * @ nscf variables
1893 */
1895 void
1897 KohnShamDFTBaseOperator<memorySpace> &kohnShamDFTEigenOperator,
1898 chebyshevOrthogonalizedSubspaceIterationSolver &subspaceIterationSolver);
1899
1900 /**
1901 * @brief compute the maximum of the residual norm of the highest occupied state among all k points
1902 */
1903 double
1905 const std::vector<std::vector<double>>
1906 &residualNormWaveFunctionsAllkPoints,
1907 const std::vector<std::vector<double>> &eigenValuesAllkPoints,
1908 const double _fermiEnergy,
1909 std::vector<double> &maxResidualsAllkPoints);
1910
1911
1912 /**
1913 * @brief compute the maximum of the residual norm of the highest state of interest among all k points
1914 */
1915 double
1917 const std::vector<std::vector<double>>
1918 &residualNormWaveFunctionsAllkPoints,
1919 const std::vector<std::vector<double>> &eigenValuesAllkPoints,
1920 const dftfe::uInt highestState,
1921 std::vector<double> &maxResidualsAllkPoints);
1922
1923
1924
1925#ifdef DFTFE_WITH_DEVICE
1926 void
1928 const dftfe::uInt s,
1929 const dftfe::uInt kPointIndex,
1931 &kohnShamDFTEigenOperator,
1932 elpaScalaManager &elpaScala,
1933 chebyshevOrthogonalizedSubspaceIterationSolverDevice
1934 &subspaceIterationSolverDevice);
1935
1936#endif
1937
1938 void
1940 const dftfe::uInt s,
1941 const dftfe::uInt kPointIndex,
1943 &kohnShamDFTEigenOperator,
1944 elpaScalaManager &elpaScala);
1945
1947 std::shared_ptr<hubbard<dataTypes::number, memorySpace>> d_hubbardClassPtr;
1949 };
1950
1951} // namespace dftfe
1952
1953#endif
Definition AuxDensityMatrix.h:41
Definition KohnShamDFTBaseOperator.h:36
This class performs the anderson mixing in a variable agnostic way This class takes can take differen...
Definition mixingClass.h:57
Definition atomCenteredPostProcessing.h:30
Definition FEBasisOperations.h:87
Concrete class implementing Chebyshev filtered orthogonalized subspace iteration solver.
Definition chebyshevOrthogonalizedSubspaceIterationSolver.h:38
dealii linear solver class wrapper
Definition dealiiLinearSolver.h:32
abstract base class for dft
Definition dftBase.h:34
dftfe::uInt getElectroQuadratureRhsId() const
void applyHomogeneousDirichletBC(const dealii::DoFHandler< 3 > &_dofHandler, const dealii::AffineConstraints< double > &onlyHangingNodeConstraints, dealii::AffineConstraints< double > &constraintMatrix)
Sets homegeneous dirichlet boundary conditions for total potential constraints on non-periodic bounda...
dftfe::utils::MemoryStorage< double, dftfe::utils::MemorySpace::HOST > d_phiOutQuadValues
Definition dft.h:1722
void updateAtomPositionsAndMoveMesh(const std::vector< dealii::Tensor< 1, 3, double > > &globalAtomsDisplacements, const double maxJacobianRatioFactor=1.25, const bool useSingleAtomSolutionsOverride=false)
Updates atom positions, remeshes/moves mesh and calls appropriate reinits.
std::vector< std::vector< double > > d_atomLocationsAutoMesh
Definition dft.h:1326
void compute_fermienergy_purestate(const std::vector< std::vector< double > > &eigenValuesInput, const double numElectronsInput)
find HOMO eigenvalue for pure state
dftfe::uInt d_phiPrimeDofHandlerIndexElectro
Definition dft.h:1500
std::shared_ptr< excManager< memorySpace > > d_excManagerPtr
Definition dft.h:1309
dftfe::uInt d_feOrderPlusOneQuadratureId
Definition dft.h:1496
dealii::AffineConstraints< double > d_constraintsForPhiPrimeElectro
Definition dft.h:1644
std::vector< std::vector< double > > d_meshSizes
Definition dft.h:1325
const std::vector< std::vector< double > > & getCell() const
Gets the current cell lattice vectors.
dealii::DoFHandler< 3 > d_dofHandlerRhoNodal
Definition dft.h:1487
std::vector< std::vector< dftfe::Int > > d_globalChargeIdToImageIdMapTrunc
globalChargeId to ImageChargeId Map generated with a truncated pspCutOff
Definition dft.h:1401
std::vector< double > d_kPointWeights
k point weights
Definition dft.h:1848
std::map< dealii::CellId, std::vector< double > > d_rhoAtomsValues
for xl-bomd
Definition dft.h:1766
void computeOutputDensityDirectionalDerivative(distributedCPUVec< double > &v, distributedCPUVec< double > &vSpin0, distributedCPUVec< double > &vSpin1, distributedCPUVec< double > &fv, distributedCPUVec< double > &fvSpin0, distributedCPUVec< double > &fvSpin1)
void localCollinearMagnetizationDensity(const std::vector< dftfe::utils::MemoryStorage< double, dftfe::utils::MemorySpace::HOST > > &densityQuadValues)
std::map< dealii::CellId, std::vector< dftfe::uInt > > d_bCellNonTrivialAtomImageIds
Definition dft.h:1439
std::vector< orbital > waveFunctionsVector
Definition dft.h:1449
void kohnShamEigenSpaceFirstOrderDensityMatResponse(const dftfe::uInt s, const dftfe::uInt kPointIndex, KohnShamDFTBaseOperator< dftfe::utils::MemorySpace::HOST > &kohnShamDFTEigenOperator, elpaScalaManager &elpaScala)
bool scfConverged
Definition dft.h:1894
std::vector< distributedCPUVec< double > > d_densityInNodalValues
Definition dft.h:1712
void computeRhoNodalFirstOrderResponseFromPSIAndPSIPrime(distributedCPUVec< double > &fv, distributedCPUVec< double > &fvSpin0, distributedCPUVec< double > &fvSpin1)
std::vector< std::vector< double > > d_localVselfs
Definition dft.h:1804
std::vector< const dealii::AffineConstraints< double > * > d_constraintsVectorElectro
Definition dft.h:1555
void saveTriaInfoAndRhoNodalData()
save triangulation information and rho quadrature data to checkpoint file for restarts
dealii::AffineConstraints< double > d_constraintsPRefinedOnlyHanging
Definition dft.h:1650
const std::vector< std::vector< double > > & getLocalVselfs() const
std::vector< std::vector< double > > d_imagePositionsTrunc
Definition dft.h:1398
std::vector< double > d_upperBoundUnwantedSpectrumValues
Definition dft.h:1884
std::shared_ptr< dftfe::linearAlgebra::BLASWrapper< dftfe::utils::MemorySpace::HOST > > d_BLASWrapperPtr
Definition dft.h:1549
const MPI_Comm & getMPIParent() const override
std::vector< dealii::Tensor< 1, 3, double > > d_gaussianMovementAtomsNetDisplacements
Definition dft.h:1473
distributedCPUVec< double > d_residualPredicted
Definition dft.h:1758
dftfe::uInt d_phiTotDofHandlerIndexElectro
Definition dft.h:1499
std::deque< distributedCPUVec< double > > d_fvcontainerVals
Definition dft.h:1753
dealii::AffineConstraints< double > constraintsNoneEigen
Definition dft.h:1639
double d_pspCutOff
distance from the domain till which periodic images will be considered
Definition dft.h:1404
void compute_localizationLength(const std::string &locLengthFileName)
compute localization length
dealii::IndexSet d_locallyRelevantDofsRhoNodal
Definition dft.h:1571
void generateImageCharges(const double pspCutOff, std::vector< dftfe::Int > &imageIds, std::vector< double > &imageCharges, std::vector< std::vector< double > > &imagePositions)
creates datastructures related to periodic image charges
dftfe::uInt d_baseDofHandlerIndexElectro
Definition dft.h:1492
dealii::IndexSet d_locallyOwnedDofsPRefined
Definition dft.h:1572
std::vector< double > d_nearestAtomDistances
nearest atom distances for all domain atoms
Definition dft.h:1366
bool isHubbardCorrectionsUsed()
Function to check if hubbard corrections is being used.
dealii::IndexSet locally_relevant_dofs
Definition dft.h:1570
const dealii::MatrixFree< 3, double > & getMatrixFreeData() const
Get matrix free data object.
void updatePRefinedConstraints()
void initHubbardOperator()
Checks if the Exc functional requires Hubbard correction and sets up the Hubbard class if required.
void moveMeshToAtoms(dealii::Triangulation< 3, 3 > &triangulationMove, dealii::Triangulation< 3, 3 > &triangulationSerial, bool reuseFlag=false, bool moveSubdivided=false)
moves the triangulation vertices using Gaussians such that the all atoms are on triangulation vertice...
void determineOrbitalFilling()
const MPI_Comm & getMPIDomain() const override
dftfe::uInt getSmearedChargeQuadratureIdElectro()
std::set< dftfe::uInt > atomTypes
Definition dft.h:1317
std::deque< distributedCPUVec< double > > d_vcontainerVals
for low rank jacobian inverse approximation
Definition dft.h:1752
std::shared_ptr< dftfe::groupSymmetryClass > groupSymmetryPtr
Definition dft.h:1579
void noRemeshRhoDataInit()
double numElectronsUp
Definition dft.h:1316
std::shared_ptr< dftfe::atomCenteredOrbitalsPostProcessing< dataTypes::number, memorySpace > > d_atomCenteredOrbitalsPostProcessingPtr
Definition dft.h:1541
double d_entropicEnergy
entropic energy
Definition dft.h:1872
dftfe::utils::MemoryStorage< double, dftfe::utils::MemorySpace::HOST > d_gradPhiInQuadValues
Definition dft.h:1724
double numElectronsDown
Definition dft.h:1316
dealii::IndexSet d_locallyOwnedDofsRhoNodal
Definition dft.h:1572
const double d_smearedChargeWidthMin
minimum smeared charge width
Definition dft.h:1447
dispersionCorrection d_dispersionCorr
Definition dft.h:1310
void initializeKohnShamDFTOperator(const bool initializeCublas=true)
void computeVselfFieldGateauxDerFD()
const std::string d_dftfeScratchFolderName
Definition dft.h:1597
dftfe::utils::MemoryStorage< double, dftfe::utils::MemorySpace::HOST > d_gradDensityTotalInValuesLpspQuad
Definition dft.h:1738
dftfe::utils::MemoryStorage< double, dftfe::utils::MemorySpace::HOST > d_gradPhiOutQuadValues
Definition dft.h:1724
std::vector< double > d_smearedChargeMoments
Definition dft.h:1747
void kohnShamEigenSpaceCompute(const dftfe::uInt s, const dftfe::uInt kPointIndex, KohnShamDFTBaseOperator< dftfe::utils::MemorySpace::HOST > &kohnShamDFTEigenOperator, elpaScalaManager &elpaScala, chebyshevOrthogonalizedSubspaceIterationSolver &subspaceIterationSolver, std::vector< double > &residualNormWaveFunctions, const bool computeResidual, const bool useMixedPrec=false, const bool isFirstScf=false)
Function that computes the Eigen space for the Kohn Sham operator.
dealii::AffineConstraints< double > d_constraintsForTotalPotentialElectro
Definition dft.h:1642
std::shared_ptr< AuxDensityMatrix< memorySpace > > d_auxDensityMatrixXCInPtr
Definition dft.h:1740
void writeDomainAndAtomCoordinates(const std::string Path) const
writes the current domain bounding vectors and atom coordinates to files for structural optimization ...
double d_freeEnergyInitial
Definition dft.h:1867
void createMasterChargeIdToImageIdMaps(const double pspCutOff, const std::vector< dftfe::Int > &imageIds, const std::vector< std::vector< double > > &imagePositions, std::vector< std::vector< dftfe::Int > > &globalChargeIdToImageIdMap)
double totalCharge(const dealii::DoFHandler< 3 > &dofHandlerOfField, const std::map< dealii::CellId, std::vector< double > > *rhoQuadValues)
std::vector< double > d_smearedChargeWidths
smeared charge widths for all domain atoms
Definition dft.h:1357
dftUtils::constraintMatrixInfo< dftfe::utils::MemorySpace::HOST > constraintsNoneEigenDataInfo
Definition dft.h:1620
double computeMaximumHighestOccupiedStateResidualNorm(const std::vector< std::vector< double > > &residualNormWaveFunctionsAllkPoints, const std::vector< std::vector< double > > &eigenValuesAllkPoints, const dftfe::uInt highestState, std::vector< double > &maxResidualsAllkPoints)
compute the maximum of the residual norm of the highest state of interest among all k points
dealii::IndexSet locally_owned_dofsEigen
Definition dft.h:1569
std::vector< double > d_gaussianConstantsAutoMesh
Definition dft.h:1347
dftParameters * d_dftParamsPtr
dftParameters object
Definition dft.h:1839
double totalCharge(const dealii::MatrixFree< 3, double > &matrixFreeDataObject, const distributedCPUVec< double > &rhoNodalField)
void compute_rhoOut(const bool isGroundState=false)
Computes output electron-density from wavefunctions.
std::vector< std::vector< dftfe::Int > > d_globalChargeIdToImageIdMap
globalChargeId to ImageChargeId Map
Definition dft.h:1386
std::vector< double > d_flatTopWidthsAutoMeshMove
Definition dft.h:1354
void run()
FIXME: legacy call, move to main.cc.
std::map< dftfe::uInt, std::map< dealii::CellId, std::vector< double > > > d_gradRhoAtomsValuesSeparate
Definition dft.h:1769
dealii::TimerOutput computing_timer
compute-time logger
Definition dft.h:1700
distributedCPUVec< double > d_preCondTotalDensityResidualVector
Definition dft.h:1728
std::shared_ptr< dftfe::linearAlgebra::BLASWrapper< dftfe::utils::MemorySpace::HOST > > getBLASWrapperHost()
const distributedCPUVec< double > & getRhoNodalSplitOut() const
void writeBands()
write the KS eigen values for given BZ sampling/path
void projectPreviousGroundStateRho()
project ground state electron density from previous mesh into the new mesh to be used as initial gues...
distributedCPUVec< double > d_rhoOutNodalValuesSplit
Definition dft.h:1727
void loadDensityFromQuadratureValues()
std::map< dealii::CellId, std::vector< double > > d_bQuadValuesAllAtoms
non-intersecting smeared charges of all atoms at quad points
Definition dft.h:1414
double computeMaximumHighestOccupiedStateResidualNorm(const std::vector< std::vector< double > > &residualNormWaveFunctionsAllkPoints, const std::vector< std::vector< double > > &eigenValuesAllkPoints, const double _fermiEnergy, std::vector< double > &maxResidualsAllkPoints)
compute the maximum of the residual norm of the highest occupied state among all k points
void writeMesh()
Writes inital density and mesh to file.
std::map< dealii::CellId, std::vector< double > > & getBQuadValuesAllAtoms()
non-intersecting smeared charges of all atoms at quad points
const std::vector< dealii::types::global_dof_index > & getLocalDofIndicesImag() const
Get local dofs global indices imag.
distributedCPUVec< double > d_rhoOutNodalValuesDistributed
Definition dft.h:1729
std::shared_ptr< AuxDensityMatrix< memorySpace > > d_auxDensityMatrixXCOutPtr
Definition dft.h:1741
void totalNonCollinearMagnetization(const std::vector< dftfe::utils::MemoryStorage< double, dftfe::utils::MemorySpace::HOST > > &densityQuadValues)
std::map< dftfe::uInt, std::map< dealii::CellId, std::vector< double > > > d_hessianRhoAtomsValuesSeparate
Definition dft.h:1770
dftUtils::constraintMatrixInfo< dftfe::utils::MemorySpace::HOST > constraintsNoneDataInfo
Definition dft.h:1630
std::vector< dftfe::utils::MemoryStorage< double, dftfe::utils::MemorySpace::HOST > > & getDensityOutValues()
distributedCPUVec< double > d_rhoInNodalValuesRead
Definition dft.h:1727
std::map< dftfe::uInt, dftfe::uInt > d_atomTypeAtributes
Definition dft.h:1321
void writeMesh(std::string meshFileName)
void resetRhoNodalIn(distributedCPUVec< double > &OutDensity)
dealii::IndexSet locally_owned_dofs
Definition dft.h:1569
std::vector< dealii::Tensor< 1, 3, double > > d_dispClosestTriaVerticesToAtoms
Definition dft.h:1852
std::map< dftfe::uInt, std::map< dftfe::uInt, std::map< dftfe::uInt, double > > > outerValues
Definition dft.h:1455
std::deque< distributedCPUVec< double > > d_fvSpin1containerVals
Definition dft.h:1757
std::shared_ptr< dftfe::oncvClass< dataTypes::number, memorySpace > > d_oncvClassPtr
Definition dft.h:1537
const MPI_Comm & getMPIInterBand() const override
const distributedCPUVec< double > & getRhoNodalOut() const
std::map< dealii::CellId, std::vector< dftfe::Int > > d_bQuadAtomIdsAllAtomsImages
Definition dft.h:1425
double computeAndPrintKE(dftfe::utils::MemoryStorage< double, dftfe::utils::MemorySpace::HOST > &kineticEnergyDensityValues)
Computes the kinetic energy.
double d_autoMeshMaxJacobianRatio
Definition dft.h:1462
double d_dftfeClassStartTime
Definition dft.h:1509
std::map< dealii::CellId, std::vector< double > > d_hessianRhoCore
Definition dft.h:1815
elpaScalaManager * d_elpaScala
Definition dft.h:1582
void totalMagnetization(const dftfe::utils::MemoryStorage< double, dftfe::utils::MemorySpace::HOST > &magQuadValues)
Computes net magnetization from the difference of local spin densities.
double d_monopole
Definition dft.h:1744
void initAtomicRho()
dftfe::uInt d_kerkerAXQuadratureIdElectro
Definition dft.h:1502
std::vector< dftfe::utils::MemoryStorage< double, dftfe::utils::MemorySpace::HOST > > d_densityInQuadValues
Definition dft.h:1710
std::map< dealii::CellId, std::vector< double > > d_gradRhoCore
Definition dft.h:1810
void createpRefinedDofHandler(dealii::parallel::distributed::Triangulation< 3 > &triangulation)
double numElectrons
Definition dft.h:1316
double computeResidualQuadData(const dftfe::utils::MemoryStorage< double, dftfe::utils::MemorySpace::HOST > &outValues, const dftfe::utils::MemoryStorage< double, dftfe::utils::MemorySpace::HOST > &inValues, dftfe::utils::MemoryStorage< double, dftfe::utils::MemorySpace::HOST > &residualValues, const dftfe::utils::MemoryStorage< double, dftfe::utils::MemorySpace::HOST > &JxW, const bool computeNorm)
Copies the residual residualValues=outValues-inValues.
dataTypes::number computeTraceXtHX(dftfe::uInt numberWaveFunctionsEstimate)
bool d_smearedChargeMomentsComputed
Definition dft.h:1748
distributedCPUVec< double > d_rhoNodalFieldRefined
Definition dft.h:1728
bool d_tolReached
Definition dft.h:1762
chebyshevOrthogonalizedSubspaceIterationSolver d_subspaceIterationSolver
Definition dft.h:1603
const MPI_Comm mpi_communicator
Definition dft.h:1560
std::vector< dealii::types::global_dof_index > local_dof_indicesImag
Definition dft.h:1574
std::vector< std::vector< double > > d_domainBoundingVectors
Definition dft.h:1325
dftfe::utils::MemoryStorage< dataTypes::number, dftfe::utils::MemorySpace::HOST > d_eigenVectorsDensityMatrixPrimeHost
Definition dft.h:1684
std::vector< double > d_imageCharges
Definition dft.h:1379
void l2ProjectionQuadToNodal(const std::shared_ptr< dftfe::basis::FEBasisOperations< double, double, dftfe::utils::MemorySpace::HOST > > &basisOperationsPtr, const dealii::AffineConstraints< double > &constraintMatrix, const dftfe::uInt dofHandlerId, const dftfe::uInt quadratureId, const dftfe::utils::MemoryStorage< double, dftfe::utils::MemorySpace::HOST > &quadratureValueData, distributedCPUVec< double > &nodalField)
l2 projection
const dftfe::utils::MemoryStorage< dataTypes::number, dftfe::utils::MemorySpace::HOST > & getEigenVectorsHost() const
Get reference the host eigen vectors.
dftfe::uInt getDensityQuadratureId()
dealii::IndexSet d_locallyRelevantDofsPRefined
Definition dft.h:1571
std::deque< distributedCPUVec< double > > d_vSpin1containerVals
Definition dft.h:1756
std::map< dealii::CellId, std::vector< dftfe::Int > > d_bQuadAtomIdsAllAtoms
non-intersecting smeared charges atom ids of all atoms at quad points
Definition dft.h:1420
void computeRhoNodalFromPSI()
computes density nodal data from wavefunctions
dftfe::uInt d_nOMPThreads
Definition dft.h:1508
std::vector< dftfe::Int > d_imageIds
Definition dft.h:1373
std::map< dealii::CellId, std::vector< double > > d_hessianRhoAtomsValues
Definition dft.h:1767
dealii::DoFHandler< 3 > dofHandler
Definition dft.h:1486
dealii::ConditionalOStream pcout
device eigenvectors
Definition dft.h:1697
std::vector< dftfe::uInt > d_nearestAtomIds
nearest atom ids for all domain atoms
Definition dft.h:1363
dftfe::uInt d_gllQuadratureId
Definition dft.h:1498
std::vector< std::vector< double > > d_densityMatDerFermiEnergy
Definition dft.h:1670
std::vector< dftfe::utils::MemoryStorage< double, dftfe::utils::MemorySpace::HOST > > d_densityOutQuadValues
Definition dft.h:1710
std::vector< dealii::Point< 3 > > d_controlPointLocationsCurrentMove
Definition dft.h:1474
std::vector< dftfe::utils::MemoryStorage< double, dftfe::utils::MemorySpace::HOST > > d_tauResidualQuadValues
Definition dft.h:1717
std::vector< double > d_smearedChargeScaling
smeared charge normalization scaling for all domain atoms
Definition dft.h:1360
dftfe::uInt getElectroQuadratureAxId() const
std::map< dealii::CellId, std::vector< double > > d_gradbQuadValuesAllAtoms
non-intersecting smeared charge gradients of all atoms at quad points
Definition dft.h:1417
void compute_fermienergy_constraintMagnetization_purestate(const std::vector< std::vector< double > > &eigenValuesInput)
Find spin-up and spin-down channel HOMO eigenvalues.
std::map< dealii::types::global_dof_index, dealii::Point< 3 > > d_supportPoints
Definition dft.h:1551
double getInternalEnergy() const
std::vector< distributedCPUVec< double > > d_densityResidualNodalValues
Definition dft.h:1713
std::map< dealii::CellId, std::vector< double > > d_pseudoVLoc
Definition dft.h:1796
void setNumElectrons(dftfe::uInt inputNumElectrons)
void applyKerkerPreconditionerToTotalDensityResidual(kerkerSolverProblemWrapperClass &kerkerPreconditionedResidualSolverProblem, dealiiLinearSolver &CGSolver, distributedCPUVec< double > &residualRho, distributedCPUVec< double > &preCondTotalDensityResidualVector)
Mixing schemes for mixing electron-density.
void compute_fermienergy(const std::vector< std::vector< double > > &eigenValuesInput, const double numElectronsInput)
Computes Fermi-energy obtained by imposing constraint on the number of electrons.
dftfe::uInt d_autoMesh
Definition dft.h:1463
void loadQuadratureData(const std::shared_ptr< dftfe::basis::FEBasisOperations< dataTypes::number, double, dftfe::utils::MemorySpace::HOST > > &basisOperationsPtr, const dftfe::uInt quadratureId, dftfe::utils::MemoryStorage< double, dftfe::utils::MemorySpace::HOST > &quadratureValueData, const dftfe::uInt fieldDimension, const std::string &fieldName, const std::string &folderPath, const MPI_Comm &mpi_comm_parent, const MPI_Comm &mpi_comm_domain, const MPI_Comm &interpoolcomm, const MPI_Comm &interBandGroupComm)
loads data from quad points of checkpoint file. Used for restart calculations, nscf and bands.
std::vector< bool > d_isFirstFilteringCall
Definition dft.h:1882
std::map< dealii::types::global_dof_index, dealii::Point< 3 > > d_supportPointsPRefined
Definition dft.h:1552
dealii::AffineConstraints< double > * getDensityConstraint()
std::vector< dftfe::Int > d_imageIdsTrunc
Definition dft.h:1390
const MPI_Comm interBandGroupComm
Definition dft.h:1566
const std::vector< std::vector< double > > & getImageAtomLocationsCart() const
Gets the current image atom Locations in cartesian form (origin at center of domain) from dftClass.
double d_wfcInitTruncation
init wfc trunctation radius
Definition dft.h:1480
dftfe::uInt d_densityQuadratureIdElectro
Definition dft.h:1506
dftClass(const MPI_Comm &mpiCommParent, const MPI_Comm &mpi_comm_domain, const MPI_Comm &interpoolcomm, const MPI_Comm &interBandGroupComm, const std::string &scratchFolderName, dftParameters &dftParams)
dftClass constructor
std::vector< double > d_kPointCoordinatesFrac
k point crystal coordinates
Definition dft.h:1845
void applyPeriodicBCHigherOrderNodes()
Computes inner Product and Y = alpha*X + Y for complex vectors used during periodic boundary conditio...
double getCellVolume() const
Gets the current cell volume.
std::vector< double > d_quadrupole
Definition dft.h:1746
double rhofieldInnerProduct(const dealii::MatrixFree< 3, double > &matrixFreeDataObject, const distributedCPUVec< double > &rhoNodalField1, const distributedCPUVec< double > &rhoNodalField2, const dftfe::uInt dofHandlerId, const dftfe::uInt quadratureId)
void reInitializeKohnShamDFTOperator()
void initUnmovedTriangulation(dealii::parallel::distributed::Triangulation< 3 > &triangulation)
dealii::AffineConstraints< double > d_noConstraints
Definition dft.h:1640
distributedCPUVec< double > d_phiExt
Definition dft.h:1791
std::shared_ptr< dftfe::basis::FEBasisOperations< dataTypes::number, double, dftfe::utils::MemorySpace::HOST > > getBasisOperationsHost()
std::vector< std::vector< double > > d_imagePositionsAutoMesh
Definition dft.h:1327
double fermiEnergy
fermi energy
Definition dft.h:1865
const dealii::MatrixFree< 3, double > & getMatrixFreeDataElectro() const
virtual void writeGSElectronDensity(const std::string Path) const
writes quadrature grid information and associated spin-up and spin-down electron-density for post-pro...
const dealii::AffineConstraints< double > * getConstraintsVectorElectro()
void compute_fermienergy_constraintMagnetization(const std::vector< std::vector< double > > &eigenValuesInput)
Computes Fermi-energy obtained by imposing separate constraints on the number of spin-up and spin-dow...
const MPI_Comm & getMPIInterPool() const override
std::vector< dftfe::utils::MemoryStorage< double, dftfe::utils::MemorySpace::HOST > > d_tauInQuadValues
Definition dft.h:1717
const std::vector< dealii::types::global_dof_index > & getLocalDofIndicesReal() const
Get local dofs global indices real.
std::vector< distributedCPUVec< double > > d_tauOutNodalValues
Definition dft.h:1714
void initpRefinedObjects(const bool recomputeBasisData, const bool meshOnlyDeformed, const bool vselfPerturbationUpdateForStress=false)
dftfe::uInt d_binsStartDofHandlerIndexElectro
Definition dft.h:1504
std::vector< dealii::Tensor< 1, 3, double > > d_atomsDisplacementsGaussianRead
Gaussian displacements of atoms read from file.
Definition dft.h:1330
std::vector< double > bLow
Definition dft.h:1878
const MPI_Comm d_mpiCommParent
Definition dft.h:1564
std::vector< bool > selectedDofsHanging
Definition dft.h:1577
std::shared_ptr< hubbard< dataTypes::number, memorySpace > > getHubbardClassPtr()
Returns the shared ptr to hubbard class.
std::map< dealii::types::global_dof_index, dealii::Point< 3 > > d_supportPointsEigen
Definition dft.h:1552
std::vector< std::vector< double > > eigenValues
Definition dft.h:1662
void writeStructureEnergyForcesDataPostProcess(const std::string Path) const
writes atomistics data for subsequent post-processing. Related to WRITE STRUCTURE ENERGY FORCES DATA ...
const std::vector< dftfe::Int > & getImageAtomIDs() const
Gets the current image atom ids from dftClass.
std::shared_ptr< dftfe::basis::FEBasisOperations< double, double, dftfe::utils::MemorySpace::HOST > > getBasisOperationsElectroHost()
distributedCPUVec< double > d_phiPrime
Definition dft.h:1788
dftfe::uInt getNumEigenValues() const
double lowrankApproxScfDielectricMatrixInv(const dftfe::uInt scfIter)
dealii::IndexSet locally_relevant_dofsEigen
Definition dft.h:1570
void computeRhoInitialGuessFromPSI(std::vector< std::vector< distributedCPUVec< double > > > eigenVectors)
const dftfe::utils::MemoryStorage< dataTypes::number, memorySpace > & getEigenVectors() const
Get reference the memorySpace templated eigen vectors.
void compute_tdos(const std::vector< std::vector< double > > &eigenValuesInput, const std::string &fileName)
compute density of states and local density of states
std::map< dftfe::uInt, std::map< dftfe::uInt, std::map< dftfe::uInt, alglib::spline1dinterpolant > > > radValues
Definition dft.h:1453
double getFermiEnergy() const
Get the value of fermi energy.
dftfe::uInt d_nonAtomicWaveFunctions
Number of random wavefunctions.
Definition dft.h:289
dftfe::uInt d_rankCurrentLRD
Definition dft.h:1759
dftfe::uInt numLevels
Definition dft.h:1315
void writeDomainAndAtomCoordinates()
writes the current domain bounding vectors and atom coordinates to files, which are required for geom...
double d_freeEnergy
Definition dft.h:1869
dealii::DoFHandler< 3 > d_dofHandlerPRefined
Definition dft.h:1486
void set()
atomic system pre-processing steps.
triangulationManager * getTriangulationManager()
void saveQuadratureData(const std::shared_ptr< dftfe::basis::FEBasisOperations< dataTypes::number, double, dftfe::utils::MemorySpace::HOST > > &basisOperationsPtr, const dftfe::uInt quadratureId, const dftfe::utils::MemoryStorage< double, dftfe::utils::MemorySpace::HOST > &quadratureValueData, const dftfe::uInt fieldDimension, const std::string &fieldName, const std::string &folderPath, const MPI_Comm &mpi_comm_parent, const MPI_Comm &mpi_comm_domain, const MPI_Comm &interpoolcomm, const MPI_Comm &interBandGroupComm)
save data of quad points to checkpoint file. Used for restart calculations, nscf and bands.
void initLocalPseudoPotential(const dealii::DoFHandler< 3 > &_dofHandler, const dftfe::uInt lpspQuadratureId, const dealii::MatrixFree< 3, double > &_matrix_free_data, const dftfe::uInt _phiExtDofHandlerIndex, const dealii::AffineConstraints< double > &phiExtConstraintMatrix, const std::map< dealii::types::global_dof_index, dealii::Point< 3 > > &supportPoints, const vselfBinsManager &vselfBinManager, distributedCPUVec< double > &phiExt, std::map< dealii::CellId, std::vector< double > > &_pseudoValues, std::map< dftfe::uInt, std::map< dealii::CellId, std::vector< double > > > &_pseudoValuesAtoms)
dealii::AffineConstraints< double > constraintsNone
Definition dft.h:1639
double d_relativeErrorJacInvApproxPrevScfLRD
Definition dft.h:1760
std::vector< std::map< dealii::CellId, std::vector< dftfe::uInt > > > d_bCellNonTrivialAtomIdsBins
Definition dft.h:1434
std::shared_ptr< dftfe::basis::FEBasisOperations< dataTypes::number, double, dftfe::utils::MemorySpace::HOST > > d_basisOperationsPtrHost
Definition dft.h:1515
std::vector< dealii::types::global_dof_index > localProc_dof_indicesReal
Definition dft.h:1575
const std::vector< dealii::types::global_dof_index > & getLocalProcDofIndicesReal() const
Get local dofs local proc indices real.
const std::vector< std::vector< double > > & getAtomLocationsFrac() const
Gets the current atom Locations in fractional form from dftClass (only applicable for periodic and se...
dftfe::utils::MemoryStorage< dataTypes::number, dftfe::utils::MemorySpace::HOST > d_eigenVectorsFlattenedHost
Definition dft.h:1680
distributedCPUVec< double > d_magYInNodalValuesRead
Definition dft.h:1732
dftfe::utils::MemoryStorage< double, dftfe::utils::MemorySpace::HOST > d_phiInQuadValues
Definition dft.h:1722
std::shared_ptr< dftfe::linearAlgebra::BLASWrapper< dftfe::utils::MemorySpace::HOST > > d_BLASWrapperPtrHost
Definition dft.h:1534
std::vector< double > d_dipole
Definition dft.h:1745
double d_domainVolume
volume of the domain
Definition dft.h:1477
void readkPointData()
void initImageChargesUpdateKPoints(bool flag=true)
generate image charges and update k point cartesian coordinates based on current lattice vectors
std::vector< double > d_generatorFlatTopWidths
composite generator flat top widths for all domain atoms
Definition dft.h:1350
expConfiningPotential d_expConfiningPot
Definition dft.h:1946
void initBoundaryConditions(const bool recomputeBasisData=true, const bool meshOnlyDeformed=false, const bool vselfPerturbationUpdateForStress=false)
distributedCPUVec< double > d_magZInNodalValuesRead
Definition dft.h:1732
void outputDensity()
write electron density solution fields
void runFunctionalTest()
double d_nlPSPCutOff
Definition dft.h:1411
std::vector< dftfe::utils::MemoryStorage< double, dftfe::utils::MemorySpace::HOST > > d_densityResidualQuadValues
Definition dft.h:1711
void solveNoSCF()
compute approximation to ground-state without solving the SCF iteration
const std::vector< dealii::types::global_dof_index > & getLocalProcDofIndicesImag() const
Get local dofs local proc indices imag.
void locatePeriodicPinnedNodes(const dealii::DoFHandler< 3 > &_dofHandler, const dealii::AffineConstraints< double > &constraintMatrixBase, dealii::AffineConstraints< double > &constraintMatrix)
Sets homogeneous dirichlet boundary conditions on a node farthest from all atoms (pinned node)....
dealii::Timer d_globalTimer
Definition dft.h:1705
double getEntropicEnergy() const
const std::vector< double > & getNearestAtomDistance() const
Gets the nearest atom distance for each atom.
std::vector< dealii::types::global_dof_index > local_dof_indicesReal
Definition dft.h:1573
poissonSolverProblemWrapperClass d_phiTotalSolverProblem
Definition dft.h:1584
dealii::MatrixFree< 3, double > matrix_free_data
Definition dft.h:1510
dealii::AffineConstraints< double > d_constraintsRhoNodalOnlyHanging
Definition dft.h:1654
double fermiEnergyDown
Definition dft.h:1865
std::vector< std::vector< double > > d_fracOccupancy
Definition dft.h:1668
std::vector< dealii::Point< 3 > > d_closestTriaVertexToAtomsLocation
closest tria vertex
Definition dft.h:1851
std::map< dftfe::uInt, std::map< dealii::CellId, std::vector< double > > > d_gradRhoCoreAtoms
Definition dft.h:1813
bool d_isAtomsGaussianDisplacementsReadFromFile
Definition dft.h:1338
poissonSolverProblemWrapperClass d_phiPrimeSolverProblem
Definition dft.h:1586
double d_residualNormPredicted
Definition dft.h:1761
static constexpr double d_tikhonovRegularizationConstantLRD
Definition dft.h:1763
std::vector< dftfe::utils::MemoryStorage< double, dftfe::utils::MemorySpace::HOST > > d_gradDensityOutQuadValues
Definition dft.h:1774
dftfe::uInt d_lpspQuadratureIdElectro
Definition dft.h:1497
dftfe::uInt d_phiExtDofHandlerIndexElectro
Definition dft.h:1488
dftfe::utils::MemoryStorage< double, dftfe::utils::MemorySpace::HOST > d_densityTotalInValuesLpspQuad
Definition dft.h:1737
std::shared_ptr< dftfe::linearAlgebra::BLASWrapper< memorySpace > > getBLASWrapperMemSpace()
dftfe::uInt getElectroDofHandlerIndex() const
const std::set< dftfe::uInt > & getAtomTypes() const
Gets the current atom types from dftClass.
std::map< dftfe::uInt, std::map< dealii::CellId, std::vector< double > > > d_rhoAtomsValuesSeparate
Definition dft.h:1769
void outputWfc(const std::string outputFileName="wfcOutput")
write wavefunction solution fields
dealii::AffineConstraints< double > d_constraintsPRefined
Definition dft.h:1648
std::map< dealii::CellId, std::vector< double > > d_gradRhoAtomsValues
Definition dft.h:1767
std::map< dealii::CellId, std::vector< dftfe::uInt > > d_bCellNonTrivialAtomIds
Definition dft.h:1429
double getFreeEnergy() const
std::vector< std::vector< double > > d_reciprocalLatticeVectors
Definition dft.h:1325
MixingScheme d_mixingScheme
Definition dft.h:1725
std::map< dealii::types::global_dof_index, double > & getAtomNodeToChargeMap()
map of atom node number and atomic weight
void initElectronicFields()
std::shared_ptr< dftfe::basis::FEBasisOperations< double, double, dftfe::utils::MemorySpace::HOST > > d_basisOperationsPtrElectroHost
Definition dft.h:1519
meshMovementAffineTransform d_affineTransformMesh
affine transformation object
Definition dft.h:1467
std::vector< double > a0
Definition dft.h:1877
void l2ProjectionQuadDensityMinusAtomicDensity(const std::shared_ptr< dftfe::basis::FEBasisOperations< double, double, dftfe::utils::MemorySpace::HOST > > &basisOperationsPtr, const dealii::AffineConstraints< double > &constraintMatrix, const dftfe::uInt dofHandlerId, const dftfe::uInt quadratureId, const dftfe::utils::MemoryStorage< double, dftfe::utils::MemorySpace::HOST > &quadratureValueData, distributedCPUVec< double > &nodalField)
l2 projection
std::vector< std::vector< double > > d_imagePositions
Definition dft.h:1383
std::shared_ptr< dftfe::basis::FEBasisOperations< dataTypes::number, double, memorySpace > > getBasisOperationsMemSpace()
triangulationManager d_mesh
Definition dft.h:1460
std::vector< dftfe::utils::MemoryStorage< double, dftfe::utils::MemorySpace::HOST > > & getDensityInValues()
std::vector< std::map< dealii::CellId, std::vector< dftfe::uInt > > > d_bCellNonTrivialAtomImageIdsBins
Definition dft.h:1444
double fermiEnergyUp
Definition dft.h:1865
std::map< dealii::types::global_dof_index, double > d_atomNodeIdToChargeMap
map of atom node number and atomic weight
Definition dft.h:1823
std::shared_ptr< hubbard< dataTypes::number, memorySpace > > d_hubbardClassPtr
Definition dft.h:1947
void applyMultipoleDirichletBC(const dealii::DoFHandler< 3 > &_dofHandler, const dealii::AffineConstraints< double > &onlyHangingNodeConstraints, dealii::AffineConstraints< double > &constraintMatrix)
Sets inhomegeneous dirichlet boundary conditions upto quadrupole for total potential constraints on n...
dftfe::uInt d_lpspQuadratureId
Definition dft.h:1495
elpaScalaManager * getElpaScalaManager() const
const std::vector< double > & getKPointWeights() const
std::vector< double > d_imageChargesTrunc
Definition dft.h:1394
void deformDomain(const dealii::Tensor< 2, 3, double > &deformationGradient, const bool vselfPerturbationUpdateForStress=false, const bool useSingleAtomSolutionsOverride=false, const bool print=true)
Deforms the domain by the given deformation gradient and reinitializes the dftClass datastructures.
void loadTriaInfoAndRhoNodalData()
load triangulation information rho quadrature data from checkpoint file for restarted run
std::deque< distributedCPUVec< double > > d_fvSpin0containerVals
Definition dft.h:1755
bool d_useHubbard
Definition dft.h:1948
dftParameters & getParametersObject() const
Get reference to dftParameters object.
dftfe::uInt d_phiTotAXQuadratureIdElectro
Definition dft.h:1501
dealii::AffineConstraints< double > d_constraintsRhoNodal
Definition dft.h:1652
double fieldGradl2Norm(const dealii::MatrixFree< 3, double > &matrixFreeDataObject, const distributedCPUVec< double > &field)
void normalizeRhoInQuadValues()
normalize the input electron density
double totalCharge(const dealii::DoFHandler< 3 > &dofHandlerOfField, const distributedCPUVec< double > &rhoNodalField)
Computes total charge by integrating the electron-density.
void compute_ldos(const std::vector< std::vector< double > > &eigenValuesInput, const std::string &fileName)
virtual void resetRhoNodalSplitIn(distributedCPUVec< double > &OutDensity)
distributedCPUVec< double > d_tempEigenVec
Definition dft.h:1886
dealii::MatrixFree< 3, double > d_matrixFreeDataPRefined
Definition dft.h:1510
const std::vector< std::vector< double > > & getAtomLocationsCart() const
Gets the current atom Locations in cartesian form (origin at center of domain) from dftClass.
std::tuple< bool, double > solve(const bool computeForces=true, const bool computestress=true, const bool restartGroundStateCalcFromChk=false)
Kohn-Sham ground-state solve using SCF iteration.
dftfe::uInt d_nonPeriodicDensityDofHandlerIndexElectro
Definition dft.h:1491
double computeTraceXtKX(dftfe::uInt numberWaveFunctionsEstimate)
void normalizeRhoOutQuadValues()
normalize the output total electron density in each scf
double d_groundStateEnergy
Definition dft.h:1865
dftfe::uInt d_smearedChargeQuadratureIdElectro
Definition dft.h:1493
double computeVolume(const dealii::DoFHandler< 3 > &_dofHandler)
Computes the volume of the domain.
const std::map< dealii::CellId, std::vector< double > > & getPseudoVLoc() const
return the pseudo potential field
dftfe::utils::MemoryStorage< double, dftfe::utils::MemorySpace::HOST > d_gradDensityTotalOutValuesLpspQuad
Definition dft.h:1738
const std::map< dealii::CellId, std::vector< dftfe::uInt > > & getbCellNonTrivialAtomIds() const
double lowrankApproxScfDielectricMatrixInvSpinPolarized(const dftfe::uInt scfIter)
dftfe::uInt d_nlpspQuadratureId
Definition dft.h:1494
double getNumElectrons() const
Get the number of electrons.
dealii::FESystem< 3 > FEEigen
Definition dft.h:1485
void recomputeKPointCoordinates()
distributedCPUVec< double > d_magXInNodalValuesRead
Definition dft.h:1733
int lowerBoundKindex
global k index of lower bound of the local k point set
Definition dft.h:1855
std::vector< distributedCPUVec< double > > d_densityOutNodalValues
Definition dft.h:1713
dealii::FESystem< 3 > FE
Definition dft.h:1485
void trivialSolveForStress()
bool d_kohnShamDFTOperatorsInitialized
Definition dft.h:1593
std::vector< const dealii::AffineConstraints< double > * > d_constraintsVector
Definition dft.h:1553
std::shared_ptr< configurationalForceClass< memorySpace > > d_configForcePtr
Definition dft.h:1580
dftfe::uInt getDensityDofHandlerIndex()
get the index of the DoF Handler corresponding to
std::vector< std::vector< double > > atomLocations
FIXME: remove atom type atributes from atomLocations.
Definition dft.h:1324
void updateAuxDensityXCMatrix(const std::vector< dftfe::utils::MemoryStorage< double, dftfe::utils::MemorySpace::HOST > > &densityQuadValues, const std::vector< dftfe::utils::MemoryStorage< double, dftfe::utils::MemorySpace::HOST > > &gradDensityQuadValues, const std::vector< dftfe::utils::MemoryStorage< double, dftfe::utils::MemorySpace::HOST > > &tauQuadValues, const std::map< dealii::CellId, std::vector< double > > &rhoCore, const std::map< dealii::CellId, std::vector< double > > &gradRhoCore, const dftfe::utils::MemoryStorage< dataTypes::number, memorySpace > &eigenVectorsFlattenedMemSpace, const std::vector< std::vector< double > > &eigenValues, const double fermiEnergy_, const double fermiEnergyUp_, const double fermiEnergyDown_, std::shared_ptr< AuxDensityMatrix< memorySpace > > auxDensityMatrixXCPtr)
std::vector< dftfe::utils::MemoryStorage< double, dftfe::utils::MemorySpace::HOST > > d_tauOutQuadValues
Definition dft.h:1717
KohnShamDFTBaseOperator< memorySpace > * getKohnShamDFTBaseOperatorClass()
get the Ptr to the operator class ( Kohn Sham Base Operator)
std::vector< dftfe::utils::MemoryStorage< double, dftfe::utils::MemorySpace::HOST > > d_gradDensityResidualQuadValues
Definition dft.h:1775
void addAtomicRhoQuadValuesGradients(dftfe::utils::MemoryStorage< double, dftfe::utils::MemorySpace::HOST > &quadratureValueData, dftfe::utils::MemoryStorage< double, dftfe::utils::MemorySpace::HOST > &quadratureGradValueData, const bool isConsiderGradData=false)
add atomic densities at quadrature points
dealii::TimerOutput computingTimerStandard
Definition dft.h:1701
void normalizeAtomicRhoQuadValues()
normalize the electron density
void solveBands()
compute bands without solving the SCF iteration
const dftfe::uInt n_mpi_processes
Definition dft.h:1567
std::map< dftfe::uInt, std::map< dealii::CellId, std::vector< double > > > d_pseudoVLocAtoms
Definition dft.h:1801
double d_atomicRhoScalingFac
Definition dft.h:988
std::vector< dftfe::utils::MemoryStorage< double, dftfe::utils::MemorySpace::HOST > > d_gradDensityInQuadValues
Definition dft.h:1774
void localNonCollinearMagnetizationDensity(const std::vector< dftfe::utils::MemoryStorage< double, dftfe::utils::MemorySpace::HOST > > &densityQuadValues)
std::vector< dealii::types::global_dof_index > localProc_dof_indicesImag
Definition dft.h:1576
std::vector< double > d_netFloatingDispSinceLastCheckForSmearedChargeOverlaps
Definition dft.h:1336
std::vector< double > d_kPointCoordinates
kPoint cartesian coordinates
Definition dft.h:1842
std::map< dealii::CellId, std::vector< double > > d_rhoCore
Definition dft.h:1808
KohnShamDFTBaseOperator< memorySpace > * d_kohnShamDFTOperatorPtr
Definition dft.h:1595
dftfe::uInt d_numEigenValues
Number of Kohn-Sham eigen values to be computed.
Definition dft.h:283
const dealii::Tensor< 2, 3, double > & getCellStress() const
Gets the current cell stress from dftClass.
const dftfe::uInt this_mpi_process
Definition dft.h:1568
void computeRhoNodalMassVector(dftfe::utils::MemoryStorage< double, dftfe::utils::MemorySpace::HOST > &massVec)
Computes the diagonal mass matrix for rho nodal grid, used for nodal mixing.
const std::vector< double > & getForceonAtoms() const
Gets the current atomic forces from dftClass.
void nscf(KohnShamDFTBaseOperator< memorySpace > &kohnShamDFTEigenOperator, chebyshevOrthogonalizedSubspaceIterationSolver &subspaceIterationSolver)
double rhofieldl2Norm(const dealii::MatrixFree< 3, double > &matrixFreeDataObject, const distributedCPUVec< double > &rhoNodalField, const dftfe::uInt dofHandlerId, const dftfe::uInt quadratureId)
double computeResidualNodalData(const distributedCPUVec< double > &outValues, const distributedCPUVec< double > &inValues, distributedCPUVec< double > &residualValues)
void finalizeKohnShamDFTOperator()
std::vector< double > d_netFloatingDispSinceLastBinsUpdate
Definition dft.h:1333
distributedCPUVec< double > d_phiTotRhoOut
Definition dft.h:1783
bool d_isRestartGroundStateCalcFromChk
Definition dft.h:1888
dftUtils::constraintMatrixInfo< dftfe::utils::MemorySpace::HOST > d_constraintsRhoNodalInfo
Definition dft.h:1657
dftfe::uInt d_densityDofHandlerIndexElectro
Definition dft.h:1490
void locateAtomCoreNodes(const dealii::DoFHandler< 3 > &_dofHandler, std::map< dealii::types::global_dof_index, double > &atomNodeIdToChargeValueMap)
Finds the global dof ids of the nodes containing atoms.
const MPI_Comm interpoolcomm
Definition dft.h:1565
double totalCharge(const dealii::DoFHandler< 3 > &dofHandlerOfField, const dftfe::utils::MemoryStorage< double, dftfe::utils::MemorySpace::HOST > &rhoQuadValues)
void initPseudoPotentialAll(const bool updateNonlocalSparsity=true)
dftfe::uInt d_helmholtzDofHandlerIndexElectro
Definition dft.h:1503
void computeMultipoleMoments(const std::shared_ptr< dftfe::basis::FEBasisOperations< double, double, dftfe::utils::MemorySpace::HOST > > &basisOperationsPtr, const dftfe::uInt densityQuadratureId, const dftfe::utils::MemoryStorage< double, dftfe::utils::MemorySpace::HOST > &rhoQuadValues, const std::map< dealii::CellId, std::vector< double > > *bQuadValues)
dftfe::uInt d_densityDofHandlerIndex
Definition dft.h:1489
const expConfiningPotential & getConfiningPotential() const
std::shared_ptr< dftfe::basis::FEBasisOperations< double, double, memorySpace > > getBasisOperationsElectroMemSpace()
~dftClass()
dftClass destructor
std::vector< std::vector< double > > atomLocationsFractional
Definition dft.h:1324
void normalizeRhoMagInInitialGuessQuadValues()
normalize input mag electron density to total magnetization for use in constraint magnetization case ...
std::vector< std::vector< double > > d_partialOccupancies
Definition dft.h:1663
dealii::DoFHandler< 3 > dofHandlerEigen
Definition dft.h:1486
chebyshevOrthogonalizedSubspaceIterationSolver * getSubspaceIterationSolverHost()
Get the Ptr to Chebyshev solver in host.
vselfBinsManager d_vselfBinsManager
vselfBinsManager object
Definition dft.h:1826
std::vector< distributedCPUVec< double > > d_vselfFieldGateauxDerStrainFDBins
Definition dft.h:1831
std::vector< double > d_gaussianConstantsForce
Definition dft.h:1343
dealii::AffineConstraints< double > d_constraintsForHelmholtzRhoNodal
Definition dft.h:1646
std::map< dftfe::uInt, std::map< dealii::CellId, std::vector< double > > > d_hessianRhoCoreAtoms
Definition dft.h:1818
dftfe::utils::MemoryStorage< double, dftfe::utils::MemorySpace::HOST > d_gradPhiResQuadValues
Definition dft.h:1724
const std::vector< std::vector< double > > & getEigenValues() const
Get reference to the eigen values.
void generateMPGrid()
const double d_pspCutOffTrunc
distance from the domain till which periodic images will be considered
Definition dft.h:1407
distributedCPUVec< double > d_phiTotRhoIn
Definition dft.h:1779
double d_minDist
Definition dft.h:1369
void init()
Does KSDFT problem pre-processing steps including mesh generation calls.
void readPSIRadialValues()
dftfe::utils::MemoryStorage< double, dftfe::utils::MemorySpace::HOST > d_densityTotalOutValuesLpspQuad
Definition dft.h:1737
std::deque< distributedCPUVec< double > > d_vSpin0containerVals
Definition dft.h:1754
std::deque< distributedCPUVec< double > > d_groundStateDensityHistory
Definition dft.h:1794
void initNoRemesh(const bool updateImagesAndKPointsAndVselfBins=true, const bool checkSmearedChargeWidthsForOverlap=true, const bool useSingleAtomSolutionOverride=false, const bool isMeshDeformed=false)
Does KSDFT problem pre-processing steps but without remeshing.
dftfe::uInt d_densityQuadratureId
Definition dft.h:1505
meshMovementGaussianClass d_gaussianMovePar
meshMovementGaussianClass object
Definition dft.h:1470
double getTotalChargeforRhoSplit()
dftfe::uInt d_sparsityPatternQuadratureId
Definition dft.h:1507
void calculateSmearedChargeWidths()
a
dftfe::uInt d_eigenDofHandlerIndex
Definition dft.h:1488
void calculateNearestAtomDistances()
a
void loadPSIFiles(dftfe::uInt Z, dftfe::uInt n, dftfe::uInt l, dftfe::uInt &flag)
Namespace which declares the input parameters and the functions to parse them from the input paramete...
Definition dftParameters.h:36
Overloads dealii's distribute and distribute_local_to_global functions associated with constraints cl...
Definition constraintMatrixInfo.h:43
Calculates dispersion correction to energy, force and stress.
Definition dftd.h:37
Manager class for ELPA and ScaLAPACK.
Definition elpaScalaManager.h:38
Definition expConfiningPotential.h:29
Definition kerkerSolverProblemWrapper.h:64
Definition BLASWrapper.h:35
Class to update triangulation under affine transformation.
Definition meshMovementAffineTransform.h:30
Class to move triangulation nodes using Gaussian functions attached to control points.
Definition meshMovementGaussian.h:30
Definition poissonSolverProblemWrapper.h:66
This class generates and stores adaptive finite element meshes for the real-space dft problem.
Definition triangulationManager.h:42
Definition MemoryStorage.h:33
Categorizes atoms into bins for efficient solution of nuclear electrostatic self-potential.
Definition vselfBinsManager.h:34
Definition FEBasisOperations.h:30
double number
Definition dftfeDataTypes.h:41
@ HOST
Definition MemorySpaceType.h:34
@ DEVICE
Definition MemorySpaceType.h:36
Definition pseudoPotentialToDftfeConverter.cc:34
std::uint32_t uInt
Definition TypeConfig.h:10
std::int32_t Int
Definition TypeConfig.h:11
dealii::LinearAlgebra::distributed::Vector< elem_type, dealii::MemorySpace::Host > distributedCPUVec
Definition headers.h:92