Commit 26449546 by Francois Gygi

merge --reintegrate of oncv branch into trunk

git-svn-id: http://qboxcode.org/svn/qb/trunk@1673 cba15fb0-1239-40c8-b417-11db7ca47a34
parent 40e559d5
......@@ -24,12 +24,12 @@
#include <iomanip>
#include <algorithm> // fill
#include <functional>
using namespace std;
////////////////////////////////////////////////////////////////////////////////
ChargeDensity::ChargeDensity(const Wavefunction& wf) : ctxt_(wf.context()),
wf_(wf), vcomm_(wf.sd(0,0)->basis().comm())
{
vbasis_ = new Basis(vcomm_, D3vector(0,0,0));
vbasis_->resize(wf.cell(),wf.refcell(),4.0*wf.ecut());
......@@ -89,6 +89,8 @@ wf_(wf), vcomm_(wf.sd(0,0)->basis().comm())
ft_[ikp] = new FourierTransform(wb,np0v,np1v,np2v);
}
}
// initialize core density ptr to null ptr
rhocore_r = 0;
}
////////////////////////////////////////////////////////////////////////////////
......@@ -171,8 +173,15 @@ void ChargeDensity::update_density(void)
tmap["charge_vft"].start();
vft_->forward(&rhotmp[0],&rhog[ispin][0]);
tmap["charge_vft"].stop();
// add core correction charge
if ( rhocore_r )
for ( int i = 0; i < rhor_size; i++ )
rhor[ispin][i] += rhocore_r[i];
}
}
////////////////////////////////////////////////////////////////////////////////
void ChargeDensity::update_rhor(void)
{
......@@ -191,9 +200,16 @@ void ChargeDensity::update_rhor(void)
const int rhor_size = rhor[ispin].size();
double *const prhor = &rhor[ispin][0];
if ( rhocore_r )
{
#pragma omp parallel for
for ( int i = 0; i < rhor_size; i++ )
prhor[i] = ( rhotmp[i].real() + rhocore_r[i] ) * omega_inv;
}
else
{
#pragma omp parallel for
for ( int i = 0; i < rhor_size; i++ )
prhor[i] = rhotmp[i].real() * omega_inv;
}
}
......
......@@ -51,7 +51,8 @@ class ChargeDensity
std::vector<std::vector<double> > rhor; // rhor[ispin][i]
std::vector<std::vector<std::complex<double> > > rhog; // rhog[ispin][ig]
// core density ptr. If non-zero, contains the real-space core density
double* rhocore_r;
void update_density(void);
void update_rhor(void);
......
......@@ -34,11 +34,12 @@
#include <iostream>
#include <iomanip>
#include <algorithm> // fill()
#include <algorithm> // fill(), copy()
#include <functional>
using namespace std;
////////////////////////////////////////////////////////////////////////////////
EnergyFunctional::EnergyFunctional( Sample& s, const ChargeDensity& cd)
EnergyFunctional::EnergyFunctional( Sample& s, ChargeDensity& cd)
: s_(s), cd_(cd)
{
const AtomSet& atoms = s_.atoms;
......@@ -121,6 +122,8 @@ EnergyFunctional::EnergyFunctional( Sample& s, const ChargeDensity& cd)
eself_ = 0.0;
// core_charge_: true if any species has a core charge density
core_charge_ = false;
for ( int is = 0; is < nsp_; is++ )
{
Species *s = atoms.species_list[is];
......@@ -133,8 +136,22 @@ EnergyFunctional::EnergyFunctional( Sample& s, const ChargeDensity& cd)
eself_ += na * s->eself();
na_[is] = na;
zv_[is] = s->zval();
zv_[is] = s->ztot();
rcps_[is] = s->rcps();
core_charge_ |= s->has_nlcc();
}
if ( core_charge_ )
{
vxc_g.resize(ngloc);
rhocore_sp_g.resize(nsp_);
for ( int is = 0; is < nsp_; is++ )
rhocore_sp_g[is].resize(ngloc);
rhocore_g.resize(ngloc);
rhocore_r.resize(vft->np012loc());
// set rhocore_r ptr in ChargeDensity
cd_.rhocore_r = &rhocore_r[0];
}
// FT for interpolation of wavefunctions on the fine grid
......@@ -229,11 +246,32 @@ void EnergyFunctional::update_vhxc(bool compute_stress)
// compute xc energy, update self-energy operator and potential
tmap["exc"].start();
for ( int ispin = 0; ispin < wf.nspin(); ispin++ )
fill(v_r[ispin].begin(),v_r[ispin].end(),0.0);
memset((void*)&v_r[ispin][0], 0, vft->np012loc()*sizeof(double));
xco->update(v_r, compute_stress);
exc_ = xco->exc();
if ( compute_stress )
xco->compute_stress(sigma_exc);
if ( core_charge_ )
{
// compute Fourier coefficients of Vxc
for ( int ispin = 0; ispin < wf.nspin(); ispin++ )
{
copy(v_r[ispin].begin(),v_r[ispin].end(),tmp_r.begin());
if ( ispin == 0 )
{
vft->forward(&tmp_r[0],&vxc_g[0]);
}
else
{
vft->forward(&tmp_r[0],&vtemp[0]);
for ( int ig = 0; ig < ngloc; ig++ )
vxc_g[ig] = 0.5 * ( vxc_g[ig] + vtemp[ig] );
}
}
}
tmap["exc"].stop();
// compute local potential energy:
......@@ -535,11 +573,22 @@ double EnergyFunctional::energy(bool compute_hpsi, Wavefunction& dwf,
{
const complex<double> sg = s[ig];
const complex<double> rg = rhogt[ig];
Species *s = s_.atoms.species_list[is];
const double g = vbasis_->g_ptr()[ig];
const double gi = vbasis_->gi_ptr()[ig];
// next line: keep only real part
const double temp = fac2 *
( rg.real() * sg.real() +
rg.imag() * sg.imag() )
// contribution of pseudocharge of ion
double temp = fac2 * (rg.real() * sg.real() + rg.imag() * sg.imag())
* rhops[is][ig] * g2i[ig];
if ( core_charge_ )
{
double rhoc_g, drhoc_g;
s->drhocoreg(g,rhoc_g,drhoc_g);
// next line: keep only real part
// contribution of core correction
temp -= (vxc_g[ig].real() * sg.real() + vxc_g[ig].imag() * sg.imag())
* drhoc_g * gi * omega_inv / fpi;
}
const double tgx = g_x[ig];
const double tgy = g_y[ig];
......@@ -670,6 +719,8 @@ double EnergyFunctional::energy(bool compute_hpsi, Wavefunction& dwf,
{
double tmp = fpi * rhops[is][ig] * g2i[ig];
vtemp[ig] = tmp * conj(rhogt[ig]) + vps[is][ig] * conj(rhoelg[ig]);
if ( core_charge_ )
vtemp[ig] += conj(vxc_g[ig]) * rhocore_sp_g[is][ig];
}
fion_nl.resize(3*na_[is]);
fion_nl = 0.0;
......@@ -863,6 +914,7 @@ double EnergyFunctional::energy(bool compute_hpsi, Wavefunction& dwf,
void EnergyFunctional::atoms_moved(void)
{
const AtomSet& atoms = s_.atoms;
const Wavefunction& wf = s_.wf;
int ngloc = vbasis_->localsize();
// fill tau0 with values in atom_list
......@@ -874,6 +926,29 @@ void EnergyFunctional::atoms_moved(void)
memset( (void*)&vion_local_g[0], 0, 2*ngloc*sizeof(double) );
memset( (void*)&dvion_local_g[0], 0, 2*ngloc*sizeof(double) );
memset( (void*)&rhopst[0], 0, 2*ngloc*sizeof(double) );
if ( core_charge_ )
{
// recalculate Fourier coefficients of the core charge
assert(rhocore_g.size()==ngloc);
memset( (void*)&rhocore_g[0], 0, 2*ngloc*sizeof(double) );
const double spin_fac = wf.cell().volume() / wf.nspin();
for ( int is = 0; is < atoms.nsp(); is++ )
{
complex<double> *s = &sf.sfac[is][0];
for ( int ig = 0; ig < ngloc; ig++ )
{
const complex<double> sg = s[ig];
// divide core charge by two if two spins
rhocore_g[ig] += spin_fac * sg * rhocore_sp_g[is][ig];
}
}
// recalculate real-space core charge
vft->backward(&rhocore_g[0],&tmp_r[0]);
for ( int i = 0; i < tmp_r.size(); i++ )
rhocore_r[i] = tmp_r[i].real();
}
for ( int is = 0; is < atoms.nsp(); is++ )
{
complex<double> *s = &sf.sfac[is][0];
......@@ -1007,13 +1082,18 @@ void EnergyFunctional::cell_moved(void)
{
Species *s = atoms.species_list[is];
const double * const g = vbasis_->g_ptr();
double v,dv;
double v,dv,rhoc;
for ( int ig = 0; ig < ngloc; ig++ )
{
rhops[is][ig] = s->rhopsg(g[ig]) * omega_inv;
s->dvlocg(g[ig],v,dv);
vps[is][ig] = v * omega_inv;
dvps[is][ig] = dv * omega_inv;
if ( core_charge_ )
{
s->rhocoreg(g[ig],rhoc);
rhocore_sp_g[is][ig] = rhoc * omega_inv;
}
}
}
......
......@@ -45,7 +45,7 @@ class EnergyFunctional
private:
Sample& s_;
const ChargeDensity& cd_;
ChargeDensity& cd_;
Basis* vbasis_;
FourierTransform *vft;
std::vector<FourierTransform*> ft;
......@@ -54,9 +54,10 @@ class EnergyFunctional
std::vector<std::vector<NonLocalPotential*> > nlp; // nlp[ispin][ikp]
std::vector<ConfinementPotential*> cfp; // cfp[ikp]
std::vector<std::vector<double> > vps, dvps, rhops;
std::vector<std::vector<double> > vps, dvps, rhops, rhocore_sp_g;
std::vector<double> rhocore_r;
std::vector<std::complex<double> > tmp_r, vion_local_g,
dvion_local_g, vlocal_g, rhopst, rhogt, rhoelg, vtemp;
dvion_local_g, vxc_g, vlocal_g, rhopst, rhogt, rhoelg, vtemp, rhocore_g;
std::vector<std::vector<double> > tau0, fion_esr;
std::vector<std::vector<double> > fext;
......@@ -68,6 +69,8 @@ class EnergyFunctional
std::valarray<double> sigma_ekin,sigma_econf,sigma_eps,sigma_ehart,sigma_exc,
sigma_enl, sigma_esr, sigma;
bool core_charge_;
public:
std::vector<std::vector<double> > v_r;
......@@ -99,7 +102,7 @@ class EnergyFunctional
void print(std::ostream& os) const;
EnergyFunctional(Sample& s, const ChargeDensity& cd);
EnergyFunctional(Sample& s, ChargeDensity& cd);
~EnergyFunctional();
};
std::ostream& operator << ( std::ostream& os, const EnergyFunctional& e );
......
......@@ -38,7 +38,7 @@ class NonLocalPotential
int nsp; // number of species
int nspnl; // number of non-local species
std::vector<int> lmax; // lmax[is]
std::vector<int> nop; // nop[is]
std::vector<int> lloc; // lloc[is]
std::vector<int> na; // na[is]
std::vector<int> npr; // npr[is]
......
......@@ -23,12 +23,29 @@
#include <string>
#include <vector>
#include <cmath>
#include <utility>
class Species
{
private:
int nlm_; // number of non-local projectors:
struct ProjectorData
{
int l, m, n;
};
enum PP_type
{
// norm conserving pseudopotential
NCPP,
// norm conserving semilocal pseudopotential
SLPP
};
PP_type type_; // identify type of the PP
int nlm_; // number of non-local projectors (including m)
int nop_; // number of non-local projectors (excluding m)
int ndft_;
std::vector<std::vector<double> > vps_spl_, vps_spl2_, phi_spl_, phi_spl2_;
......@@ -36,6 +53,11 @@ class Species
std::vector<std::vector<double> > vnlg_spl_, vnlg_spl2_;
std::vector<double> wsg_; // wsg_[l] Kleinman-Bylander weight
// 1/<phi|delta_V|phi>
std::vector<double> nlcc_spl_, nlcc_spl2_;
std::vector<double> nlccg_spl_, nlccg_spl2_;
// compensation charge
std::vector<std::vector<std::vector<double> > > q_spl_, q_spl2_;
std::vector<std::vector<std::vector<double> > > qg_spl_, qg_spl2_;
std::vector<double> rps_spl_; // radial linear mesh (same for all l)
......@@ -48,6 +70,8 @@ class Species
std::string description_; // description of the pseudopotential
int zval_; // valence charge
double zcore_; // core charge
double ztot_; // total charge
int lmax_; // largest angular momentum
int llocal_; // angular momentum taken as local
int nquad_; // number of semi-local quadrature points
......@@ -56,6 +80,39 @@ class Species
std::vector<std::vector<double> > vps_; // potentials for each l (input)
std::vector<std::vector<double> > phi_; // atomic wf for each l (input)
double rcps_; // cutoff radius of gaussian pseudocharge
std::vector<double> nlcc_; // non linear core correction
// for USPP
// indices: l - angular moment of projector
// n,m - differentiate projector with same l
// i,j - differentiate projector of all l
// max_i = sum_l max_n(l)
// map index of projector -> angular moment
std::vector<int> lmap_;
// local potential in radial representation
std::vector<double> vlocal_;
// projector in radial representation, indices proj_[l,n,r]
std::vector<std::vector<std::vector<double> > > proj_;
// compensation charge, indices q_[i,j,r]
std::vector<std::vector<std::vector<double> > > q_;
// matrix D in block diagonal storage, indices d_[l,n,m]
std::vector<std::vector<std::vector<double> > > d_;
// for norm conserving semilocal PP
// potential associated with a projector in radial representation,
// indices pot_[l,n,r]
std::vector<std::vector<std::vector<double> > > pot_;
// initialize a norm conserving PP
bool initialize_ncpp();
// initialize a semilocal potential
bool initialize_slpp();
// non linear core correction
void initialize_nlcc();
// helper function that extracts l, m and n from projector index
ProjectorData get_proj_data(int ipr);
public:
......@@ -67,6 +124,8 @@ class Species
const std::string& uri(void) const { return uri_; }
int atomic_number(void) const { return atomic_number_; }
int zval(void) const { return zval_; }
double zcore(void) const { return zcore_; }
double ztot(void) const { return ztot_; }
double mass(void) const { return mass_; }
int lmax(void) const { return lmax_; }
int llocal(void) const { return llocal_; }
......@@ -74,23 +133,35 @@ class Species
double rquad(void) const { return rquad_; }
double deltar(void) const { return deltar_; }
double rcps(void) const { return rcps_; }
bool has_nlcc(void) const { return nlcc_.size() > 0; }
bool has_dmatrix(void) const { return d_.size() > 0; }
// number of non-local projectors sum_(l!=llocal) (2*l+1)
int nlm(void) { return nlm_; }
bool non_local(void) { return lmax_ > 0; };
// number of non-local projectors w/o m degeneracy
int nop(void) { return nop_; }
// angular moment of projector with index iop
int l(int iop) { return lmap_[iop]; }
// extract D matrix
double dmatrix(int ipr, int jpr);
bool non_local(void) { return nop_ > 0; };
double eself(void)
{ return zval_ * zval_ / ( sqrt ( 2.0 * M_PI ) * rcps_ ); };
{ return ztot_ * ztot_ / ( sqrt ( 2.0 * M_PI ) * rcps_ ); };
void phi(int l, double r, double &val); // phi(l,r) in r space
void vpsr(int l, double r, double &v); // Vps(l,r) in r space
void dvpsr(int l, double r, double &v, double &dv); // Vps and dVps/dr
void vlocg(double q, double &v); // Vloc(g)
void dvlocg(double q, double &v, double &dv); // Vloc(g) and dVloc/dg
void vnlg(int l, double q, double &v); // Vnl(l,g)
void dvnlg(int l, double q, double &v, double &dv); // Vnl(l,g) and dVnl/dg
void vnlg(int iop, double q, double &v); // Vnl(iop,g)
void dvnlg(int iop, double q, double &v, double &dv);// Vnl(iop,g) and dVnl/dg
double rhopsg(double q); // pseudocharge in g space
double wsg(int l) { return wsg_[l]; };
void rhocoreg(double q, double &rho); // core correction in g space
// core correction and its derivative in g space
void drhocoreg(double q, double &rho, double &drho);
double wsg(int iop) { return wsg_[iop]; };
double rcut_loc(double epsilon); // radius beyond which potential is local
const std::vector<std::vector<double> >& vps(void) const { return vps_spl_; }
......@@ -99,6 +170,8 @@ class Species
bool initialize(double rcps);
void info(std::ostream& os);
void print(std::ostream& os, bool expanded_form);
void print_nlcc(std::ostream& os);
void print_radial_function(std::ostream& os, const std::vector<double>& rad_func);
friend class SpeciesReader;
friend class SpeciesHandler;
......
......@@ -26,17 +26,63 @@ using namespace xercesc;
using namespace std;
////////////////////////////////////////////////////////////////////////////////
SpeciesHandler::SpeciesHandler(Species& sp) : sp_(sp) {}
SpeciesHandler::SpeciesHandler(Species& sp) :
sp_(sp), d_ij_alloc(false) {}
////////////////////////////////////////////////////////////////////////////////
SpeciesHandler::~SpeciesHandler() {}
////////////////////////////////////////////////////////////////////////////////
// read attributes
void SpeciesHandler::read(const Attributes& attributes)
{
unsigned int len = attributes.getLength();
for ( unsigned int index = 0; index < len; index++ )
{
string attrname(XMLString::transcode(attributes.getLocalName(index)));
if ( attrname == "l" )
{
current_l = atoi(StrX(attributes.getValue(index)).localForm());
}
else if ( attrname == "size" )
{
current_size = atoi(StrX(attributes.getValue(index)).localForm());
}
else if ( attrname == "i" )
{
current_i = atoi(StrX(attributes.getValue(index)).localForm()) - 1;
}
else if ( attrname == "j" )
{
current_j = atoi(StrX(attributes.getValue(index)).localForm()) - 1;
}
}
}
////////////////////////////////////////////////////////////////////////////////
// allocate array for d_ij matrix
void SpeciesHandler::alloc_d_ij()
{
// allocate for every l
sp_.d_.resize(sp_.proj_.size());
for ( int l = 0; l < sp_.d_.size(); l++ )
{
const int size = sp_.proj_[l].size();
sp_.d_[l].resize(size);
for ( int i = 0; i < size; i++ )
{
sp_.d_[l][i].resize(size);
}
}
d_ij_alloc = true;
}
////////////////////////////////////////////////////////////////////////////////
void SpeciesHandler::startElement(const XMLCh* const uri,
const XMLCh* const localname, const XMLCh* const qname,
const Attributes& attributes)
{
// cout << " SpeciesHandler::startElement " << StrX(qname) << endl;
// cout << " SpeciesHandler::startElement " << StrX(qname) << endl;
string locname(XMLString::transcode(localname));
......@@ -44,10 +90,10 @@ void SpeciesHandler::startElement(const XMLCh* const uri,
{
// check for the case where the species is a link to another uri
unsigned int len = attributes.getLength();
for (unsigned int index = 0; index < len; index++)
for ( unsigned int index = 0; index < len; index++ )
{
string attrname(XMLString::transcode(attributes.getLocalName(index)));
if ( attrname == "name")
if ( attrname == "name" )
{
current_name = StrX(attributes.getValue(index)).localForm();
sp_.name_ = current_name;
......@@ -59,21 +105,29 @@ void SpeciesHandler::startElement(const XMLCh* const uri,
}
}
}
else if ( locname == "projector" )
else if ( locname == "norm_conserving_pseudopotential" )
{
unsigned int len = attributes.getLength();
for (unsigned int index = 0; index < len; index++)
sp_.type_ = Species::NCPP;
}
else if ( locname == "norm_conserving_semilocal_pseudopotential" )
{
string attrname(XMLString::transcode(attributes.getLocalName(index)));
if ( attrname == "l")
sp_.type_ = Species::SLPP;
}
else if ( locname == "core_density" )
{
current_l = atoi(StrX(attributes.getValue(index)).localForm());
read(attributes);
}
else if ( attrname == "size" )
else if ( locname == "local_potential" )
{
current_size = atoi(StrX(attributes.getValue(index)).localForm());
read(attributes);
}
else if ( locname == "projector" )
{
read(attributes);
}
else if ( locname == "d_ij" )
{
read(attributes);
}
}
......@@ -83,18 +137,17 @@ void SpeciesHandler::endElement(const XMLCh* const uri,
{
string locname(XMLString::transcode(localname));
istringstream stst(content);
// cout << " SpeciesHandler::endElement " << StrX(qname)
// << " content=" << string(content,0,20) << endl;
// cout << " SpeciesHandler::endElement " << StrX(qname) << " content="
// << string(content,0,20) << endl;
if ( locname == "description")
if ( locname == "description" )
{
// reject ambiguous case where both the href and the definition are given
if ( current_href != "" )
{
cout << " SpeciesHandler: ambiguous definition: uri="
<< StrX(uri) << endl
<< " using local definition (href: " << current_href << " ignored)"
<< endl;
cout << " SpeciesHandler: ambiguous definition: uri=" << StrX(uri)
<< endl << " using local definition (href: " << current_href
<< " ignored)" << endl;
}
sp_.description_ = content;
}
......@@ -136,21 +189,58 @@ void SpeciesHandler::endElement(const XMLCh* const uri,
}
else if ( locname == "radial_potential" )
{
if ( current_l+1 > sp_.vps_.size() )
assert(sp_.type_ == Species::NCPP);
if ( current_l + 1 > sp_.vps_.size() )
{
sp_.vps_.resize(current_l+1);
sp_.phi_.resize(current_l+1);
sp_.vps_.resize(current_l + 1);
sp_.phi_.resize(current_l + 1);
}
sp_.vps_[current_l].resize(current_size);
for ( int i = 0; i < current_size; i++ )
stst >> sp_.vps_[current_l][i];
}
else if ( locname == "local_potential" )
{
assert(sp_.type_ == Species::SLPP);
sp_.vlocal_.resize(current_size);
for ( int i = 0; i < current_size; i++ )
stst >> sp_.vlocal_[i];
}
else if ( locname == "radial_function" )
{
assert(sp_.type_ == Species::NCPP);
sp_.phi_[current_l].resize(current_size);
for ( int i = 0; i < current_size; i++ )
stst >> sp_.phi_[current_l][i];
}
else if ( locname == "core_density" )
{
// read charge for nonlinear core correction
sp_.nlcc_.resize(current_size);
for ( int i = 0; i < current_size; i++ )
stst >> sp_.nlcc_[i];
}
else if ( locname == "projector" and sp_.type_ == Species::SLPP )
{
// read one of the projector with this angular momentum
// resize vector if necessary
if ( current_l >= sp_.proj_.size() ) sp_.proj_.resize(current_l + 1);
if ( current_i >= sp_.proj_[current_l].size() )
sp_.proj_[current_l].resize(current_i + 1);
sp_.proj_[current_l][current_i].resize(current_size);
// read the projector
for ( int i = 0; i < current_size; i++ )
stst >> sp_.proj_[current_l][current_i][i];
}
else if ( locname == "d_ij" )
{
assert(sp_.type_ == Species::SLPP);
if ( not d_ij_alloc ) alloc_d_ij();
assert(current_l < sp_.d_.size());
assert(current_i < sp_.d_[current_l].size());
assert(current_j < sp_.d_[current_l][current_i].size());
stst >> sp_.d_[current_l][current_i][current_j];
}
}
////////////////////////////////////////////////////////////////////////////////
......
......@@ -28,9 +28,16 @@ class SpeciesHandler : public StructureHandler
private:
Species& sp_;
int current_l, current_size;
int current_l, current_size, current_i, current_j;
std::string current_name, current_href;
double current_interval;
bool d_ij_alloc;
// read attributes
void read(const Attributes& attributes);
// allocate array for d_ij matrix
void alloc_d_ij();
public:
......
--------------------------------------------------------------------------------
Planned merging of oncv branch into trunk. Modifications.
Species.h: Must check use of namespace to define proj_data object.
Species.h: ztot_ zcore_ are double. Related to the definition given
in the species that may not be an integer charge?
Species.C: use of "or" keyword in test is a C++-11 feature. Revert to || for
consistency with the rest of the code.
EnergyFunctional.h: ChargeDensity cd_ not const. Due to core charge update?
ChargeDensity.C::update_rhor: use of transform algorithm. Add and subtract
core charge can be replaced by use of temporary array, simpler loop.
sample.xsd: remove definition of doubleListType which appears in species.xsd
merging trunk into oncv branch before reintegrate.
--------------------------------------------------------------------------------
rel1_60_9
r1656 ExchangeOperator.C: changed default bisection tolerance to 1.0
......
......@@ -19,5 +19,5 @@
#include "release.h"
std::string release(void)
{
return std::string("wrk");
return std::string("oncv");
}
......@@ -148,10 +148,6 @@
</simpleContent>
</complexType>
<simpleType name="doubleListType">
<list itemType="double"/>
</simpleType>
<complexType name="density_matrixType">
<simpleContent>
<extension base="fpmd:doubleListType">
......
......@@ -4,10 +4,10 @@
xmlns:fpmd="http://www.quantum-simulation.org/ns/fpmd/fpmd-1.0">
<annotation>
<documentation> $Id: species.xsd,v 1.6 2008-03-07 20:05:31 fgygi Exp $
<documentation>
http://www.quantum-simulation.org
FPMD atomic species XML Schema specification.
Copyright (c) 2006-2008 The Regents of the University of California.
Copyright (c) 2006-2014 The Regents of the University of California.
</documentation>
</annotation>
......@@ -19,8 +19,12 @@
<element name="symbol" type="NMTOKEN"/>
<element name="atomic_number" type="nonNegativeInteger"/>
<element name="mass" type="fpmd:positiveDouble"/>
<element name="norm_conserving_pseudopotential" minOccurs="0" maxOccurs="1"
<choice>
<element name="norm_conserving_pseudopotential"
type="fpmd:norm_conserving_pseudopotentialType"/>
<element name="norm_conserving_semilocal_pseudopotential"
type="fpmd:norm_conserving_semilocal_pseudopotentialType"/>