QMCPACK
BsplineReader Struct Referenceabstract

Each SplineC2X needs a reader derived from BsplineReader. More...

+ Inheritance diagram for BsplineReader:
+ Collaboration diagram for BsplineReader:

Public Member Functions

 BsplineReader (EinsplineSetBuilder *e)
 
virtual ~BsplineReader ()
 
std::string getSplineDumpFileName (const BandInfoGroup &bandgroup) const
 
template<typename GT , typename BCT >
bool set_grid (const TinyVector< int, 3 > &halfg, GT *xyz_grid, BCT *xyz_bc) const
 read gvectors and set the mesh, and prepare for einspline More...
 
template<typename SPE >
void check_twists (SPE &bspline, const BandInfoGroup &bandgroup) const
 initialize twist-related data for N orbitals More...
 
std::string psi_g_path (int ti, int spin, int ib) const
 return the path name in hdf5 More...
 
std::string psi_r_path (int ti, int spin, int ib) const
 return the path name in hdf5 More...
 
virtual std::unique_ptr< SPOSetcreate_spline_set (const std::string &my_name, int spin, const BandInfoGroup &bandgroup)=0
 create the actual spline sets More...
 
void setCommon (xmlNodePtr cur)
 setting common parameters More...
 
std::unique_ptr< SPOSetcreate_spline_set (int spin, xmlNodePtr cur, SPOSetInputInfo &input_info)
 create the spline after one of the kind is created More...
 
std::unique_ptr< SPOSetcreate_spline_set (int spin, xmlNodePtr cur)
 create the spline set More...
 
void setCheckNorm (bool new_checknorm)
 Set the checkNorm variable. More...
 
void setRotate (bool new_rotate)
 Set the orbital rotation flag. More...
 
void initialize_spo2band (int spin, const std::vector< BandInfo > &bigspace, SPOSetInfo &sposet, std::vector< int > &band2spo)
 build index tables to map a state to band with k-point folidng More...
 

Public Attributes

EinsplineSetBuildermybuilder
 pointer to the EinsplineSetBuilder More...
 
CommunicatemyComm
 communicator More...
 
bool checkNorm
 mesh size check the norm of orbitals More...
 
bool saveSplineCoefs
 save spline coefficients to storage More...
 
bool rotate
 apply orbital rotations More...
 
std::vector< std::vector< int > > spo2band
 map from spo index to band index More...
 

Detailed Description

Each SplineC2X needs a reader derived from BsplineReader.

This base class handles common chores

  • check_twists : read gvectors, set twists for folded bands if needed, and set the phase for the special K
  • set_grid : create the basic grid and boundary conditions for einspline Note that template is abused but it works.

Definition at line 39 of file BsplineReader.h.

Constructor & Destructor Documentation

◆ BsplineReader()

Definition at line 30 of file BsplineReader.cpp.

References MPIObjectBase::getCommunicator(), BsplineReader::mybuilder, and BsplineReader::myComm.

31  : mybuilder(e), checkNorm(true), saveSplineCoefs(false), rotate(true)
32 {
34 }
Communicate * myComm
communicator
Definition: BsplineReader.h:44
Communicate * getCommunicator() const
return myComm
Definition: MPIObjectBase.h:41
EinsplineSetBuilder * mybuilder
pointer to the EinsplineSetBuilder
Definition: BsplineReader.h:42
bool rotate
apply orbital rotations
Definition: BsplineReader.h:51
bool checkNorm
mesh size check the norm of orbitals
Definition: BsplineReader.h:47
bool saveSplineCoefs
save spline coefficients to storage
Definition: BsplineReader.h:49

◆ ~BsplineReader()

~BsplineReader ( )
virtualdefault

Member Function Documentation

◆ check_twists()

void check_twists ( SPE &  bspline,
const BandInfoGroup bandgroup 
) const
inline

initialize twist-related data for N orbitals

Definition at line 101 of file BsplineReader.h.

References qmcplusplus::abs(), qmcplusplus::app_log(), bspline(), qmcplusplus::dot(), qmcplusplus::Units::charge::e, BandInfoGroup::getFirstSPO(), BandInfoGroup::getLastSPO(), ParticleSet::getLattice(), BandInfoGroup::getNumDistinctOrbitals(), BandInfoGroup::getNumSPOs(), CrystalLattice< T, D >::k_cart(), BandInfoGroup::myBands, BsplineReader::mybuilder, qmcplusplus::Units::force::N, EinsplineSetBuilder::PrimCell, EinsplineSetBuilder::primcell_kpoints, EinsplineSetBuilder::TargetPtcl, qmcplusplus::transpose(), and BandInfoGroup::TwistIndex.

102  {
103  //init(orbitalSet,bspline);
104  bspline.PrimLattice = mybuilder->PrimCell;
105  bspline.GGt = dot(transpose(bspline.PrimLattice.G), bspline.PrimLattice.G);
106 
107  int N = bandgroup.getNumDistinctOrbitals();
108  int numOrbs = bandgroup.getNumSPOs();
109 
110  bspline.setOrbitalSetSize(numOrbs);
111  bspline.resizeStorage(N, N);
112 
113  bspline.first_spo = bandgroup.getFirstSPO();
114  bspline.last_spo = bandgroup.getLastSPO();
115 
116  int num = 0;
117  const std::vector<BandInfo>& cur_bands = bandgroup.myBands;
118  for (int iorb = 0; iorb < N; iorb++)
119  {
120  int ti = cur_bands[iorb].TwistIndex;
121  bspline.kPoints[iorb] = mybuilder->PrimCell.k_cart(-mybuilder->primcell_kpoints[ti]);
122  bspline.MakeTwoCopies[iorb] = (num < (numOrbs - 1)) && cur_bands[iorb].MakeTwoCopies;
123  num += bspline.MakeTwoCopies[iorb] ? 2 : 1;
124  }
125 
126  app_log() << "NumDistinctOrbitals " << N << " numOrbs = " << numOrbs << std::endl;
127 
128  bspline.HalfG = 0;
129  TinyVector<int, 3> bconds = mybuilder->TargetPtcl.getLattice().BoxBConds;
130  if (!bspline.isComplex())
131  {
132  //no k-point folding, single special k point (G, L ...)
133  TinyVector<double, 3> twist0 = mybuilder->primcell_kpoints[bandgroup.TwistIndex];
134  for (int i = 0; i < 3; i++)
135  if (bconds[i] && ((std::abs(std::abs(twist0[i]) - 0.5) < 1.0e-8)))
136  bspline.HalfG[i] = 1;
137  else
138  bspline.HalfG[i] = 0;
139  app_log() << " TwistIndex = " << cur_bands[0].TwistIndex << " TwistAngle " << twist0 << std::endl;
140  app_log() << " HalfG = " << bspline.HalfG << std::endl;
141  }
142  app_log().flush();
143  }
double bspline(double x, const std::vector< double > &t, const std::vector< double > &c, int k)
MakeReturn< UnaryNode< FnFabs, typename CreateLeaf< Vector< T1, C1 > >::Leaf_t > >::Expression_t abs(const Vector< T1, C1 > &l)
std::ostream & app_log()
Definition: OutputManager.h:65
ParticleSet & TargetPtcl
quantum particle set
SingleParticlePos k_cart(const SingleParticlePos &kin) const
conversion of a reciprocal-vector
AntiSymTensor< T, D > transpose(const AntiSymTensor< T, D > &rhs)
const auto & getLattice() const
Definition: ParticleSet.h:251
std::vector< TinyVector< double, OHMMS_DIM > > primcell_kpoints
Tensor< typename BinaryReturn< T1, T2, OpMultiply >::Type_t, D > dot(const AntiSymTensor< T1, D > &lhs, const AntiSymTensor< T2, D > &rhs)
EinsplineSetBuilder * mybuilder
pointer to the EinsplineSetBuilder
Definition: BsplineReader.h:42

◆ create_spline_set() [1/3]

virtual std::unique_ptr<SPOSet> create_spline_set ( const std::string &  my_name,
int  spin,
const BandInfoGroup bandgroup 
)
pure virtual

◆ create_spline_set() [2/3]

std::unique_ptr< SPOSet > create_spline_set ( int  spin,
xmlNodePtr  cur,
SPOSetInputInfo input_info 
)

create the spline after one of the kind is created

Definition at line 123 of file BsplineReader.cpp.

References OhmmsAttributeSet::add(), SPOSetBuilder::clear_states(), BsplineReader::create_spline_set(), EinsplineSetBuilder::FullBands, MPIObjectBase::getName(), BandInfoGroup::GroupID, BsplineReader::initialize_spo2band(), qmcplusplus::make_bandgroup_name(), SPOSetInputInfo::max_index(), SPOSetInputInfo::min_index(), BsplineReader::mybuilder, BandInfoGroup::myName, OhmmsAttributeSet::put(), BandInfoGroup::selectBands(), BsplineReader::spo2band, SPOSetBuilder::states, EinsplineSetBuilder::TileMatrix, EinsplineSetBuilder::twist_num_, and BandInfoGroup::TwistIndex.

124 {
125  std::string spo_object_name;
127  a.add(spo_object_name, "name");
128  a.add(spo_object_name, "id");
129  a.put(cur);
130 
131  if (spo2band.empty())
132  spo2band.resize(mybuilder->states.size());
133 
134  std::vector<BandInfo>& fullband = (*(mybuilder->FullBands[spin]));
135 
136  if (spo2band[spin].empty())
137  {
138  spo2band[spin].reserve(fullband.size());
139  if (!mybuilder->states[spin])
140  mybuilder->states[spin] = std::make_unique<SPOSetInfo>();
141  mybuilder->clear_states(spin);
142  initialize_spo2band(spin, fullband, *mybuilder->states[spin], spo2band[spin]);
143  }
144 
145  BandInfoGroup vals;
146  vals.TwistIndex = fullband[0].TwistIndex;
147  vals.GroupID = 0;
149  input_info.min_index(), input_info.max_index());
150  vals.selectBands(fullband, spo2band[spin][input_info.min_index()], input_info.max_index() - input_info.min_index(),
151  false);
152 
153  return create_spline_set(spo_object_name, spin, vals);
154 }
Tensor< int, OHMMS_DIM > TileMatrix
bool put(xmlNodePtr cur)
assign attributes to the set
Definition: AttributeSet.h:55
std::vector< std::unique_ptr< SPOSetInfo > > states
state info of all possible states available in the basis
Definition: SPOSetBuilder.h:57
std::vector< std::unique_ptr< std::vector< BandInfo > > > FullBands
Helper vector for sorting bands.
std::string make_bandgroup_name(const std::string &root, int spin, int twist, const Tensor< int, 3 > &tilematrix, int first, int last)
std::vector< std::vector< int > > spo2band
map from spo index to band index
Definition: BsplineReader.h:53
void clear_states(int index=0)
clear state information
Definition: SPOSetBuilder.h:69
class to handle a set of attributes of an xmlNode
Definition: AttributeSet.h:24
virtual std::unique_ptr< SPOSet > create_spline_set(const std::string &my_name, int spin, const BandInfoGroup &bandgroup)=0
create the actual spline sets
void initialize_spo2band(int spin, const std::vector< BandInfo > &bigspace, SPOSetInfo &sposet, std::vector< int > &band2spo)
build index tables to map a state to band with k-point folidng
const std::string & getName() const
return the name
Definition: MPIObjectBase.h:54
EinsplineSetBuilder * mybuilder
pointer to the EinsplineSetBuilder
Definition: BsplineReader.h:42
void add(PDT &aparam, const std::string &aname, std::vector< PDT > candidate_values={}, TagStatus status=TagStatus::OPTIONAL)
add a new attribute
Definition: AttributeSet.h:42

◆ create_spline_set() [3/3]

std::unique_ptr< SPOSet > create_spline_set ( int  spin,
xmlNodePtr  cur 
)

create the spline set

Definition at line 87 of file BsplineReader.cpp.

References OhmmsAttributeSet::add(), APP_ABORT_TRACE, SPOSetBuilder::clear_states(), BsplineReader::create_spline_set(), EinsplineSetBuilder::FullBands, MPIObjectBase::getName(), BandInfoGroup::GroupID, BsplineReader::initialize_spo2band(), qmcplusplus::make_bandgroup_name(), BsplineReader::mybuilder, BandInfoGroup::myName, qmcplusplus::Units::time::ns, OhmmsAttributeSet::put(), BandInfoGroup::selectBands(), BsplineReader::spo2band, SPOSetBuilder::states, EinsplineSetBuilder::TileMatrix, EinsplineSetBuilder::twist_num_, and BandInfoGroup::TwistIndex.

88 {
89  int ns(0);
90  std::string spo_object_name;
92  a.add(ns, "size");
93  a.add(spo_object_name, "name");
94  a.add(spo_object_name, "id");
95  a.put(cur);
96 
97  if (ns == 0)
98  APP_ABORT_TRACE(__FILE__, __LINE__, "parameter/@size missing");
99 
100  if (spo2band.empty())
101  spo2band.resize(mybuilder->states.size());
102 
103  std::vector<BandInfo>& fullband = (*(mybuilder->FullBands[spin]));
104 
105  if (spo2band[spin].empty())
106  {
107  spo2band[spin].reserve(fullband.size());
108  if (!mybuilder->states[spin])
109  mybuilder->states[spin] = std::make_unique<SPOSetInfo>();
110  mybuilder->clear_states(spin);
111  initialize_spo2band(spin, fullband, *mybuilder->states[spin], spo2band[spin]);
112  }
113 
114  BandInfoGroup vals;
115  vals.TwistIndex = fullband[0].TwistIndex;
116  vals.GroupID = 0;
118  vals.selectBands(fullband, 0, ns, false);
119 
120  return create_spline_set(spo_object_name, spin, vals);
121 }
Tensor< int, OHMMS_DIM > TileMatrix
bool put(xmlNodePtr cur)
assign attributes to the set
Definition: AttributeSet.h:55
std::vector< std::unique_ptr< SPOSetInfo > > states
state info of all possible states available in the basis
Definition: SPOSetBuilder.h:57
std::vector< std::unique_ptr< std::vector< BandInfo > > > FullBands
Helper vector for sorting bands.
std::string make_bandgroup_name(const std::string &root, int spin, int twist, const Tensor< int, 3 > &tilematrix, int first, int last)
std::vector< std::vector< int > > spo2band
map from spo index to band index
Definition: BsplineReader.h:53
void clear_states(int index=0)
clear state information
Definition: SPOSetBuilder.h:69
class to handle a set of attributes of an xmlNode
Definition: AttributeSet.h:24
virtual std::unique_ptr< SPOSet > create_spline_set(const std::string &my_name, int spin, const BandInfoGroup &bandgroup)=0
create the actual spline sets
void initialize_spo2band(int spin, const std::vector< BandInfo > &bigspace, SPOSetInfo &sposet, std::vector< int > &band2spo)
build index tables to map a state to band with k-point folidng
const std::string & getName() const
return the name
Definition: MPIObjectBase.h:54
EinsplineSetBuilder * mybuilder
pointer to the EinsplineSetBuilder
Definition: BsplineReader.h:42
#define APP_ABORT_TRACE(f, l, msg)
Definition: AppAbort.h:34
void add(PDT &aparam, const std::string &aname, std::vector< PDT > candidate_values={}, TagStatus status=TagStatus::OPTIONAL)
add a new attribute
Definition: AttributeSet.h:42

◆ getSplineDumpFileName()

std::string getSplineDumpFileName ( const BandInfoGroup bandgroup) const
inline

Definition at line 59 of file BsplineReader.h.

References EinsplineSetBuilder::MeshSize, BsplineReader::mybuilder, and BandInfoGroup::myName.

60  {
61  auto& MeshSize = mybuilder->MeshSize;
62  std::ostringstream oo;
63  oo << bandgroup.myName << ".g" << MeshSize[0] << "x" << MeshSize[1] << "x" << MeshSize[2] << ".h5";
64  return oo.str();
65  }
EinsplineSetBuilder * mybuilder
pointer to the EinsplineSetBuilder
Definition: BsplineReader.h:42

◆ initialize_spo2band()

void initialize_spo2band ( int  spin,
const std::vector< BandInfo > &  bigspace,
SPOSetInfo sposet,
std::vector< int > &  spo2band 
)

build index tables to map a state to band with k-point folidng

Parameters
bigspacefull BandInfo constructed by EinsplineSetBuilder
sposetSPOSetInfo owned by someone, most likely EinsplinseSetBuilder
spo2bandspo2band[i] is the index in bigspace

At gamma or arbitrary kpoints with complex wavefunctions, spo2band[i]==i

Definition at line 163 of file BsplineReader.cpp.

References SPOSetInfo::add(), qmcplusplus::comm, qmcplusplus::Units::charge::e, Communicate::getGroupID(), MPIObjectBase::getName(), CrystalLattice< T, D >::k_cart(), qmcplusplus::make_bandinfo_filename(), BsplineReader::mybuilder, BsplineReader::myComm, qmcplusplus::Units::time::ns, EinsplineSetBuilder::PrimCell, EinsplineSetBuilder::primcell_kpoints, Communicate::rank(), qmcplusplus::Units::time::s, BsplineReader::spo2band, EinsplineSetBuilder::TileMatrix, and EinsplineSetBuilder::twist_num_.

Referenced by BsplineReader::create_spline_set().

167 {
168  spo2band.reserve(bigspace.size());
169  int ns = 0;
170  for (int i = 0; i < bigspace.size(); ++i)
171  {
172  spo2band.push_back(i);
173  SPOInfo a(ns, bigspace[i].Energy);
174  sposet.add(a);
175  ns++;
176  if (bigspace[i].MakeTwoCopies)
177  {
178  spo2band.push_back(i);
179  SPOInfo b(ns, bigspace[i].Energy);
180  sposet.add(b);
181  ns++;
182  }
183  }
184 
185  //write to a file
186  const Communicate* comm = myComm;
187  if (comm->rank())
188  return;
189 
190  std::filesystem::path aname = make_bandinfo_filename(mybuilder->getName(), spin, mybuilder->twist_num_,
192  aname += ".bandinfo.dat";
193 
194  std::ofstream o(aname.c_str());
195  std::array<char, 1024> s;
196  ns = 0;
197  using PosType = QMCTraits::PosType;
198  o << "# Band State TwistIndex BandIndex Energy Kx Ky Kz K1 K2 K3 KmK "
199  << std::endl;
200  for (int i = 0; i < bigspace.size(); ++i)
201  {
202  int ti = bigspace[i].TwistIndex;
203  int bi = bigspace[i].BandIndex;
204  double e = bigspace[i].Energy;
205  int nd = (bigspace[i].MakeTwoCopies) ? 2 : 1;
207  int s_size = std::snprintf(s.data(), s.size(), "%8d %8d %8d %8d %12.6f %7.4f %7.4f %7.4f %7.4f %7.4f %7.4f %6d\n",
208  i, ns, ti, bi, e, k[0], k[1], k[2], mybuilder->primcell_kpoints[ti][0],
210  if (s_size < 0)
211  throw std::runtime_error("Error generating bandinfo");
212  o << s.data();
213  ns += nd;
214  }
215 }
Tensor< int, OHMMS_DIM > TileMatrix
int rank() const
return the rank
Definition: Communicate.h:116
Communicate * myComm
communicator
Definition: BsplineReader.h:44
std::string make_bandinfo_filename(const std::string &root, int spin, int twist, const Tensor< int, 3 > &tilematrix, int gid)
QMCTraits::PosType PosType
Wrapping information on parallelism.
Definition: Communicate.h:68
SingleParticlePos k_cart(const SingleParticlePos &kin) const
conversion of a reciprocal-vector
int getGroupID() const
return the group id
Definition: Communicate.h:121
std::vector< std::vector< int > > spo2band
map from spo index to band index
Definition: BsplineReader.h:53
QTBase::PosType PosType
Definition: Configuration.h:61
const std::string & getName() const
return the name
Definition: MPIObjectBase.h:54
std::vector< TinyVector< double, OHMMS_DIM > > primcell_kpoints
EinsplineSetBuilder * mybuilder
pointer to the EinsplineSetBuilder
Definition: BsplineReader.h:42

◆ psi_g_path()

std::string psi_g_path ( int  ti,
int  spin,
int  ib 
) const
inline

return the path name in hdf5

Parameters
titwist index
spinspin index
ibband index

Definition at line 150 of file BsplineReader.h.

151  {
152  std::ostringstream path;
153  path << "/electrons/kpoint_" << ti << "/spin_" << spin << "/state_" << ib << "/psi_g";
154  return path.str();
155  }

◆ psi_r_path()

std::string psi_r_path ( int  ti,
int  spin,
int  ib 
) const
inline

return the path name in hdf5

Parameters
titwist index
spinspin index
ibband index

Definition at line 162 of file BsplineReader.h.

163  {
164  std::ostringstream path;
165  path << "/electrons/kpoint_" << ti << "/spin_" << spin << "/state_" << ib << "/psi_r";
166  return path.str();
167  }

◆ set_grid()

bool set_grid ( const TinyVector< int, 3 > &  halfg,
GT *  xyz_grid,
BCT *  xyz_bc 
) const
inline

read gvectors and set the mesh, and prepare for einspline

Definition at line 70 of file BsplineReader.h.

References EinsplineSetBuilder::MeshSize, BsplineReader::mybuilder, and EinsplineSetBuilder::ReadGvectors_ESHDF().

71  {
72  //This sets MeshSize from the input file
73  bool havePsig = mybuilder->ReadGvectors_ESHDF();
74 
75  for (int j = 0; j < 3; ++j)
76  {
77  xyz_grid[j].start = 0.0;
78  xyz_grid[j].end = 1.0;
79  xyz_grid[j].num = mybuilder->MeshSize[j];
80 
81  if (halfg[j])
82  {
83  xyz_bc[j].lCode = ANTIPERIODIC;
84  xyz_bc[j].rCode = ANTIPERIODIC;
85  }
86  else
87  {
88  xyz_bc[j].lCode = PERIODIC;
89  xyz_bc[j].rCode = PERIODIC;
90  }
91 
92  xyz_bc[j].lVal = 0.0;
93  xyz_bc[j].rVal = 0.0;
94  }
95  return havePsig;
96  }
bool ReadGvectors_ESHDF()
read gvectors for each twist
EinsplineSetBuilder * mybuilder
pointer to the EinsplineSetBuilder
Definition: BsplineReader.h:42

◆ setCheckNorm()

void setCheckNorm ( bool  new_checknorm)
inline

Set the checkNorm variable.

Definition at line 186 of file BsplineReader.h.

References BsplineReader::checkNorm.

186 { checkNorm = new_checknorm; };
bool checkNorm
mesh size check the norm of orbitals
Definition: BsplineReader.h:47

◆ setCommon()

void setCommon ( xmlNodePtr  cur)

setting common parameters

Definition at line 68 of file BsplineReader.cpp.

References OhmmsAttributeSet::add(), qmcplusplus::app_log(), BsplineReader::checkNorm, OhmmsAttributeSet::put(), and BsplineReader::saveSplineCoefs.

69 {
70  // check orbital normalization by default
71  std::string checkOrbNorm("yes");
72  std::string saveCoefs("no");
74  a.add(checkOrbNorm, "check_orb_norm");
75  a.add(saveCoefs, "save_coefs");
76  a.put(cur);
77 
78  // allow user to turn off norm check with a warning
79  if (checkOrbNorm == "no")
80  {
81  app_log() << "WARNING: disable orbital normalization check!" << std::endl;
82  checkNorm = false;
83  }
84  saveSplineCoefs = saveCoefs == "yes";
85 }
std::ostream & app_log()
Definition: OutputManager.h:65
bool put(xmlNodePtr cur)
assign attributes to the set
Definition: AttributeSet.h:55
class to handle a set of attributes of an xmlNode
Definition: AttributeSet.h:24
void add(PDT &aparam, const std::string &aname, std::vector< PDT > candidate_values={}, TagStatus status=TagStatus::OPTIONAL)
add a new attribute
Definition: AttributeSet.h:42
bool checkNorm
mesh size check the norm of orbitals
Definition: BsplineReader.h:47
bool saveSplineCoefs
save spline coefficients to storage
Definition: BsplineReader.h:49

◆ setRotate()

void setRotate ( bool  new_rotate)
inline

Set the orbital rotation flag.

Rotations are applied to balance the real/imaginary components.

Definition at line 189 of file BsplineReader.h.

References BsplineReader::rotate.

189 { rotate = new_rotate; };
bool rotate
apply orbital rotations
Definition: BsplineReader.h:51

Member Data Documentation

◆ checkNorm

bool checkNorm

mesh size check the norm of orbitals

Definition at line 47 of file BsplineReader.h.

Referenced by BsplineReader::setCheckNorm(), and BsplineReader::setCommon().

◆ mybuilder

◆ myComm

Communicate* myComm

communicator

Definition at line 44 of file BsplineReader.h.

Referenced by BsplineReader::BsplineReader(), and BsplineReader::initialize_spo2band().

◆ rotate

bool rotate

apply orbital rotations

Definition at line 51 of file BsplineReader.h.

Referenced by BsplineReader::setRotate().

◆ saveSplineCoefs

bool saveSplineCoefs

save spline coefficients to storage

Definition at line 49 of file BsplineReader.h.

Referenced by BsplineReader::setCommon().

◆ spo2band

std::vector<std::vector<int> > spo2band

map from spo index to band index

Definition at line 53 of file BsplineReader.h.

Referenced by BsplineReader::create_spline_set(), and BsplineReader::initialize_spo2band().


The documentation for this struct was generated from the following files: