FOX/ObjCryst++  2022
ScatteringPower.cpp
1 /* ObjCryst++ Object-Oriented Crystallographic Library
2  (c) 2000-2002 Vincent Favre-Nicolin vincefn@users.sourceforge.net
3  2000-2001 University of Geneva (Switzerland)
4 
5  This program is free software; you can redistribute it and/or modify
6  it under the terms of the GNU General Public License as published by
7  the Free Software Foundation; either version 2 of the License, or
8  (at your option) any later version.
9 
10  This program is distributed in the hope that it will be useful,
11  but WITHOUT ANY WARRANTY; without even the implied warranty of
12  MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
13  GNU General Public License for more details.
14 
15  You should have received a copy of the GNU General Public License
16  along with this program; if not, write to the Free Software
17  Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
18 */
19 #include <cmath>
20 #include <typeinfo>
21 #include <iomanip>
22 #include <fstream>
23 #include <cstring>
24 #include <stdio.h> //for sprintf()
25 
26 
27 #include "cctbx/eltbx/xray_scattering.h"
28 #include "cctbx/eltbx/tiny_pse.h"
29 #include "cctbx/eltbx/icsd_radii.h"
30 #include "cctbx/eltbx/covalent_radii.h"
31 #include "cctbx/eltbx/henke.h"
32 #include "cctbx/eltbx/neutron.h"
33 
34 #include "ObjCryst/ObjCryst/ScatteringPower.h"
35 #include "ObjCryst/Quirks/VFNStreamFormat.h"
36 #include "ObjCryst/Quirks/VFNDebug.h"
37 #include "ObjCryst/ObjCryst/Colours.h"
38 
39 #ifdef __WX__CRYST__
40  #include "ObjCryst/wxCryst/wxScatteringPower.h"
41 #endif
42 
43 namespace ObjCryst
44 {
45 
46 const RefParType *gpRefParTypeScattPow=0;
47 const RefParType *gpRefParTypeScattPowResonant=0;
48 const RefParType *gpRefParTypeScattPowTemperature=0;
49 const RefParType *gpRefParTypeScattPowTemperatureIso=0;
50 const RefParType *gpRefParTypeScattPowTemperatureAniso=0;
51 long NiftyStaticGlobalObjectsInitializer_ScatteringPower::mCount=0;
52 //######################################################################
53 //
54 // Bij to Betaij conversion
55 //
56 //######################################################################
57 CrystMatrix_REAL Bij2Betaij(const CrystVector_REAL &Bij, const UnitCell &cell)
58 {
59  // Willis & Pryor,p 101: (betaij) = 2*pi^2 * transpose(cell.Bmatrix) * (Bij) * cell.Bmatrix
60  // :TODO: this needs to be checked before being used
61  const REAL B11=Bij(0);
62  const REAL B22=Bij(1);
63  const REAL B33=Bij(2);
64  const REAL B12=Bij(3);
65  const REAL B13=Bij(4);
66  const REAL B23=Bij(5);
67  CrystMatrix_REAL B(3,3);
68  B(0,0)=B11;
69  B(0,1)=B12;
70  B(0,1)=B13;
71  B(1,0)=B12;
72  B(1,1)=B22;
73  B(1,1)=B23;
74  B(2,0)=B13;
75  B(2,1)=B23;
76  B(2,1)=B33;
77  CrystMatrix_REAL b(3,3);
78  b=cell.GetBMatrix().transpose().Mult(B.Mult(cell.GetBMatrix()));
79  b*=2*M_PI*M_PI;
80  return b;
81 }
82 
83 //######################################################################
84 //
85 // SCATTERING POWER
86 //
87 //######################################################################
88 ObjRegistry<ScatteringPower> gScatteringPowerRegistry("Global ScatteringPower Registry");
89 
90 ScatteringPower::ScatteringPower():mDynPopCorrIndex(0),mBiso(1.0),mIsIsotropic(true),
91 mMaximumLikelihoodNbGhost(0),mFormalCharge(0.0)
92 {
93  VFN_DEBUG_MESSAGE("ScatteringPower::ScatteringPower():"<<mName,5)
94  mBeta.resize(6);
95  mBeta = 0;
96  mB.resize(6);
97  mB = 0;
98  gScatteringPowerRegistry.Register(*this);
99  this->Init();
100  mClockMaster.AddChild(mClock);
101  mClockMaster.AddChild(mMaximumLikelihoodParClock);
102 }
103 ScatteringPower::ScatteringPower(const ScatteringPower& old):
104 mDynPopCorrIndex(old.mDynPopCorrIndex),mBiso(old.mBiso),mIsIsotropic(old.mIsIsotropic),
105 mBeta(old.mBeta),mB(old.mB),
106 mFormalCharge(old.mFormalCharge)
107 {
108  VFN_DEBUG_MESSAGE("ScatteringPower::ScatteringPower(&old):"<<mName,5)
109  gScatteringPowerRegistry.Register(*this);
110  this->Init();
111  mMaximumLikelihoodPositionError=old.mMaximumLikelihoodPositionError;
112  mMaximumLikelihoodNbGhost=old.mMaximumLikelihoodNbGhost;
113  mClockMaster.AddChild(mClock);
114  mClockMaster.AddChild(mMaximumLikelihoodParClock);
115 }
116 ScatteringPower::~ScatteringPower()
117 {
118  VFN_DEBUG_MESSAGE("ScatteringPower::~ScatteringPower():"<<mName,5)
119  gScatteringPowerRegistry.DeRegister(*this);
120 }
121 
122 const string& ScatteringPower::GetClassName()const
123 {
124  const static string className="ScatteringPower";
125  return className;
126 }
127 
128 void ScatteringPower::operator=(const ScatteringPower& rhs)
129 {
130  VFN_DEBUG_MESSAGE("ScatteringPower::operator=():"<<mName,2)
131  mDynPopCorrIndex=rhs.mDynPopCorrIndex;
132  mBiso=rhs.mBiso;
133  mIsIsotropic=rhs.mIsIsotropic;
134  mBeta=rhs.mBeta;
135  mB=rhs.mB;
136 }
137 
138 bool ScatteringPower::operator==(const ScatteringPower& rhs) const
139 {
140  if(mBiso!=rhs.GetBiso()) return false;
141  for(unsigned int i=0;i<6;i++)
142  if(this->GetBij(i) != rhs.GetBij(i)) return false;
143  if(this->GetClassName() != rhs.GetClassName()) return false;
144  if(this->GetSymbol() != rhs.GetSymbol()) return false;
145  return true;
146 }
147 
148 bool ScatteringPower::operator!=(const ScatteringPower& rhs) const
149 {
150  return !(*this == rhs);
151 }
152 
153 bool ScatteringPower::IsScatteringFactorAnisotropic()const{return false;}
154 bool ScatteringPower::IsTemperatureFactorAnisotropic()const{return false;}
155 bool ScatteringPower::IsResonantScatteringAnisotropic()const{return false;}
156 
157 const string& ScatteringPower::GetSymbol() const {return this->GetName();}
158 REAL ScatteringPower::GetBiso() const {return mBiso;}
159 REAL& ScatteringPower::GetBiso() {mClock.Click();return mBiso;}
160 void ScatteringPower::SetBiso(const REAL newB) { mClock.Click();mBiso=newB;mIsIsotropic=true;}
161 REAL ScatteringPower::GetBij(const size_t &i, const size_t &j) const
162 {
163  size_t idx = 0;
164  if(i == j)
165  {
166  idx = i - 1;
167  }
168  else
169  {
170  idx = i + j;
171  }
172  return this->GetBij(idx);
173 }
174 REAL ScatteringPower::GetBij(const size_t &idx) const
175 {
176  return mB(idx);
177 }
178 void ScatteringPower::SetBij(const size_t &i, const size_t &j, const REAL newB)
179 {
180  size_t idx = 0;
181  if(i == j)
182  {
183  idx = i - 1;
184  }
185  else
186  {
187  idx = i + j;
188  }
189  this->SetBij(idx, newB);
190 }
191 void ScatteringPower::SetBij(const size_t &idx, const REAL newB)
192 {
193  mClock.Click();
194  mIsIsotropic=false;
195  mB(idx) = newB;
196 }
197 bool ScatteringPower::IsIsotropic() const {return mIsIsotropic;}
198 long ScatteringPower::GetDynPopCorrIndex() const {return mDynPopCorrIndex;}
199 long ScatteringPower::GetNbScatteringPower()const {return gScatteringPowerRegistry.GetNb();}
200 const RefinableObjClock& ScatteringPower::GetLastChangeClock()const {return mClock;}
201 
202 const string& ScatteringPower::GetColourName()const{ return mColourName;}
203 const float* ScatteringPower::GetColourRGB()const{ return mColourRGB;}
204 void ScatteringPower::SetColour(const string& colourName)
205 {
206  mColourName=colourName;
207  this->InitRGBColour();
208 }
209 void ScatteringPower::SetColour(const float r,const float g,const float b)
210 {
211  mColourRGB[0]=r;
212  mColourRGB[1]=g;
213  mColourRGB[2]=b;
214 }
215 void ScatteringPower::GetGeneGroup(const RefinableObj &obj,
216  CrystVector_uint & groupIndex,
217  unsigned int &first) const
218 {
219  // One group for all parameters
220  unsigned int index=0;
221  VFN_DEBUG_MESSAGE("ScatteringPower::GetGeneGroup()",4)
222  for(long i=0;i<obj.GetNbPar();i++)
223  for(long j=0;j<this->GetNbPar();j++)
224  if(&(obj.GetPar(i)) == &(this->GetPar(j)))
225  {
226  if(index==0) index=first++;
227  groupIndex(i)=index;
228  }
229 }
230 
231 REAL ScatteringPower::GetMaximumLikelihoodPositionError()const
232 {return mMaximumLikelihoodPositionError;}
233 
234 const RefinableObjClock& ScatteringPower::GetMaximumLikelihoodParClock()const
235 {return mMaximumLikelihoodParClock;}
236 
237 void ScatteringPower::SetMaximumLikelihoodPositionError(const REAL mle)
238 {
239  if(mle!=mMaximumLikelihoodPositionError)
240  {
241  mMaximumLikelihoodPositionError=mle;
242  mMaximumLikelihoodParClock.Click();
243  }
244 }
245 
246 REAL ScatteringPower::GetMaximumLikelihoodNbGhostAtom()const
247 {return mMaximumLikelihoodNbGhost;}
248 
249 void ScatteringPower::SetMaximumLikelihoodNbGhostAtom(const REAL nb)
250 {
251  if(nb!=mMaximumLikelihoodNbGhost)
252  {
253  mMaximumLikelihoodNbGhost=nb;
254  mMaximumLikelihoodParClock.Click();
255  }
256 }
257 
258 REAL ScatteringPower::GetFormalCharge()const{return mFormalCharge;}
259 void ScatteringPower::SetFormalCharge(const REAL charge)
260 {mFormalCharge=charge;}
261 
262 void ScatteringPower::Init()
263 {
264  VFN_DEBUG_MESSAGE("ScatteringPower::Init():"<<mName,2)
265  mColourName="White";
266  mMaximumLikelihoodPositionError=0;
267  mMaximumLikelihoodNbGhost=0;
268  VFN_DEBUG_MESSAGE("ScatteringPower::Init():End",2)
269 }
270 void ScatteringPower::InitRGBColour()
271 {
272  VFN_DEBUG_MESSAGE("ScatteringPower::InitRGBColour()",2)
273  for(long i=0;;)
274  {
275  if(gPOVRayColours[i].mName==mColourName)
276  {
277  mColourRGB[0]=gPOVRayColours[i].mRGB[0];
278  mColourRGB[1]=gPOVRayColours[i].mRGB[1];
279  mColourRGB[2]=gPOVRayColours[i].mRGB[2];
280  break;
281  }
282  i++;
283  if(strncmp(gPOVRayColours[i].mName,"",3)==0)
284  {//could not find colour !
285  cout << "Could not find colour:"<<mColourName<<" for ScatteringPower "<<mName<<endl;
286  mColourRGB[0]=1;
287  mColourRGB[1]=1;
288  mColourRGB[2]=1;
289  break;
290  }
291  }
292  VFN_DEBUG_MESSAGE("->RGBColour:"<<mColourName<<mColourRGB[0]<<" "<<mColourRGB[1]<<" "<<mColourRGB[2],2)
293 }
294 
295 //######################################################################
296 //
297 // SCATTERING POWER ATOM
298 //
299 //######################################################################
301  gScatteringPowerAtomRegistry("Global ScatteringPowerAtom Registry");
302 
303 ScatteringPowerAtom::ScatteringPowerAtom():
304 ScatteringPower(),mSymbol(""),mAtomicNumber(0),mpGaussian(0)
305 {
306  VFN_DEBUG_MESSAGE("ScatteringPowerAtom::ScatteringPowerAtom():"<<mName,5)
307  gScatteringPowerAtomRegistry.Register(*this);
308  this->InitRefParList();
309 }
310 
312  const string &symbol,
313  const REAL bIso):
314 mpGaussian(0)
315 {
316  VFN_DEBUG_MESSAGE("ScatteringPowerAtom::ScatteringPowerAtom(n,s,B):"<<name,5)
317  gScatteringPowerAtomRegistry.Register(*this);
318  this->InitRefParList();
319  this->Init(name,symbol,bIso);
320 }
321 
322 ScatteringPowerAtom::ScatteringPowerAtom(const ScatteringPowerAtom& old):
323 mpGaussian(0)
324 {
325  VFN_DEBUG_MESSAGE("ScatteringPowerAtom::ScatteringPowerAtom(&old):"<<old.mName,5)
326  gScatteringPowerAtomRegistry.Register(*this);
327  this->Init(old.GetName(),old.mSymbol,old.mBiso);
328  //this->InitRefParList(); //?? :TODO: Check
329 }
330 
332 {
333  VFN_DEBUG_MESSAGE("ScatteringPowerAtom::~ScatteringPowerAtom():"<<mName,5)
334  gScatteringPowerAtomRegistry.DeRegister(*this);
335  delete mpGaussian;
336 }
337 
338 const string& ScatteringPowerAtom::GetClassName() const
339 {
340  const static string className="ScatteringPowerAtom";
341  return className;
342 }
343 
344 // Disable the base-class function.
346 {
347  // This should be never called.
348  abort();
349 }
350 
351 void ScatteringPowerAtom::Init(const string &name,const string &symbol,const REAL bIso)
352 {
353  VFN_DEBUG_MESSAGE("ScatteringPowerAtom::Init(n,s,b)"<<mName,4)
354  this->ScatteringPower::Init();
355  this->SetName(name);
356  mSymbol=symbol;
357  mBiso=bIso;
358  mIsIsotropic=true;
359  if(mpGaussian!=0) delete mpGaussian;
360  try
361  {
362  cctbx::eltbx::xray_scattering::wk1995 wk95t(mSymbol);
363  mpGaussian=new cctbx::eltbx::xray_scattering::gaussian(wk95t.fetch());
364 
365  this->InitAtNeutronScattCoeffs();
366 
367  cctbx::eltbx::tiny_pse::table tpse(mSymbol);
368  mAtomicNumber=tpse.atomic_number();
369  mAtomicWeight=tpse.weight();
370 
371  cctbx::eltbx::icsd_radii::table ticsd(mSymbol);
372  mRadius= ticsd.radius();
373  cctbx::eltbx::covalent_radii::table tcov(mSymbol);
374  mCovalentRadius=tcov.radius();
375  }
376  catch(exception &err)
377  {
378  cout << "WARNING: could not interpret Symbol name !"<<mSymbol<<endl
379  << " Reverting to H !"<<endl;
380  (*fpObjCrystInformUser)("Symbol not understood:"+mSymbol);
381  this->Init(name,"H",bIso);
382  }
383 
384 
385  VFN_DEBUG_MESSAGE("ScatteringPowerAtom::Init():/Name="<<this->GetName() \
386  <<" /Symbol="<<mSymbol<<" /Atomic Number=" << mAtomicNumber,4)
387 
389 
390  //Init default atom colours for POVRay/GUI
391  cctbx::eltbx::tiny_pse::table tpse(mSymbol);
392  mColourName= tpse.symbol();
393  this->InitRGBColour();
394  switch(mAtomicNumber)
395  {// Values from OpenBabel element.txt
396  case 0: mMaxCovBonds=0;break;break;break;//Xx
397  case 1: mMaxCovBonds=1;break;break;break;//H
398  case 2: mMaxCovBonds=0;break;break;break;//He
399  case 3: mMaxCovBonds=1;break;break;break;//Li
400  case 4: mMaxCovBonds=2;break;//Be
401  case 5: mMaxCovBonds=4;break;//B
402  case 6: mMaxCovBonds=4;break;//C
403  case 7: mMaxCovBonds=4;break;//N
404  case 8: mMaxCovBonds=2;break;//O
405  case 9: mMaxCovBonds=1;break;//F
406  case 10: mMaxCovBonds=0;break;//Ne
407  case 11: mMaxCovBonds=1;break;//Na
408  case 12: mMaxCovBonds=2;break;//Mg
409  case 13: mMaxCovBonds=6;break;//Al
410  case 14: mMaxCovBonds=6;break;//Si
411  case 15: mMaxCovBonds=6;break;//P
412  case 16: mMaxCovBonds=6;break;//S
413  case 17: mMaxCovBonds=1;break;//Cl
414  case 18: mMaxCovBonds=0;break;//Ar
415  case 19: mMaxCovBonds=1;break;//K
416  case 20: mMaxCovBonds=2;break;//Ca
417  case 21: mMaxCovBonds=6;break;//Sc
418  case 22: mMaxCovBonds=6;break;//Ti
419  case 23: mMaxCovBonds=6;break;//V
420  case 24: mMaxCovBonds=6;break;//Cr
421  case 25: mMaxCovBonds=8;break;//Mn
422  case 26: mMaxCovBonds=6;break;//Fe
423  case 27: mMaxCovBonds=6;break;//Co
424  case 28: mMaxCovBonds=6;break;//Ni
425  case 29: mMaxCovBonds=6;break;//Cu
426  case 30: mMaxCovBonds=6;break;//Zn
427  case 31: mMaxCovBonds=3;break;//Ga
428  case 32: mMaxCovBonds=4;break;//Ge
429  case 33: mMaxCovBonds=3;break;//As
430  case 34: mMaxCovBonds=2;break;//Se
431  case 35: mMaxCovBonds=1;break;//Br
432  case 36: mMaxCovBonds=0;break;//Kr
433  case 37: mMaxCovBonds=1;break;//Rb
434  case 38: mMaxCovBonds=2;break;//Sr
435  case 39: mMaxCovBonds=6;break;//Y
436  case 40: mMaxCovBonds=6;break;//Zr
437  case 41: mMaxCovBonds=6;break;//Nb
438  case 42: mMaxCovBonds=6;break;//Mo
439  case 43: mMaxCovBonds=6;break;//Tc
440  case 44: mMaxCovBonds=6;break;//Ru
441  case 45: mMaxCovBonds=6;break;//Rh
442  case 46: mMaxCovBonds=6;break;//Pd
443  case 47: mMaxCovBonds=6;break;//Ag
444  case 48: mMaxCovBonds=6;break;//Cd
445  case 49: mMaxCovBonds=3;break;//In
446  case 50: mMaxCovBonds=4;break;//Sn
447  case 51: mMaxCovBonds=3;break;//Sb
448  case 52: mMaxCovBonds=2;break;//Te
449  case 53: mMaxCovBonds=1;break;//I
450  case 54: mMaxCovBonds=0;break;//Xe
451  case 55: mMaxCovBonds=1;break;//Cs
452  case 56: mMaxCovBonds=2;break;//Ba
453  case 57: mMaxCovBonds=12;break;//La
454  case 58: mMaxCovBonds=6;break;//Ce
455  case 59: mMaxCovBonds=6;break;//Pr
456  case 60: mMaxCovBonds=6;break;//Nd
457  case 61: mMaxCovBonds=6;break;//Pm
458  case 62: mMaxCovBonds=6;break;//Sm
459  case 63: mMaxCovBonds=6;break;//Eu
460  case 64: mMaxCovBonds=6;break;//Gd
461  case 65: mMaxCovBonds=6;break;//Tb
462  case 66: mMaxCovBonds=6;break;//Dy
463  case 67: mMaxCovBonds=6;break;//Ho
464  case 68: mMaxCovBonds=6;break;//Er
465  case 69: mMaxCovBonds=6;break;//Tm
466  case 70: mMaxCovBonds=6;break;//Yb
467  case 71: mMaxCovBonds=6;break;//Lu
468  case 72: mMaxCovBonds=6;break;//Hf
469  case 73: mMaxCovBonds=6;break;//Ta
470  case 74: mMaxCovBonds=6;break;//W
471  case 75: mMaxCovBonds=6;break;//Re
472  case 76: mMaxCovBonds=6;break;//Os
473  case 77: mMaxCovBonds=6;break;//Ir
474  case 78: mMaxCovBonds=6;break;//Pt
475  case 79: mMaxCovBonds=6;break;//Au
476  case 80: mMaxCovBonds=6;break;//Hg
477  case 81: mMaxCovBonds=3;break;//Tl
478  case 82: mMaxCovBonds=4;break;//Pb
479  case 83: mMaxCovBonds=3;break;//Bi
480  case 84: mMaxCovBonds=2;break;//Po
481  case 85: mMaxCovBonds=1;break;//At
482  case 86: mMaxCovBonds=0;break;//Rn
483  case 87: mMaxCovBonds=1;break;//Fr
484  case 88: mMaxCovBonds=2;break;//Ra
485  default: mMaxCovBonds=6;break;
486  }
487 
488  VFN_DEBUG_MESSAGE("ScatteringPowerAtom::Init(n,s,b):End",3)
489 }
490 
492  const int spgSymPosIndex) const
493 {
494  VFN_DEBUG_MESSAGE("ScatteringPower::GetScatteringFactor(&data):"<<mName,3)
495  CrystVector_REAL sf(data.GetNbRefl());
496  switch(data.GetRadiationType())
497  {
498  case(RAD_NEUTRON):
499  {
500  VFN_DEBUG_MESSAGE("ScatteringPower::GetScatteringFactor():NEUTRON:"<<mName,3)
502  break;
503  }
504  case(RAD_XRAY):
505  {
506  VFN_DEBUG_MESSAGE("ScatteringPower::GetScatteringFactor():XRAY:"<<mName,3)
507  if(mpGaussian!=0)
508  {
509  const long nb=data.GetSinThetaOverLambda().numElements();
510  const REAL *pstol=data.GetSinThetaOverLambda().data();
511  for(long i=0;i<nb;i++)
512  sf(i)=mpGaussian->at_stol(*pstol++);
513  }
514  else sf=1.0;//:KLUDGE: Should never happen
515  break;
516  }
517  case(RAD_ELECTRON):
518  {
519  VFN_DEBUG_MESSAGE("ScatteringPower::GetScatteringFactor():ELECTRON:"<<mName,3)
520  if(mpGaussian!=0)
521  {
522  const REAL z=this->GetAtomicNumber();
523  const long nb=data.GetSinThetaOverLambda().numElements();
524  const REAL *pstol=data.GetSinThetaOverLambda().data();
525  for(long i=0;i<nb;i++)
526  {
527  sf(i)=(z-mpGaussian->at_stol(*pstol))/(*pstol * *pstol);
528  pstol++;
529  }
530  }
531  else sf=1.0;//:KLUDGE: Should never happen
532  break;
533  }
534  }
535  VFN_DEBUG_MESSAGE("ScatteringPower::GetScatteringFactor(&data):End",3)
536  return sf;
537 }
538 
540 {
541  REAL sf = 0;
542  switch(type)
543  {
544  case(RAD_NEUTRON):
545  {
547  break;
548  }
549  case(RAD_XRAY):
550  {
551  if(mpGaussian!=0)
552  sf=mpGaussian->at_stol(0);
553  else sf=1.0;
554  break;
555  }
556  case(RAD_ELECTRON):
557  {
558  const REAL z=this->GetAtomicNumber();
559  sf=(z-mpGaussian->at_stol(0.0001))/(.0001 * .0001);
560  }
561  }
562  VFN_DEBUG_MESSAGE("ScatteringPower::GetScatteringFactor(&data):End",3)
563  return sf;
564 }
565 
566 static bool warnADP=true;
568  const int spgSymPosIndex) const
569 {
570  VFN_DEBUG_MESSAGE("ScatteringPower::GetTemperatureFactor(&data):"<<mName,3)
571  CrystVector_REAL sf(data.GetNbRefl());
572  if((mIsIsotropic==false) && warnADP)
573  { // Warn once
574  cout<<"========================== WARNING ========================="<<endl
575  <<" In ScatteringPowerAtom::GetTemperatureFactor():"<<endl
576  <<" Anisotropic Displacement Parameters are not currently properly handled"<<endl
577  <<" for Debye-Waller calculations (no symmetry handling for ADPs)."<<endl
578  <<" =>The Debye-Waller calculations will instead use only isotropic DPs"<<endl<<endl;
579  warnADP=false;
580  }
581 
582  if(true)//(mIsIsotropic)
583  {
584  CrystVector_REAL stolsq(data.GetNbRefl());
585  const CrystVector_REAL stol=data.GetSinThetaOverLambda();
586  stolsq=stol;
587  stolsq*=stol;
588 
589  #ifdef __VFN_VECTOR_USE_BLITZ__
590  #define SF sf
591  #define STOLSQ stolsq
592  #else
593  #define SF (*ssf)
594  #define STOLSQ (*sstolsq)
595 
596  REAL *ssf=sf.data();
597  const REAL *sstolsq=stolsq.data();
598 
599  for(long ii=0;ii<sf.numElements();ii++)
600  {
601  #endif
602 
603  SF=exp(-mBiso*STOLSQ);
604 
605  #ifdef __VFN_VECTOR_USE_BLITZ__
606 
607  #else
608  ssf++;
609  sstolsq++;
610  }
611  #endif
612 
613  #undef SF
614  #undef STOLSQ
615  }
616  else
617  {// :TODO: handle ADP - requires taking into account symmetries...
618  const REAL b11=mBeta(0);
619  const REAL b22=mBeta(1);
620  const REAL b33=mBeta(2);
621  const REAL b12=mBeta(3);
622  const REAL b13=mBeta(4);
623  const REAL b23=mBeta(5);
624  #ifdef __VFN_VECTOR_USE_BLITZ__
625  #define HH data.H()
626  #define KK data.K()
627  #define LL data.L()
628  #define SF sf
629  #else
630  #define HH (*hh)
631  #define KK (*kk)
632  #define LL (*ll)
633  #define SF (*ssf)
634 
635  const REAL *hh=(data.GetH()).data();
636  const REAL *kk=(data.GetK()).data();
637  const REAL *ll=(data.GetL()).data();
638  REAL *ssf=sf.data();
639 
640  for(long ii=0;ii<sf.numElements();ii++)
641  {
642  #endif
643 
644  SF= exp( -b11*pow(HH,2)
645  -b22*pow(KK,2)
646  -b33*pow(LL,2)
647  -2*b12*HH*KK
648  -2*b13*HH*LL
649  -2*b23*KK*LL);
650 
651  #ifdef __VFN_VECTOR_USE_BLITZ__
652 
653  #else
654  hh++;
655  kk++;
656  ll++;
657  ssf++;
658  }
659  #endif
660 
661  #undef HH
662  #undef KK
663  #undef LL
664  #undef SF
665  }
666  return sf;
667 }
668 
669 CrystMatrix_REAL ScatteringPowerAtom::
671  const int spgSymPosIndex) const
672 {
673  VFN_DEBUG_MESSAGE("ScatteringPower::GetResonantScattFactReal(&data):"<<mName,3)
674  CrystMatrix_REAL fprime(1,1);//:TODO: More than one wavelength
675  CrystMatrix_REAL fsecond(1,1);
676  switch(data.GetRadiationType())
677  {
678  case(RAD_NEUTRON):
679  {
680  fprime=0;
681  fsecond=mNeutronScattLengthImag;
682  break;
683  }
684  case(RAD_XRAY):
685  {
686  try
687  {
688  cctbx::eltbx::henke::table thenke(mSymbol);
689  cctbx::eltbx::fp_fdp f=thenke.at_angstrom(data.GetWavelength()(0));
690 
691  if(f.is_valid_fp()) fprime(0)=f.fp();
692  else fprime(0)=0;
693  if(f.is_valid_fdp()) fsecond(0)=f.fdp();
694  else fsecond(0)=0;
695  }
696  catch(cctbx::error)
697  {
698  fprime(0)=0;
699  fsecond(0)=0;
700  }
701  break;
702  }
703  case(RAD_ELECTRON):
704  {
705  fprime=0;
706  fsecond=0;
707  break;
708  }
709  }
710  return fprime;
711 }
712 
713 CrystMatrix_REAL ScatteringPowerAtom::
715  const int spgSymPosIndex) const
716 {
717  VFN_DEBUG_MESSAGE("ScatteringPower::GetResonantScattFactImag():"<<mName,3)
718  CrystMatrix_REAL fprime(1,1);//:TODO: More than one wavelength
719  CrystMatrix_REAL fsecond(1,1);
720  switch(data.GetRadiationType())
721  {
722  case(RAD_NEUTRON):
723  {
724  fprime=0;
725  fsecond=mNeutronScattLengthImag;
726  break;
727  }
728  case(RAD_XRAY):
729  {
730  try
731  {
732  cctbx::eltbx::henke::table thenke(mSymbol);
733  cctbx::eltbx::fp_fdp f=thenke.at_angstrom(data.GetWavelength()(0));
734 
735  if(f.is_valid_fp()) fprime(0)=f.fp();
736  else fprime(0)=0;
737  if(f.is_valid_fdp()) fsecond(0)=f.fdp();
738  else fsecond(0)=0;
739  }
740  catch(cctbx::error)
741  {
742  fprime(0)=0;
743  fsecond(0)=0;
744  }
745  break;
746  }
747  case(RAD_ELECTRON):
748  {
749  fprime=0;
750  fsecond=0;
751  break;
752  }
753  }
754  return fsecond;
755 }
756 
757 
758 void ScatteringPowerAtom::SetSymbol(const string& symbol)
759 {
760  VFN_DEBUG_MESSAGE("ScatteringPowerAtom::SetSymbol():"<<mName,5)
761  this->Init(this->GetName(),symbol,this->GetBiso());
762 }
763 const string& ScatteringPowerAtom::GetSymbol() const
764 {
765  VFN_DEBUG_MESSAGE("ScatteringPowerAtom::GetSymbol():"<<mName,5)
766  return mSymbol;
767 }
768 
770 {
771  VFN_DEBUG_MESSAGE("ScatteringPowerAtom::GetElementName():"<<mName,2)
772  try
773  {
774  cctbx::eltbx::tiny_pse::table tpse(mSymbol);
775  return tpse.name();
776  }
777  catch(cctbx::error)
778  {
779  cout << "WARNING: could not interpret Symbol:"<<mSymbol<<endl;
780  }
781  return "Unknown";
782 }
783 
789 
790 void ScatteringPowerAtom::Print()const
791 {
792  VFN_DEBUG_MESSAGE("ScatteringPowerAtom::Print()",1)
793  cout << "ScatteringPowerAtom ("<<this->GetName()<<","
794  << FormatString(this->GetSymbol(),4) << ") :"
795  << FormatFloat(this->GetBiso());
796  VFN_DEBUG_MESSAGE_SHORT("at "<<this,10)
797  cout << endl;
798 }
799 
801 {
802  VFN_DEBUG_MESSAGE("ScatteringPowerAtom::InitAtNeutronScattCoeffs():"<<mName,3)
803  mClock.Click();
804  try
805  {
806  cctbx::eltbx::neutron::neutron_news_1992_table nn92t(mSymbol);
807  mNeutronScattLengthReal=nn92t.bound_coh_scatt_length_real();
808  mNeutronScattLengthImag=nn92t.bound_coh_scatt_length_imag();
809  }
810  catch(cctbx::error)
811  {
812  cout << "WARNING: could not interpret symbol for neutron coeefs:"<<mSymbol<<endl;
813  }
814 
815  VFN_DEBUG_MESSAGE("ScatteringPowerAtom::InitAtNeutronScattCoeffs():End",3)
816 }
817 
818 void ScatteringPowerAtom::InitRefParList()
819 {
820  VFN_DEBUG_MESSAGE("ScatteringPowerAtom::InitRefParList():"<<mName,5)
821  {
822  RefinablePar tmp("Biso",&mBiso,0.1,5.,
823  gpRefParTypeScattPowTemperatureIso,REFPAR_DERIV_STEP_ABSOLUTE,
824  true,true,true,false);
825  tmp.SetDerivStep(1e-3);
826  tmp.SetGlobalOptimStep(.5);
827  tmp.AssignClock(mClock);
828  this->AddPar(tmp);
829  }
830  {
831  REAL* bdata = (REAL*) mB.data();
832 
833  RefinablePar B11("B11",&bdata[0],0.1,5.,
834  gpRefParTypeScattPowTemperatureAniso,REFPAR_DERIV_STEP_ABSOLUTE,
835  true,true,false,false);
836  B11.SetDerivStep(1e-3);
837  B11.SetGlobalOptimStep(.5);
838  B11.AssignClock(mClock);
839  this->AddPar(B11);
840 
841  RefinablePar B22("B22",&bdata[1],0.1,5.,
842  gpRefParTypeScattPowTemperatureAniso,REFPAR_DERIV_STEP_ABSOLUTE,
843  true,true,false,false);
844  B22.SetDerivStep(1e-3);
845  B22.SetGlobalOptimStep(.5);
846  B22.AssignClock(mClock);
847  this->AddPar(B22);
848 
849  RefinablePar B33("B33",&bdata[2],0.1,5.,
850  gpRefParTypeScattPowTemperatureAniso,REFPAR_DERIV_STEP_ABSOLUTE,
851  true,true,false,false);
852  B33.SetDerivStep(1e-3);
853  B33.SetGlobalOptimStep(.5);
854  B33.AssignClock(mClock);
855  this->AddPar(B33);
856 
857  RefinablePar B12("B12",&bdata[3],-5,5.,
858  gpRefParTypeScattPowTemperatureAniso,REFPAR_DERIV_STEP_ABSOLUTE,
859  true,true,false,false);
860  B12.SetDerivStep(1e-3);
861  B12.SetGlobalOptimStep(.5);
862  B12.AssignClock(mClock);
863  this->AddPar(B12);
864 
865  RefinablePar B13("B13",&bdata[4],-5,5.,
866  gpRefParTypeScattPowTemperatureAniso,REFPAR_DERIV_STEP_ABSOLUTE,
867  true,true,false,false);
868  B13.SetDerivStep(1e-3);
869  B13.SetGlobalOptimStep(.5);
870  B13.AssignClock(mClock);
871  this->AddPar(B13);
872 
873  RefinablePar B23("B23",&bdata[5],-5,5.,
874  gpRefParTypeScattPowTemperatureAniso,REFPAR_DERIV_STEP_ABSOLUTE,
875  true,true,false,false);
876  B23.SetDerivStep(1e-3);
877  B23.SetGlobalOptimStep(.5);
878  B23.AssignClock(mClock);
879  this->AddPar(B23);
880  }
881  {
882  RefinablePar tmp("ML Error",&mMaximumLikelihoodPositionError,0.,1.,
883  gpRefParTypeScattPow,REFPAR_DERIV_STEP_ABSOLUTE,
884  false,true,true,false);
885  tmp.SetDerivStep(1e-4);
886  tmp.SetGlobalOptimStep(.001);
887  tmp.AssignClock(mMaximumLikelihoodParClock);
888  this->AddPar(tmp);
889  }
890  {
891  RefinablePar tmp("ML-Nb Ghost Atoms",&mMaximumLikelihoodNbGhost,0.,10.,
892  gpRefParTypeScattPow,REFPAR_DERIV_STEP_ABSOLUTE,
893  true,true,true,false);
894  tmp.SetDerivStep(1e-3);
895  tmp.SetGlobalOptimStep(.05);
896  tmp.AssignClock(mMaximumLikelihoodParClock);
897  this->AddPar(tmp);
898  }
899  {
900  RefinablePar tmp("Formal Charge",&mFormalCharge,-10.,10.,
901  gpRefParTypeScattPow,REFPAR_DERIV_STEP_ABSOLUTE,
902  true,true,true,false);
903  tmp.SetDerivStep(1e-3);
904  tmp.SetGlobalOptimStep(.05);
905  tmp.AssignClock(mClock);
906  this->AddPar(tmp);
907  }
908 }
909 #ifdef __WX__CRYST__
910 WXCrystObjBasic* ScatteringPowerAtom::WXCreate(wxWindow* parent)
911 {
912  //:TODO: Check mpWXCrystObj==0
913  mpWXCrystObj=new WXScatteringPowerAtom(parent,this);
914  return mpWXCrystObj;
915 }
916 #endif
917 
918 //######################################################################
919 //
920 // SCATTERING COMPONENT
921 //
922 //######################################################################
923 ScatteringComponent::ScatteringComponent():
924 mX(0),mY(0),mZ(0),mOccupancy(0),mpScattPow(0),mDynPopCorr(0)
925 {}
926 bool ScatteringComponent::operator==(const ScatteringComponent& rhs) const
927 {
928  return ((mX==rhs.mX) && (mY==rhs.mY) && (mZ==rhs.mZ) &&
929  (mOccupancy==rhs.mOccupancy) && (mpScattPow==rhs.mpScattPow));
930 }
931 bool ScatteringComponent::operator!=(const ScatteringComponent& rhs) const
932 {
933  return ((mX!=rhs.mX) || (mY!=rhs.mY) || (mZ!=rhs.mZ) ||
934  (mOccupancy!=rhs.mOccupancy) || (mpScattPow!=rhs.mpScattPow));
935 }
937 {
938  cout <<mX<<" "<<mY<<" "<<mZ<<" "<<mOccupancy<<" "<<mDynPopCorr<<" "<<mpScattPow;
939  if(0!=mpScattPow) cout <<" "<<mpScattPow->GetName();
940  cout<<endl;
941 }
942 
943 //######################################################################
944 //
945 // SCATTERING COMPONENT LIST
946 //
947 //######################################################################
948 
949 ScatteringComponentList::ScatteringComponentList()
950 {
951 }
952 
953 ScatteringComponentList::ScatteringComponentList(const long nbComponent)
954 {
955  mvScattComp.resize(nbComponent);
956 }
957 
958 ScatteringComponentList::ScatteringComponentList(const ScatteringComponentList &old):
959 mvScattComp(old.mvScattComp)
960 {
961 }
962 
963 ScatteringComponentList::~ScatteringComponentList()
964 {
965 }
967 {
968  mvScattComp.clear();
969 }
970 
972 {
973  VFN_DEBUG_MESSAGE("ScatteringComponentList::operator()("<<i<<")",1)
974  if(i>=this->GetNbComponent())
975  {
976  this->Print();
977  throw ObjCrystException("ScatteringComponentList::operator()(i)::i>mNbComponent!!");
978  }
979  if(i<0) throw ObjCrystException("ScatteringComponentList::operator()&(i)::i<0!!");
980  return mvScattComp[i];
981 }
982 
984 {
985  VFN_DEBUG_MESSAGE("ScatteringComponentList::operator()&("<<i<<")",1)
986  if(i>=this->GetNbComponent())
987  {
988  this->Print();
989  throw ObjCrystException("ScatteringComponentList::operator()&(i)::i>mNbComponent!!");
990  }
991  if(i<0) throw ObjCrystException("ScatteringComponentList::operator()&(i):: i<0!!");
992  return mvScattComp[i];
993 }
994 
996 
998 {
999  VFN_DEBUG_MESSAGE("ScatteringComponentList::operator=()",1)
1001  VFN_DEBUG_MESSAGE("ScatteringComponentList::operator=():End",0)
1002 }
1004 {
1005  if(rhs.GetNbComponent() != this->GetNbComponent()) return false;
1006  for(long i=0;i<this->GetNbComponent();i++)
1007  if( (*this)(i) != rhs(i) ) return false;
1008  return true;
1009 }
1010 
1012 {
1013  for(long i=0;i<rhs.GetNbComponent();i++)
1014  mvScattComp.push_back(rhs(i));
1015 }
1017 {
1018  VFN_DEBUG_MESSAGE("ScatteringComponentList::operator+=()",1)
1019  mvScattComp.push_back(rhs);
1020 }
1021 
1023 {
1024  VFN_DEBUG_MESSAGE("ScatteringComponentList::operator++()",1)
1025  mvScattComp.resize(this->GetNbComponent()+1);
1026 }
1027 
1029 {
1030  VFN_DEBUG_MESSAGE("ScatteringComponentList::operator--()",1)
1031  if(this->GetNbComponent()>0)
1032  mvScattComp.resize(this->GetNbComponent()-1);
1033 }
1034 
1036 {
1037  VFN_DEBUG_ENTRY("ScatteringComponentList::Print()",5)
1038  cout<<"Number of Scattering components:"<<this->GetNbComponent()<<endl;
1039  for(long i=0;i<this->GetNbComponent();i++)
1040  {
1041  cout << i<<":";
1042  (*this)(i).Print();
1043  }
1044  VFN_DEBUG_EXIT("ScatteringComponentList::Print()",5)
1045 }
1046 
1047 }//namespace
The namespace which includes all objects (crystallographic and algorithmic) in ObjCryst++.
Definition: doc-main.h:25
RadiationType
Type of radiation used.
Definition: General.h:96
ObjRegistry< ScatteringPower > gScatteringPowerRegistry("Global ScatteringPower Registry")
Global registry for all ScatteringPower objects.
ObjRegistry< ScatteringPowerAtom > gScatteringPowerAtomRegistry("Global ScatteringPowerAtom Registry")
Global registry for all ScatteringPowerAtom objects.
Exception class for ObjCryst++ library.
Definition: General.h:122
Class to compute structure factors for a set of reflections and a Crystal.
const CrystVector_REAL & GetSinThetaOverLambda() const
Return an array with for all reflections.
const CrystVector_REAL & GetK() const
Return the 1D array of K coordinates for all reflections.
CrystVector_REAL GetWavelength() const
wavelength of the experiment (in Angstroems)
const CrystVector_REAL & GetH() const
Return the 1D array of H coordinates for all reflections.
RadiationType GetRadiationType() const
Neutron or x-ray experiment ? Wavelength ?
const CrystVector_REAL & GetL() const
Return the 1D array of L coordinates for all reflections.
long GetNbRefl() const
Return the number of reflections in this experiment.
Abstract Base Class to describe the scattering power of any Scatterer component in a crystal.
long mDynPopCorrIndex
number identifying this kind of scatterer, for the dynamical occupancy correction.
CrystVector_REAL mBeta
Anisotropic Beta(ij)
REAL mFormalCharge
Formal Charge.
string mColourName
Colour for this ScatteringPower (from POVRay)
REAL mMaximumLikelihoodPositionError
estimated error (sigma) on the positions for this type of element.
virtual void Init()
Initialization of the object, used by all constructors, and operator=.
bool mIsIsotropic
Is the scattering isotropic ?
REAL GetBiso() const
Returns the isotropic temperature B factor.
REAL mMaximumLikelihoodNbGhost
Number of ghost atoms in the asymmetric unit.
CrystVector_REAL mB
Anisotropic B(ij)
RefinableObjClock mClock
Clock.
REAL mBiso
Temperature isotropic B factor.
virtual void InitRGBColour()
Get RGB Colour coordinates from Colour Name.
The Scattering Power for an Atom.
string GetElementName() const
Returns the standard name of the element (ie "hydrogen", "tantalum",..).
virtual const string & GetSymbol() const
Returns the symbol ('Ta', 'O2-',...) of the atom.
REAL mCovalentRadius
Covalent Radius for this atom, in Angstroems (from cctbx)
virtual CrystVector_REAL GetScatteringFactor(const ScatteringData &data, const int spgSymPosIndex=0) const
Get the Scattering factor for all reflections of a given ScatteringData object.
virtual CrystMatrix_REAL GetResonantScattFactReal(const ScatteringData &data, const int spgSymPosIndex=0) const
Get the real part of the resonant scattering factor.
REAL GetCovalentRadius() const
Covalent Radius for this atom, in Angstroems (from cctbx)
virtual CrystMatrix_REAL GetResonantScattFactImag(const ScatteringData &data, const int spgSymPosIndex=0) const
Get the imaginary part of the resonant scattering factor.
REAL mRadius
Radius of the atom or ion, in Angstroems (ICSD table from cctbx)
REAL GetAtomicWeight() const
Atomic weight (g/mol) for this atom.
void SetSymbol(const string &symbol)
Set the symbol for this atom.
REAL mNeutronScattLengthReal
Neutron Bond Coherent Scattering lengths.
string mSymbol
Symbol of this atom.
unsigned int GetMaxCovBonds() const
Maximum number of covalent bonds (from openbabel element.txt)
virtual const string & GetClassName() const
Name for this class ("RefinableObj", "Crystal",...).
REAL GetRadius() const
Atomic radius for this atom or ion, in Angstroems (ICSD table from cctbx)
void Init()
Initialization of the object, used by all constructors, and operator=.
virtual REAL GetForwardScatteringFactor(const RadiationType) const
Get the scattering factor at (0,0,0).
int mAtomicNumber
atomic number (Z) for the atom
int GetAtomicNumber() const
Atomic number for this atom.
virtual CrystVector_REAL GetTemperatureFactor(const ScatteringData &data, const int spgSymPosIndex=0) const
Get the temperature factor for all reflections of a given ScatteringData object.
REAL mAtomicWeight
atomic weight (g/mol) for the atom
unsigned int mMaxCovBonds
Maximum number of covalent bonds.
cctbx::eltbx::xray_scattering::gaussian * mpGaussian
Pointer to cctbx's gaussian describing the thomson x-ray scattering factor.
A scattering position in a crystal, associated with the corresponding occupancy and a pointer to the ...
const ScatteringPower * mpScattPow
The ScatteringPower associated with this position.
REAL mDynPopCorr
Dynamical Population Correction.
REAL mX
Coordinates of scattering positions i the crystal with the corresponding occupancy.
void Print() const
Print one line oabout this component.
list of scattering positions in a crystal, associated with the corresponding occupancy and a pointer ...
bool operator==(const ScatteringComponentList &rhs) const
Compare two lists.
vector< ScatteringComponent > mvScattComp
The vector of components.
const ScatteringComponent & operator()(const long i) const
Access to a component.
void operator+=(const ScatteringComponentList &rhs)
Add another list of components.
void operator++()
Add component (the whole list should be updated after that)
void operator=(const ScatteringComponentList &rhs)
Assignement operator.
long GetNbComponent() const
Number of components.
void Print() const
Print the list of Scattering components. For debugging.
void operator--()
Remove component (the whole list should be updated after that)
We need to record exactly when refinable objects have been modified for the last time (to avoid re-co...
Definition: RefinableObj.h:140
void Click()
Record an event for this clock (generally, the 'time' an object has been modified,...
Generic class for parameters of refinable objects.
Definition: RefinableObj.h:225
Object Registry.
Definition: RefinableObj.h:645
Generic Refinable Object.
Definition: RefinableObj.h:784
void AddPar(const RefinablePar &newRefPar)
Add a refinable parameter.
virtual void SetName(const string &name)
Name of the object.
RefinablePar & GetPar(const long i)
Access all parameters in the order they were inputted.
virtual const string & GetName() const
Name of the object.
long GetNbPar() const
Total number of refinable parameter in the object.
string mName
Name for this RefinableObject. Should be unique, at least in the same scope.+.
output a number as a formatted float:
output a string with a fixed length (adding necessary space or removing excess characters) :