FOX/ObjCryst++  1.10.X (development)
ScatteringData.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 /*
20 * source file for LibCryst++ ScatteringData class
21 *
22 */
23 
24 #include <cmath>
25 
26 #include <typeinfo>
27 
28 #include "cctbx/sgtbx/space_group.h"
29 #include "cctbx/miller/index_generator.h"
30 #include "cctbx/miller/sym_equiv.h"
31 #include "cctbx/eltbx/wavelengths.h"
32 
33 #include "ObjCryst/ObjCryst/ScatteringData.h"
34 #include "ObjCryst/Quirks/VFNDebug.h"
35 #include "ObjCryst/Quirks/VFNStreamFormat.h"
36 #include "ObjCryst/Quirks/Chronometer.h"
37 
38 #ifdef __WX__CRYST__
39  #include "ObjCryst/wxCryst/wxPowderPattern.h"
40  #include "ObjCryst/wxCryst/wxRadiation.h"
41 #endif
42 
43 #include <fstream>
44 #include <iomanip>
45 #include <stdio.h> //for sprintf()
46 
47 #ifdef HAVE_SSE_MATHFUN
48 #include "ObjCryst/Quirks/sse_mathfun.h"
49 #endif
50 
51 #define POSSIBLY_UNUSED(expr) (void)(expr);
52 
53 namespace ObjCryst
54 {
72 
73 const RefParType *gpRefParTypeRadiation=0;
74 const RefParType *gpRefParTypeRadiationWavelength=0;
75 
76 long NiftyStaticGlobalObjectsInitializer_ScatteringData::mCount=0;
77 
78 #ifndef HAVE_SSE_MATHFUN
79 //######################################################################
80 // Tabulated math functions for faster (&less precise) F(hkl) calculation
81 //These function are defined and used in cristallo-spacegroup.cpp
82 //Currently tabulating sine and cosine only
83 //######################################################################
84 
85 //static bool sLibCrystTabulCosineIsInit=false;
86 
87 //conversion value
88 static REAL sLibCrystTabulCosineRatio;
89 // Number of tabulated values of cosine between [0;2pi]
90 // 100 000 is far enough for a model search, yielding a maximum
91 // error less than .05%... 10000 should be enough, too, with (probably) a higher cache hit
92 #define sLibCrystNbTabulSine 8192
93 #define sLibCrystNbTabulSineMASK 8191
94 //storage of tabulated values of cosine, and a table with interlaced csoine/sine values
95 static REAL *spLibCrystTabulCosine;
96 static REAL *spLibCrystTabulCosineSine;
97 
98 void InitLibCrystTabulCosine()
99 {
100  VFN_DEBUG_MESSAGE("InitLibCrystTabulCosine()",10)
101  spLibCrystTabulCosine=new REAL[sLibCrystNbTabulSine];
102  spLibCrystTabulCosineSine=new REAL[sLibCrystNbTabulSine*2];
103  REAL *tmp=spLibCrystTabulCosine;
104  sLibCrystTabulCosineRatio=sLibCrystNbTabulSine/2./M_PI;
105  for(REAL i=0;i<sLibCrystNbTabulSine;i++) *tmp++ = cos(i/sLibCrystTabulCosineRatio);
106  tmp=spLibCrystTabulCosineSine;
107  for(REAL i=0;i<sLibCrystNbTabulSine;i++)
108  {
109  *tmp++ = cos(i/sLibCrystTabulCosineRatio);
110  *tmp++ = sin(i/sLibCrystTabulCosineRatio);
111  }
112 }
113 
114 void DeleteLibCrystTabulCosine()
115 {
116  delete[] spLibCrystTabulCosine;
117  delete[] spLibCrystTabulCosineSine;
118 }
119 
120 // Same for exponential calculations (used for global temperature factors)
121 static bool sLibCrystTabulExpIsInit=false;
122 static const long sLibCrystNbTabulExp=10000;
123 static const REAL sLibCrystMinTabulExp=-5.;
124 static const REAL sLibCrystMaxTabulExp=10.;
125 static REAL *spLibCrystTabulExp;
126 void InitLibCrystTabulExp()
127 {
128  VFN_DEBUG_MESSAGE("InitLibCrystTabulExp()",10)
129  spLibCrystTabulExp=new REAL[sLibCrystNbTabulExp];
130  REAL *tmp=spLibCrystTabulExp;
131  for(REAL i=0;i<sLibCrystNbTabulExp;i++)
132  *tmp++ = exp(sLibCrystMinTabulExp+i*(sLibCrystMaxTabulExp-sLibCrystMinTabulExp)/sLibCrystNbTabulExp);
133  sLibCrystTabulExpIsInit=true;
134 }
135 
136 void DeleteLibCrystTabulExp() { delete[] spLibCrystTabulExp;}
137 
138 //:KLUDGE: The allocated memory for cos and sin table is never freed...
139 // This should be done after the last ScatteringData object is deleted.
140 #endif
141 //
143 // Radiation
144 //
147 mWavelength(1),mXRayTubeName(""),mXRayTubeDeltaLambda(0.),
148 mXRayTubeAlpha2Alpha1Ratio(0.5),mLinearPolarRate(0)
149 {
150  mWavelength=1;
151  this->InitOptions();
152  mRadiationType.SetChoice(RAD_XRAY);
153  mWavelengthType.SetChoice(WAVELENGTH_MONOCHROMATIC);
154  mClockMaster.AddChild(mClockWavelength);
155  mClockMaster.AddChild(mClockRadiation);
156 }
157 
158 Radiation::Radiation(const RadiationType rad,const REAL wavelength)
159 {
160  this->InitOptions();
161  mRadiationType.SetChoice(rad);
162  mWavelengthType.SetChoice(WAVELENGTH_MONOCHROMATIC);
163  mWavelength=wavelength;
164  mXRayTubeName="";
165  mXRayTubeDeltaLambda=0.;//useless here
166  mXRayTubeAlpha2Alpha1Ratio=0.5;//useless here
167  mLinearPolarRate=0.95;//assume it's synchrotron ?
168  mClockMaster.AddChild(mClockWavelength);
169  mClockMaster.AddChild(mClockRadiation);
170 }
171 
172 Radiation::Radiation(const string &XRayTubeElementName,const REAL alpha2Alpha2ratio)
173 {
174  this->InitOptions();
175  this->SetWavelength(XRayTubeElementName,alpha2Alpha2ratio);
176  mClockMaster.AddChild(mClockWavelength);
177  mClockMaster.AddChild(mClockRadiation);
178 }
179 
181 mRadiationType(old.mRadiationType),
182 mWavelengthType(old.mWavelengthType),
183 mWavelength(old.mWavelength),
184 mXRayTubeName(old.mXRayTubeName),
185 mXRayTubeDeltaLambda(old.mXRayTubeDeltaLambda),
186 mXRayTubeAlpha2Alpha1Ratio(old.mXRayTubeAlpha2Alpha1Ratio),
187 mLinearPolarRate(old.mLinearPolarRate)
188 {
189  mClockWavelength.Click();
190  mClockMaster.AddChild(mClockWavelength);
191  mClockMaster.AddChild(mClockRadiation);
192 }
193 
194 Radiation::~Radiation()
195 {}
196 
197 const string& Radiation::GetClassName() const
198 {
199  const static string className="Radiation";
200  return className;
201 }
202 
203 void Radiation::operator=(const Radiation &old)
204 {
211  mClockWavelength.Click();
212  mRadiationType.SetChoice(old.mRadiationType.GetChoice());
213 }
214 
216 {return (RadiationType) mRadiationType.GetChoice();}
217 
219 {
220  mRadiationType.SetChoice(rad);
221  if(rad == RAD_NEUTRON) mLinearPolarRate=0;
222  if(rad == RAD_ELECTRON) mLinearPolarRate=0;
223 }
224 
226 {
227  mWavelengthType.SetChoice((unsigned long) type);
228  if(type==WAVELENGTH_TOF) this->SetRadiationType(RAD_NEUTRON);
229  if(type==WAVELENGTH_ALPHA12) this->SetRadiationType(RAD_XRAY);
230 }
231 
233 {return (WavelengthType) mWavelengthType.GetChoice();}
234 
235 const CrystVector_REAL& Radiation::GetWavelength() const {return mWavelength;}
236 void Radiation::SetWavelength(const REAL l)
237 {
238  mWavelength.resize(1);
239  mWavelength=l;
240  mClockWavelength.Click();
241 }
242 void Radiation::SetWavelength(const string &XRayTubeElementName,
243  const REAL alpha2Alpha2ratio)
244 {
245  VFN_DEBUG_MESSAGE("Radiation::SetWavelength(tubeName,ratio):",5)
246  mXRayTubeName=XRayTubeElementName;
247  this->SetRadiationType(RAD_XRAY);
248  mWavelength.resize(1);
250 
251  if(XRayTubeElementName.length() >=3) //:KLUDGE:
252  {
253  mWavelengthType.SetChoice(WAVELENGTH_MONOCHROMATIC);
254  if(XRayTubeElementName=="CoA1")
255  {
256  mWavelength=1.78901;
257  }
258  else
259  {
260  try
261  {
262  cctbx::eltbx::wavelengths::characteristic ch(mXRayTubeName);
263  if(!ch.is_valid())
264  {
265  cout << "WARNING: could not interpret X-Ray tube name:"<<XRayTubeElementName<<endl
266  << " not modifying wavelength !"<<endl;
267  return;
268  }
269  mWavelength=ch.as_angstrom();
270  }
271  catch(cctbx::error)
272  {
273  cout << "WARNING: could not interpret X-Ray tube name:"<<XRayTubeElementName<<endl
274  << " not modifying wavelength !"<<endl;
275  }
276  }
277  }
278  else
279  {
280  mWavelengthType.SetChoice(WAVELENGTH_ALPHA12);
281  mXRayTubeAlpha2Alpha1Ratio=alpha2Alpha2ratio;
282  REAL lambda1 = 0, lambda2 = 0;
283  if(XRayTubeElementName=="Co")
284  {
285  lambda1=1.78901;
286  lambda2=1.79290;
287  }
288  else
289  {
290  try
291  {
292  cctbx::eltbx::wavelengths::characteristic ch(mXRayTubeName+"A1");
293  if(!ch.is_valid())
294  {
295  cout << "WARNING: could not interpret X-Ray tube name:"<<XRayTubeElementName<<endl
296  << " not modifying wavelength !"<<endl;
297  return;
298  }
299  lambda1=ch.as_angstrom();
300  cctbx::eltbx::wavelengths::characteristic ch2(mXRayTubeName+"A2");
301  if(!ch2.is_valid())
302  {
303  cout << "WARNING: could not interpret X-Ray tube name:"<<XRayTubeElementName<<endl
304  << " not modifying wavelength !"<<endl;
305  return;
306  }
307  lambda2=ch2.as_angstrom();
308  }
309  catch(cctbx::error)
310  {
311  cout << "WARNING: could not interpret X-Ray tube name:"<<XRayTubeElementName<<endl
312  << " not modifying wavelength !"<<endl;
313  }
314  }
315  mXRayTubeDeltaLambda=lambda2-lambda1;
316  mWavelength=lambda1
318  }
319  mClockWavelength.Click();
320 }
321 
323 
325 
326 
327 const RefinableObjClock& Radiation::GetClockWavelength() const {return mClockWavelength;}
329 
330 void Radiation::Print()const
331 {
332  VFN_DEBUG_MESSAGE("Radiation::Print():"<<this->GetName(),5)
333  cout << "Radiation:" << " " ;
334 
335  switch(mRadiationType.GetChoice())
336  {
337  case RAD_NEUTRON: cout<< "Neutron,";break;
338  case RAD_XRAY: cout<< "X-Ray,";break;
339  case RAD_ELECTRON: cout<< "Electron,";break;
340  }
341 
342  cout << "Wavelength=" <<" ";
343  switch(mWavelengthType.GetChoice())
344  {
345  case WAVELENGTH_MONOCHROMATIC: cout<< "monochromatic:"<<" "<<mWavelength(0) <<endl;break;
346  case WAVELENGTH_ALPHA12: cout << "tube:"<<" "<<mXRayTubeName<<", Alpha1/Alpha2= "
347  << mXRayTubeAlpha2Alpha1Ratio<<endl;break;
348  case WAVELENGTH_MAD: cout<< "mad"<<" "<<endl;break;
349  case WAVELENGTH_DAFS: cout<< "dafs"<<" "<<endl;break;
350  case WAVELENGTH_LAUE: cout<< "laue"<<" "<<endl;break;
351  case WAVELENGTH_TOF: cout<< "Time Of Flight"<<" "<<endl;break;
352  }
353 }
354 
355 REAL Radiation::GetLinearPolarRate()const{return mLinearPolarRate;}
356 
357 void Radiation::SetLinearPolarRate(const REAL f){mLinearPolarRate=f;}
358 
359 void Radiation::InitOptions()
360 {
361  static string RadiationTypeName;
362  static string RadiationTypeChoices[3];
363  static string WavelengthTypeName;
364  static string WavelengthTypeChoices[3];
365 
366  static bool needInitNames=true;
367  if(true==needInitNames)
368  {
369  RadiationTypeName="Radiation";
370  RadiationTypeChoices[0]="Neutron";
371  RadiationTypeChoices[1]="X-Ray";
372  RadiationTypeChoices[2]="Electron";
373 
374  WavelengthTypeName="Spectrum";
375  WavelengthTypeChoices[0]="Monochromatic";
376  WavelengthTypeChoices[1]="X-Ray Tube";
377  WavelengthTypeChoices[2]="Time Of Flight";
378  //WavelengthTypeChoices[2]="MAD";
379  //WavelengthTypeChoices[3]="DAFS";
380  //WavelengthTypeChoices[4]="LAUE";
381 
382  needInitNames=false;//Only once for the class
383  }
384  mRadiationType.Init(3,&RadiationTypeName,RadiationTypeChoices);
385  mWavelengthType.Init(3,&WavelengthTypeName,WavelengthTypeChoices);
386  this->AddOption(&mRadiationType);
387  this->AddOption(&mWavelengthType);
388 
389  {//Fixed by default
390  RefinablePar tmp("Wavelength",mWavelength.data(),0.05,20.,
391  gpRefParTypeRadiationWavelength,REFPAR_DERIV_STEP_ABSOLUTE,
392  true,true,true,false,1.0);
393  tmp.SetDerivStep(1e-4);
394  tmp.AssignClock(mClockWavelength);
395  this->AddPar(tmp);
396  }
397 }
398 
399 #ifdef __WX__CRYST__
400 WXCrystObjBasic* Radiation::WXCreate(wxWindow* parent)
401 {
402  //:TODO: Check mpWXCrystObj==0
403  mpWXCrystObj=new WXRadiation(parent,this);
404  return mpWXCrystObj;
405 }
406 #endif
407 
409 //
410 // ScatteringData
411 //
413 
414 ScatteringData::ScatteringData():
415 mNbRefl(0),
416 mpCrystal(0),mGlobalBiso(0),mUseFastLessPreciseFunc(false),
417 mIgnoreImagScattFact(false),mMaxSinThetaOvLambda(10)
418 {
419  VFN_DEBUG_MESSAGE("ScatteringData::ScatteringData()",10)
420  {//This should be done elsewhere...
421  RefinablePar tmp("Global Biso",&mGlobalBiso,-1.,1.,
422  gpRefParTypeScattPowTemperatureIso,REFPAR_DERIV_STEP_ABSOLUTE,
423  true,true,true,false,1.0);
424  tmp.SetDerivStep(1e-4);
425  tmp.AssignClock(mClockGlobalBiso);
426  this->AddPar(tmp);
427  }
428  mClockMaster.AddChild(mClockHKL);
429  mClockMaster.AddChild(mClockGlobalBiso);
430  mClockMaster.AddChild(mClockNbReflUsed);
431  mClockMaster.AddChild(mClockFhklObsSq);
432 }
433 
434 ScatteringData::ScatteringData(const ScatteringData &old):
435 mNbRefl(old.mNbRefl),
436 mpCrystal(old.mpCrystal),mUseFastLessPreciseFunc(old.mUseFastLessPreciseFunc),
437 //Do not copy temporary arrays
438 mClockHKL(old.mClockHKL),
439 mIgnoreImagScattFact(old.mIgnoreImagScattFact),
440 mMaxSinThetaOvLambda(old.mMaxSinThetaOvLambda)
441 {
442  VFN_DEBUG_MESSAGE("ScatteringData::ScatteringData(&old)",10)
443  mClockStructFactor.Reset();
444  mClockTheta.Reset();
445  mClockScattFactor.Reset();
446  mClockScattFactorResonant.Reset();
447  mClockThermicFact.Reset();
448  this->SetHKL(old.GetH(),old.GetK(),old.GetL());
449  VFN_DEBUG_MESSAGE("ScatteringData::ScatteringData(&old):End",5)
450  {//This should be done elsewhere...
451  RefinablePar tmp("Global Biso",&mGlobalBiso,-1.,1.,
452  gpRefParTypeScattPowTemperatureIso,REFPAR_DERIV_STEP_ABSOLUTE,
453  true,true,true,false,1.0);
454  tmp.SetDerivStep(1e-4);
455  tmp.AssignClock(mClockGlobalBiso);
456  this->AddPar(tmp);
457  }
458  mClockMaster.AddChild(mClockHKL);
459  mClockMaster.AddChild(mClockGlobalBiso);
460  mClockMaster.AddChild(mClockNbReflUsed);
461  mClockMaster.AddChild(mClockFhklObsSq);
462 }
463 
464 ScatteringData::~ScatteringData()
465 {
466  VFN_DEBUG_MESSAGE("ScatteringData::~ScatteringData()",10)
467 }
468 
469 void ScatteringData::SetHKL(const CrystVector_REAL &h,
470  const CrystVector_REAL &k,
471  const CrystVector_REAL &l)
472 {
473  VFN_DEBUG_ENTRY("ScatteringData::SetHKL(h,k,l)",5)
474  mNbRefl=h.numElements();
475  mH=h;
476  mK=k;
477  mL=l;
478  mClockHKL.Click();
479  this->PrepareHKLarrays();
480  VFN_DEBUG_EXIT("ScatteringData::SetHKL(h,k,l):End",5)
481 }
482 
483 void ScatteringData::GenHKLFullSpace2(const REAL maxSTOL,const bool unique)
484 {
485  //(*fpObjCrystInformUser)("Generating Full HKL list...");
486  VFN_DEBUG_ENTRY("ScatteringData::GenHKLFullSpace2()",5)
487  TAU_PROFILE("ScatteringData::GenHKLFullSpace2()","void (REAL,bool)",TAU_DEFAULT);
488  if(0==mpCrystal)
489  {
490  throw ObjCrystException("ScatteringData::GenHKLFullSpace2() \
491  no crystal assigned yet to this ScatteringData object.");
492  }
493  cctbx::uctbx::unit_cell uc=cctbx::uctbx::unit_cell(scitbx::af::double6(mpCrystal->GetLatticePar(0),
496  mpCrystal->GetLatticePar(3)*RAD2DEG,
497  mpCrystal->GetLatticePar(4)*RAD2DEG,
498  mpCrystal->GetLatticePar(5)*RAD2DEG));
499  cctbx::miller::index_generator igen(uc,
500  this->GetCrystal().GetSpaceGroup().GetCCTbxSpg().type(),
501  !(this->IsIgnoringImagScattFact()),
502  1/(2*maxSTOL));
503  if(unique)
504  {
505  mNbRefl=0;
506  CrystVector_long H(mNbRefl);
507  CrystVector_long K(mNbRefl);
508  CrystVector_long L(mNbRefl);
509  mMultiplicity.resize(mNbRefl);
510  for(;;)
511  {
512  if(mNbRefl==H.numElements())
513  {
514  H.resizeAndPreserve(mNbRefl+100);
515  K.resizeAndPreserve(mNbRefl+100);
516  L.resizeAndPreserve(mNbRefl+100);
517  mMultiplicity.resizeAndPreserve(mNbRefl+100);
518  }
519  cctbx::miller::index<> h = igen.next();
520  if (h.is_zero()) break;
521  H(mNbRefl)=h[0];
522  K(mNbRefl)=h[1];
523  L(mNbRefl)=h[2];
524  cctbx::miller::sym_equiv_indices sei(this->GetCrystal().GetSpaceGroup().GetCCTbxSpg(),h);
525  mMultiplicity(mNbRefl)=sei.multiplicity(!(this->IsIgnoringImagScattFact()));
526  mNbRefl++;
527  }
528  H.resizeAndPreserve(mNbRefl);
529  K.resizeAndPreserve(mNbRefl);
530  L.resizeAndPreserve(mNbRefl);
531  mMultiplicity.resizeAndPreserve(mNbRefl);
532  this->SetHKL(H,K,L);
533  this->SortReflectionBySinThetaOverLambda(maxSTOL);
534  }
535  else
536  {
537  mNbRefl=0;
538  CrystVector_long H(mNbRefl);
539  CrystVector_long K(mNbRefl);
540  CrystVector_long L(mNbRefl);
541  mMultiplicity.resize(mNbRefl);
542  for(;;)
543  {
544  cctbx::miller::index<> h = igen.next();
545  if (h.is_zero()) break;
546  cctbx::miller::sym_equiv_indices sei(this->GetCrystal().GetSpaceGroup().GetCCTbxSpg(),h);
547  for(int i=0;i<sei.multiplicity(true);i++)
548  {
549  cctbx::miller::index<> k = sei(i).h();
550  if(mNbRefl==H.numElements())
551  {
552  H.resizeAndPreserve(mNbRefl+100);
553  K.resizeAndPreserve(mNbRefl+100);
554  L.resizeAndPreserve(mNbRefl+100);
555  mMultiplicity.resizeAndPreserve(mNbRefl+100);
556  }
557  mMultiplicity(mNbRefl)=sei.multiplicity(!(this->IsIgnoringImagScattFact()));
558  H(mNbRefl)=k[0];
559  K(mNbRefl)=k[1];
560  L(mNbRefl++)=k[2];
561  }
562  }
563  H.resizeAndPreserve(mNbRefl);
564  K.resizeAndPreserve(mNbRefl);
565  L.resizeAndPreserve(mNbRefl);
566  mMultiplicity.resizeAndPreserve(mNbRefl);
567  this->SetHKL(H,K,L);
568  this->SortReflectionBySinThetaOverLambda(maxSTOL);
569  }
570  mClockHKL.Click();
571  /*{
572  char buf [200];
573  sprintf(buf,"Generating Full HKL list...Done (kept %d reflections)",(int)mNbRefl);
574  (*fpObjCrystInformUser)((string)buf);
575  }*/
576  VFN_DEBUG_EXIT("ScatteringData::GenHKLFullSpace2():End",5)
577 }
578 
579 void ScatteringData::GenHKLFullSpace(const REAL maxTheta,const bool useMultiplicity)
580 {
581  VFN_DEBUG_ENTRY("ScatteringData::GenHKLFullSpace()",5)
582  if(this->GetRadiation().GetWavelength()(0) <=.01)
583  {
584  throw ObjCrystException("ScatteringData::GenHKLFullSpace() \
585  no wavelength assigned yet to this ScatteringData object.");;
586  }
587  this->GenHKLFullSpace2(sin(maxTheta)/this->GetRadiation().GetWavelength()(0),useMultiplicity);
588  VFN_DEBUG_EXIT("ScatteringData::GenHKLFullSpace()",5)
589 }
590 
592 
594 {
595  VFN_DEBUG_MESSAGE("ScatteringData::SetCrystal()",5)
596  if(mpCrystal!=0) mpCrystal->DeRegisterClient(*this);
597  mpCrystal=&crystal;
598  this->AddSubRefObj(crystal);
599  crystal.RegisterClient(*this);
603 }
605 
607 
608 bool ScatteringData::HasCrystal()const {return mpCrystal!=0;}
609 
610 long ScatteringData::GetNbRefl() const {return mNbRefl;}
611 
612 const CrystVector_REAL& ScatteringData::GetH() const {return mH;}
613 const CrystVector_REAL& ScatteringData::GetK() const {return mK;}
614 const CrystVector_REAL& ScatteringData::GetL() const {return mL;}
615 
616 const CrystVector_REAL& ScatteringData::GetH2Pi() const {return mH2Pi;}
617 const CrystVector_REAL& ScatteringData::GetK2Pi() const {return mK2Pi;}
618 const CrystVector_REAL& ScatteringData::GetL2Pi() const {return mH2Pi;}
619 
620 const CrystVector_REAL& ScatteringData::GetReflX() const
621 {
622  VFN_DEBUG_ENTRY("ScatteringData::GetReflX()",1)
623  this->CalcSinThetaLambda();
624  VFN_DEBUG_EXIT("ScatteringData::GetReflX()",1)
625  return mX;
626 }
627 const CrystVector_REAL& ScatteringData::GetReflY() const
628 {
629  VFN_DEBUG_ENTRY("ScatteringData::GetReflY()",1)
630  this->CalcSinThetaLambda();
631  VFN_DEBUG_EXIT("ScatteringData::GetReflY()",1)
632  return mY;
633 }
634 const CrystVector_REAL& ScatteringData::GetReflZ() const
635 {
636  VFN_DEBUG_ENTRY("ScatteringData::GetReflZ()",1)
637  this->CalcSinThetaLambda();
638  VFN_DEBUG_EXIT("ScatteringData::GetReflZ()",1)
639  return mZ;
640 }
641 
642 const CrystVector_REAL& ScatteringData::GetSinThetaOverLambda()const
643 {
644  VFN_DEBUG_ENTRY("ScatteringData::GetSinThetaOverLambda()",1)
645  this->CalcSinThetaLambda();
646  VFN_DEBUG_EXIT("ScatteringData::GetSinThetaOverLambda()",1)
647  return mSinThetaLambda;
648 }
649 
650 const CrystVector_REAL& ScatteringData::GetTheta()const
651 {
652  VFN_DEBUG_ENTRY("ScatteringData::GetTheta()",1)
653  this->CalcSinThetaLambda();
654  VFN_DEBUG_EXIT("ScatteringData::GetTheta()",1)
655  return mTheta;
656 }
657 
659 {
660  return mClockTheta;
661 }
662 
663 const CrystVector_REAL& ScatteringData::GetFhklCalcSq() const
664 {
665  VFN_DEBUG_ENTRY("ScatteringData::GetFhklCalcSq()",2)
666  this->CalcStructFactor();
668  #ifdef __LIBCRYST_VECTOR_USE_BLITZ__
669  mFhklCalcSq=pow2(mFhklCalcReal)+pow2(mFhklCalcImag);
670  #else
671  const REAL *pr,*pi;
672  REAL *p;
673  pr=mFhklCalcReal.data();
674  pi=mFhklCalcImag.data();
675  p=mFhklCalcSq.data();
676  for(long i=0;i<mNbReflUsed;i++)
677  {
678  *p++ = *pr * *pr + *pi * *pi;
679  pr++;
680  pi++;
681  }
682  for(long i=mNbReflUsed;i<mNbRefl;i++) *p++ = 0;
683  #endif
685  VFN_DEBUG_EXIT("ScatteringData::GetFhklCalcSq()",2)
686  return mFhklCalcSq;
687 }
688 
689 std::map<RefinablePar*, CrystVector_REAL>& ScatteringData::GetFhklCalcSq_FullDeriv(std::set<RefinablePar *> &vPar)
690 {
691  TAU_PROFILE("ScatteringData::GetFhklCalcSq_FullDeriv()","void ()",TAU_DEFAULT);
692  VFN_DEBUG_ENTRY("ScatteringData::GetFhklCalcSq()",2)
693  this->CalcStructFactor_FullDeriv(vPar);
694  mFhklCalcSq_FullDeriv[0]=this->GetFhklCalcSq();
695  mFhklCalcSq_FullDeriv.clear();// :TODO: avoid complete clear
696  const REAL *pr,*pi,*prd,*pid;
697  REAL *p;
698  for(std::set<RefinablePar *>::iterator par=vPar.begin();par!=vPar.end();par++)
699  {
700  if((*par)==0) continue;
701  if(mFhklCalcReal_FullDeriv[*par].size()==0)
702  {
703  mFhklCalcSq_FullDeriv[*par].resize(0);
704  continue;
705  }
706  mFhklCalcSq_FullDeriv[*par].resize(mNbRefl);//Should use mNbRefleUsed instead ?
707  pr=mFhklCalcReal.data();
708  pi=mFhklCalcImag.data();
709  prd=mFhklCalcReal_FullDeriv[*par].data();
710  pid=mFhklCalcImag_FullDeriv[*par].data();
711  p=mFhklCalcSq_FullDeriv[*par].data();
712  for(long i=0;i<mNbReflUsed;i++)
713  *p++ = 2*(*pr++ * *prd++ + *pi++ * *pid++);
714  for(long i=mNbReflUsed;i<mNbRefl;i++) *p++ = 0;
715  }
716  VFN_DEBUG_EXIT("ScatteringData::GetFhklCalcSq()",2)
717  return mFhklCalcSq_FullDeriv;
718 }
719 
720 const CrystVector_REAL& ScatteringData::GetFhklCalcReal() const
721 {
722  VFN_DEBUG_ENTRY("ScatteringData::GetFhklCalcReal()",2)
723  this->CalcStructFactor();
724  VFN_DEBUG_EXIT("ScatteringData::GetFhklCalcReal()",2)
725  return mFhklCalcReal;
726 }
727 
728 const CrystVector_REAL& ScatteringData::GetFhklCalcImag() const
729 {
730  VFN_DEBUG_ENTRY("ScatteringData::GetFhklCalcImag()",2)
731  this->CalcStructFactor();
732  VFN_DEBUG_EXIT("ScatteringData::GetFhklCalcImag()",2)
733  return mFhklCalcImag;
734 }
735 
736 const CrystVector_REAL& ScatteringData::GetFhklObsSq() const
737 {
738  return mFhklObsSq;
739 }
740 
741 const map<const ScatteringPower*,CrystVector_REAL>& ScatteringData::GetScatteringFactor() const
742 {
743  this->CalcScattFactor();
744  return mvScatteringFactor;
745 }
746 
747 CrystVector_REAL ScatteringData::GetWavelength()const {return this->GetRadiation().GetWavelength();}
748 #if 0
749 void ScatteringData::SetUseFastLessPreciseFunc(const bool useItOrNot)
750 {
751  mUseFastLessPreciseFunc=useItOrNot;
754 }
755 #endif
757 {
758  VFN_DEBUG_MESSAGE("ScatteringData::SetIsIgnoringImagScattFact():"<<b,10)
762 }
764 
765 void ScatteringData::PrintFhklCalc(ostream &os)const
766 {
767  VFN_DEBUG_ENTRY("ScatteringData::PrintFhklCalc()",5)
768  this->GetFhklCalcSq();
769  CrystVector_REAL theta;
770  theta=mTheta;
771  theta *= RAD2DEG;
772  os <<" Number of reflections:"<<mNbRefl<<endl;
773  os <<" H K L F(hkl)^2 Re(F) Im(F)";
774  os <<" Theta 1/2d"<<endl;
775  os << FormatVertVectorHKLFloats<REAL>
776  (mH,mK,mL,mFhklCalcSq,mFhklCalcReal,mFhklCalcImag,theta,mSinThetaLambda,12,4,mNbReflUsed);
777  VFN_DEBUG_EXIT("ScatteringData::PrintFhklCalc()",5)
778 }
779 
781 {
782  VFN_DEBUG_ENTRY("ScatteringData::PrintFhklCalcDetail()",5)
783  this->GetFhklCalcSq();
784  CrystVector_REAL theta;
785  theta=mTheta;
786  theta *= RAD2DEG;
787  vector<const CrystVector_REAL *> v;
788  v.push_back(&mH);
789  v.push_back(&mK);
790  v.push_back(&mL);
791  v.push_back(&mSinThetaLambda);
792  v.push_back(&theta);
793  v.push_back(&mFhklCalcSq);
794  v.push_back(&mFhklCalcReal);
795  v.push_back(&mFhklCalcImag);
796  os <<" Number of reflections:"<<mNbRefl<<endl;
797  os <<" H K L 1/2d Theta F(hkl)^2";
798  os <<" Re(F) Im(F) ";
799  vector<CrystVector_REAL> sf;
800  sf.resize(mvRealGeomSF.size()*2);
801  long i=0;
802  for(map<const ScatteringPower*,CrystVector_REAL>::const_iterator
803  pos=mvRealGeomSF.begin();pos!=mvRealGeomSF.end();++pos)
804  {
805  os << FormatString("Re(F)_"+pos->first->GetName(),14)
806  << FormatString("Im(F)_"+pos->first->GetName(),14);
807  cout<<pos->first->GetName()<<":"<<pos->first->GetForwardScatteringFactor(RAD_XRAY)<<endl;
808  sf[2*i] = mvRealGeomSF[pos->first];
809  sf[2*i] *= mvScatteringFactor[pos->first];
810  sf[2*i] *= mvTemperatureFactor[pos->first];
811  sf[2*i+1] = mvImagGeomSF[pos->first];
812  sf[2*i+1] *= mvScatteringFactor[pos->first];
813  sf[2*i+1] *= mvTemperatureFactor[pos->first];
814  v.push_back(&(sf[2*i]));
815  v.push_back(&(sf[2*i+1]));
816  //v.push_back(mvRealGeomSF[pos->first]);
817  //v.push_back(mvImagGeomSF[pos->first]);
818  //v.push_back(mvScatteringFactor[pos->first]);
819  i++;
820  }
821  os<<endl;
822  os << FormatVertVectorHKLFloats<REAL>(v,12,4,mNbReflUsed);
823  VFN_DEBUG_EXIT("ScatteringData::PrintFhklCalcDetail()",5)
824 }
825 
826 void ScatteringData::BeginOptimization(const bool allowApproximations,
827  const bool enableRestraints)
828 {
829  if(mUseFastLessPreciseFunc!=allowApproximations)
830  {
834  }
835  mUseFastLessPreciseFunc=allowApproximations;
836  this->RefinableObj::BeginOptimization(allowApproximations,enableRestraints);
837 }
839 {
840  if(mOptimizationDepth==1)
841  {
842  if(mUseFastLessPreciseFunc==true)
843  {
847  }
849  }
851 }
852 
854 {
855  if(mUseFastLessPreciseFunc!=allow)
856  {
860  }
863 }
864 
866 {
867  VFN_DEBUG_ENTRY("ScatteringData::PrepareHKLarrays()"<<mNbRefl<<" reflections",5)
868  mFhklCalcReal.resize(mNbRefl);
869  mFhklCalcImag.resize(mNbRefl);
870  mFhklCalcSq.resize(mNbRefl);
871  mFhklCalcReal=0;
872  mFhklCalcImag=0;
873  mFhklCalcSq=0;
874 
875  mIntH=mH;
876  mIntK=mK;
877  mIntL=mL;
878 
879  mH2Pi=mH;
880  mK2Pi=mK;
881  mL2Pi=mL;
882  mH2Pi*=(2*M_PI);
883  mK2Pi*=(2*M_PI);
884  mL2Pi*=(2*M_PI);
885 
886  // If we extracted some intensities, try to keep them.
887  // Do not do this if the number of reflections changed too much, if there is no crystal
888  // structure associated, or if the spacegroup changed.
889  bool noSpgChange=false;
891  if( (mFhklObsSq.numElements()>0) && (abs(mFhklObsSq.numElements()-mNbRefl)<(0.1*mNbRefl)) && (noSpgChange) )
892  mFhklObsSq.resizeAndPreserve(mNbRefl);
893  else mFhklObsSq.resize(0);
895 
896  mNbReflUsed=mNbRefl;
897 
899  for(long i=0;i<mNbRefl;i++)
900  {
903  }
904  /*
905  {
906  mpCrystal->GetSpaceGroup().Print();
907  for(long i=0;i<mNbRefl;i++)
908  {
909  cout<<mIntH(i)<<" "<<mIntK(i)<<" "<<mIntL(i)<<" "<<mExpectedIntensityFactor(i)<<endl;
910  }
911  }
912  */
913 
914  mClockHKL.Click();
915  VFN_DEBUG_EXIT("ScatteringData::PrepareHKLarrays()"<<mNbRefl<<" reflections",5)
916 }
917 
921 {
922  if(this->IsBeingRefined()) return mNbReflUsed;
923  VFN_DEBUG_MESSAGE("ScatteringData::GetNbReflBelowMaxSinThetaOvLambda()",4)
924  this->CalcSinThetaLambda();
925  if((mNbReflUsed>0)&&(mNbReflUsed<mNbRefl))
926  {
927  if( (mSinThetaLambda(mNbReflUsed )>mMaxSinThetaOvLambda)
928  &&(mSinThetaLambda(mNbReflUsed-1)<=mMaxSinThetaOvLambda)) return mNbReflUsed;
929  }
930 
931  if((mNbReflUsed==mNbRefl)&&(mSinThetaLambda(mNbRefl-1)<=mMaxSinThetaOvLambda))
932  return mNbReflUsed;
933  long i;
934  for(i=0;i<mNbRefl;i++) if(mSinThetaLambda(i)>mMaxSinThetaOvLambda) break;
935  if(i!=mNbReflUsed)
936  {
937  mNbReflUsed=i;
939  VFN_DEBUG_MESSAGE("->Changed Max sin(theta)/lambda="<<mMaxSinThetaOvLambda\
940  <<" nb refl="<<mNbReflUsed,4)
941  }
942  return mNbReflUsed;
943 }
945 {return mClockNbReflUsed;}
946 
947 CrystVector_long ScatteringData::SortReflectionBySinThetaOverLambda(const REAL maxSTOL)
948 {
949  TAU_PROFILE("ScatteringData::SortReflectionBySinThetaOverLambda()","void ()",TAU_DEFAULT);
950  VFN_DEBUG_ENTRY("ScatteringData::SortReflectionBySinThetaOverLambda()",5)
951  this->CalcSinThetaLambda();
952  CrystVector_long sortedSubs;
953  sortedSubs=SortSubs(mSinThetaLambda);
954  CrystVector_long oldH,oldK,oldL,oldMult;
955  oldH=mH;
956  oldK=mK;
957  oldL=mL;
958  oldMult=mMultiplicity;
959  long subs;
960  long shift=0;
961 
962  //get rid of [0,0,0] reflection
963  VFN_DEBUG_MESSAGE("ScatteringData::SortReflectionBySinThetaOverLambda() 1",2)
964  if(0==mSinThetaLambda(sortedSubs(0)))
965  {
966  shift=1;
967  mNbRefl -= 1;
968  mH.resize(mNbRefl);
969  mK.resize(mNbRefl);
970  mL.resize(mNbRefl);
971  mMultiplicity.resize(mNbRefl);
972  }
973  VFN_DEBUG_MESSAGE("ScatteringData::SortReflectionBySinThetaOverLambda() 2",2)
974  for(long i=0;i<mNbRefl;i++)
975  {
976  subs=sortedSubs(i+shift);
977  mH(i)=oldH(subs);
978  mK(i)=oldK(subs);
979  mL(i)=oldL(subs);
980  mMultiplicity(i)=oldMult(subs);
981  }
982  mClockHKL.Click();
983  VFN_DEBUG_MESSAGE("ScatteringData::SortReflectionBySinThetaOverLambda() 3",2)
984  this->PrepareHKLarrays();
985  this->CalcSinThetaLambda();
986 
987  VFN_DEBUG_MESSAGE("ScatteringData::SortReflectionBySinThetaOverLambda() 4",2)
988  if(0<maxSTOL)
989  {
990  VFN_DEBUG_MESSAGE("ScatteringData::SortReflectionBySinThetaOverLambda() 5"<<maxSTOL,2)
991  long maxSubs=0;
992  VFN_DEBUG_MESSAGE(" "<< mIntH(maxSubs)<<" "<< mIntK(maxSubs)<<" "<< mIntL(maxSubs)<<" "<<mSinThetaLambda(maxSubs),1)
993  for(maxSubs=0;mSinThetaLambda(maxSubs)<maxSTOL;maxSubs++)
994  {
995  VFN_DEBUG_MESSAGE(" "<< mIntH(maxSubs)<<" "<< mIntK(maxSubs)<<" "<< mIntL(maxSubs)<<" "<<mSinThetaLambda(maxSubs),1)
996  if(maxSubs==(mNbRefl-1))
997  {
998  maxSubs=mNbRefl;
999  break;
1000  }
1001  }
1002  if(maxSubs==mNbRefl)
1003  {
1004  VFN_DEBUG_EXIT("ScatteringData::SortReflectionBySinThetaOverLambda():"<<mNbRefl<<" reflections",5)
1005  return sortedSubs;
1006  }
1007  mNbRefl=maxSubs;
1008  mH.resizeAndPreserve(mNbRefl);
1009  mK.resizeAndPreserve(mNbRefl);
1010  mL.resizeAndPreserve(mNbRefl);
1011  mMultiplicity.resizeAndPreserve(mNbRefl);
1012  sortedSubs.resizeAndPreserve(mNbRefl);
1013  mClockHKL.Click();
1014  this->PrepareHKLarrays();
1015  }
1016  VFN_DEBUG_EXIT("ScatteringData::SortReflectionBySinThetaOverLambda():"<<mNbRefl<<" reflections",5)
1017  return sortedSubs;
1018 }
1019 
1021 {
1022  TAU_PROFILE("ScatteringData::EliminateExtinctReflections()","void ()",TAU_DEFAULT);
1023  VFN_DEBUG_ENTRY("ScatteringData::EliminateExtinctReflections()",7)
1024 
1025  long nbKeptRefl=0;
1026  CrystVector_long subscriptKeptRefl(mNbRefl);
1027  subscriptKeptRefl=0;
1028  for(long j=0;j<mNbRefl;j++)
1029  {
1030  if( this->GetCrystal().GetSpaceGroup().IsReflSystematicAbsent(mH(j),mK(j),mL(j))==false )
1031  subscriptKeptRefl(nbKeptRefl++)=j;
1032  }
1033  VFN_DEBUG_MESSAGE("ScatteringData::EliminateExtinctReflections():4",5)
1034  //Keep only the elected reflections
1035  mNbRefl=nbKeptRefl;
1036  {
1037  CrystVector_long oldH,oldK,oldL;
1038  CrystVector_int oldMulti;
1039  long subs;
1040 
1041  oldH=mH;
1042  oldK=mK;
1043  oldL=mL;
1044  oldMulti=mMultiplicity;
1045 
1046  mMultiplicity.resize(mNbRefl);
1047  mH.resize(mNbRefl);
1048  mK.resize(mNbRefl);
1049  mL.resize(mNbRefl);
1050  for(long i=0;i<mNbRefl;i++)
1051  {
1052  subs=subscriptKeptRefl(i);
1053  mH(i)=oldH(subs);
1054  mK(i)=oldK(subs);
1055  mL(i)=oldL(subs);
1056  mMultiplicity(i)=oldMulti(subs);
1057  }
1058  }
1059  this->PrepareHKLarrays();
1060  VFN_DEBUG_EXIT("ScatteringData::EliminateExtinctReflections():End",7)
1061  return subscriptKeptRefl;
1062 }
1063 
1065 {
1066  if(mClockTheta>mClockMaster) return;
1067  if( 0 == mpCrystal) throw ObjCrystException("ScatteringData::CalcSinThetaLambda() \
1068  Cannot compute sin(theta)/lambda : there is no crystal affected to this \
1069  ScatteringData object yet.");
1070 
1071  if( 0 == this->GetNbRefl()) throw ObjCrystException("ScatteringData::CalcSinThetaLambda() \
1072  Cannot compute sin(theta)/lambda : there are no reflections !");
1073 
1074  if( (mClockTheta>this->GetRadiation().GetClockWavelength())
1078 
1079  VFN_DEBUG_ENTRY("ScatteringData::CalcSinThetaLambda()",3)
1080  TAU_PROFILE("ScatteringData::CalcSinThetaLambda()","void (bool)",TAU_DEFAULT);
1081  mSinThetaLambda.resize(mNbRefl);
1082 
1083  const CrystMatrix_REAL bMatrix= this->GetBMatrix();
1084  mX.resize(this->GetNbRefl());
1085  mY.resize(this->GetNbRefl());
1086  mZ.resize(this->GetNbRefl());
1087  for(int i=0;i<this->GetNbRefl();i++)
1088  { //:TODO: faster,nicer
1089  mX(i)=bMatrix(0,0)*mH(i)+bMatrix(0,1)*mK(i)+bMatrix(0,2)*mL(i);
1090  mY(i)=bMatrix(1,0)*mH(i)+bMatrix(1,1)*mK(i)+bMatrix(1,2)*mL(i);
1091  mZ(i)=bMatrix(2,0)*mH(i)+bMatrix(2,1)*mK(i)+bMatrix(2,2)*mL(i);
1092  }
1093  //cout << bMatrix << endl << xyz<<endl;
1094  for(int i=0;i< (this->GetNbRefl());i++)
1095  mSinThetaLambda(i)=sqrt(pow(mX(i),2)+pow(mY(i),2)+pow(mZ(i),2))/2;
1096 
1097  #if 0
1098  // Direct calculation from a,b,c,alpha,beta,gamma
1099  const REAL a=mpCrystal->GetLatticePar(0);
1100  const REAL b=mpCrystal->GetLatticePar(1);
1101  const REAL c=mpCrystal->GetLatticePar(2);
1102  const REAL ca=cos(mpCrystal->GetLatticePar(3));
1103  const REAL sa=sin(mpCrystal->GetLatticePar(3));
1104  const REAL cb=cos(mpCrystal->GetLatticePar(4));
1105  const REAL sb=sin(mpCrystal->GetLatticePar(4));
1106  const REAL cg=cos(mpCrystal->GetLatticePar(5));
1107  const REAL sg=sin(mpCrystal->GetLatticePar(5));
1108  for(int i=0;i< (this->GetNbRefl());i++)
1109  {
1110  const REAL h=mH(i),k=mK(i),l=mL(i);
1111  mSinThetaLambda(i)=0.5*sqrt((h*h/(a*a)*sa*sa+k*k/(b*b)*sb*sb+l*l/(c*c)*sg*sg+2*k*l/(b*c)*(cb*cg-ca)+2*l*h/(c*a)*(cg*ca-cb)+2*h*k/(a*b)*(ca*cb-cg))/(1-ca*ca-cb*cb-cg*cg+2*ca*cb*cg));
1112  }
1113  #endif
1114 
1115 
1116  if(this->GetRadiation().GetWavelengthType()!=WAVELENGTH_TOF)
1117  {
1118  if(this->GetRadiation().GetWavelength()(0) > 0)
1119  {
1120  mTheta.resize(mNbRefl);
1121  for(int i=0;i< (this->GetNbRefl());i++)
1122  {
1123  if( (mSinThetaLambda(i)*this->GetRadiation().GetWavelength()(0))>1)
1124  {
1125  //:KLUDGE: :TODO:
1126  mTheta(i)=M_PI;
1127  /*
1128  ofstream out("log.txt");
1129  out << "Error when computing Sin(theta) :"
1130  << "i="<<i<<" ,mSinThetaLambda(i)="<<mSinThetaLambda(i)
1131  << " ,this->GetRadiation().GetWavelength()(0)="
1132  << this->GetRadiation().GetWavelength()(0)
1133  << " ,H="<<mH(i)
1134  << " ,K="<<mK(i)
1135  << " ,L="<<mL(i)
1136  <<endl;
1137  out.close();
1138  abort();
1139  */
1140  }
1141  else
1142  {
1143  mTheta(i)=asin(mSinThetaLambda(i)*this->GetRadiation().GetWavelength()(0));
1144  }
1145  }
1146  } else
1147  {
1148  cout << "Wavelength not given in ScatteringData::CalcSinThetaLambda() !" <<endl;
1149  throw 0;
1150  }
1151  }
1152  else mTheta.resize(0);
1153 
1154  mClockTheta.Click();
1155  VFN_DEBUG_EXIT("ScatteringData::CalcSinThetaLambda()",3)
1156 }
1157 
1158 REAL ScatteringData::CalcSinThetaLambda(REAL h, REAL k, REAL l)const
1159 {
1160  const CrystMatrix_REAL bMatrix= this->GetBMatrix();
1161  const REAL x=bMatrix(0,0)*h+bMatrix(0,1)*k+bMatrix(0,2)*l;
1162  const REAL y=bMatrix(1,0)*h+bMatrix(1,1)*k+bMatrix(1,2)*l;
1163  const REAL z=bMatrix(2,0)*h+bMatrix(2,1)*k+bMatrix(2,2)*l;
1164  return sqrt(x*x+y*y+z*z)/2;
1165 }
1166 
1167 const CrystMatrix_REAL& ScatteringData::GetBMatrix() const
1168 {
1169  return this->GetCrystal().GetBMatrix();
1170 }
1171 
1173 {
1174  //if(mClockScattFactor>mClockMaster) return;
1175  if( (mClockScattFactor>this->GetRadiation().GetClockWavelength())
1180  TAU_PROFILE("ScatteringData::CalcScattFactor()","void (bool)",TAU_DEFAULT);
1181  VFN_DEBUG_ENTRY("ScatteringData::CalcScattFactor()",4)
1182  this->CalcResonantScattFactor();
1183  mvScatteringFactor.clear();
1184  for(int i=mpCrystal->GetScatteringPowerRegistry().GetNb()-1;i>=0;i--)
1185  {
1186  const ScatteringPower *pScattPow=&(mpCrystal->GetScatteringPowerRegistry().GetObj(i));
1187  mvScatteringFactor[pScattPow]=pScattPow->GetScatteringFactor(*this);
1188  //Directly add Fprime
1189  mvScatteringFactor[pScattPow]+= this->mvFprime[pScattPow];
1190  VFN_DEBUG_MESSAGE("-> H K L sin(t/l) f0+f'"
1192  mvScatteringFactor[pScattPow],10,4,mNbReflUsed),1);
1193  }
1195  VFN_DEBUG_EXIT("ScatteringData::CalcScattFactor()",4)
1196 }
1197 
1199 {
1200  //if(mClockThermicFact>mClockMaster) return;
1201  if( (mClockThermicFact>this->GetRadiation().GetClockWavelength())
1206  TAU_PROFILE("ScatteringData::CalcTemperatureFactor()","void (bool)",TAU_DEFAULT);
1207  VFN_DEBUG_ENTRY("ScatteringData::CalcTemperatureFactor()",4)
1208  mvTemperatureFactor.clear();
1209  for(int i=mpCrystal->GetScatteringPowerRegistry().GetNb()-1;i>=0;i--)
1210  {
1211  const ScatteringPower *pScattPow=&(mpCrystal->GetScatteringPowerRegistry().GetObj(i));
1212  mvTemperatureFactor[pScattPow]=pScattPow->GetTemperatureFactor(*this);
1213  VFN_DEBUG_MESSAGE("-> H K L sin(t/l) DebyeWaller"<<endl
1215  mvTemperatureFactor[pScattPow],10,4,mNbReflUsed),1);
1216  }
1218  VFN_DEBUG_EXIT("ScatteringData::CalcTemperatureFactor()",4)
1219 }
1220 
1222 {
1225  VFN_DEBUG_ENTRY("ScatteringData::CalcResonantScattFactor()",4)
1226  TAU_PROFILE("ScatteringData::CalcResonantScattFactor()","void (bool)",TAU_DEFAULT);
1227 
1228  mvFprime.clear();
1229  mvFsecond.clear();
1230  if(this->GetRadiation().GetWavelength()(0) == 0)
1231  {
1232  VFN_DEBUG_EXIT("ScatteringData::CalcResonantScattFactor()->Lambda=0. fprime=fsecond=0",4)
1233  return;
1234  }
1235  else
1236  {
1237  for(int i=mpCrystal->GetScatteringPowerRegistry().GetNb()-1;i>=0;i--)
1238  {
1239  const ScatteringPower *pScattPow=&(mpCrystal->GetScatteringPowerRegistry().GetObj(i));
1240  mvFprime [pScattPow]=pScattPow->GetResonantScattFactReal(*this)(0);
1241  mvFsecond[pScattPow]=pScattPow->GetResonantScattFactImag(*this)(0);
1242  }
1243  }
1245  VFN_DEBUG_EXIT("ScatteringData::CalcResonantScattFactor()",4)
1246 }
1247 
1249 {
1250  this->GetNbReflBelowMaxSinThetaOvLambda();//update mNbReflUsed, also recalc sin(theta)/lambda
1256  VFN_DEBUG_MESSAGE("ScatteringData::CalcGlobalTemperatureFactor()",2)
1257  TAU_PROFILE("ScatteringData::CalcGlobalTemperatureFactor()","void ()",TAU_DEFAULT);
1258 
1260  //if(true==mUseFastLessPreciseFunc) //:TODO:
1261  {
1262  }
1263  //else
1264  {
1265  const REAL *stol=this->GetSinThetaOverLambda().data();
1266  REAL *fact=mGlobalTemperatureFactor.data();
1267  for(long i=0;i<mNbReflUsed;i++) {*fact++ = exp(-mGlobalBiso * *stol * *stol);stol++;}
1268  }
1270 }
1271 
1273 {
1274  this->GetNbReflBelowMaxSinThetaOvLambda();//check mNbReflUsed, also recalc sin(theta)/lambda
1275  if(mClockStructFactor>mClockMaster) return;
1276 
1277  //:TODO: Anisotropic Thermic factors
1278  //TAU_PROFILE_TIMER(timer1,"ScatteringData::CalcStructFactor1:Prepare","", TAU_FIELD);
1279  //TAU_PROFILE_TIMER(timer2,"ScatteringData::CalcStructFactor2:GeomStructFact","", TAU_FIELD);
1280  //TAU_PROFILE_TIMER(timer3,"ScatteringData::CalcStructFactor3:Scatt.Factors","", TAU_FIELD);
1281  //TAU_PROFILE_TIMER(timer4,"ScatteringData::CalcStructFactor4:Finish,DynCorr","", TAU_FIELD);
1282 
1283  //TAU_PROFILE_START(timer1);
1284 
1285  const long nbRefl=this->GetNbRefl();
1286  this->CalcSinThetaLambda();
1287  //TAU_PROFILE_STOP(timer1);
1288  //TAU_PROFILE_START(timer2);
1289  this->CalcGeomStructFactor();
1290  //TAU_PROFILE_STOP(timer2);
1291  //TAU_PROFILE_START(timer3);
1292  this->CalcScattFactor();
1293  this->CalcResonantScattFactor();
1294  this->CalcTemperatureFactor();
1296  this->CalcLuzzatiFactor();
1297  this->CalcStructFactVariance();
1298  //TAU_PROFILE_STOP(timer3);
1299 
1300  //OK, really must recompute SFs?
1301  VFN_DEBUG_MESSAGE("ScatteringData::CalcStructFactor():Fhkl Recalc ?"<<endl
1302  <<"mClockStructFactor<mClockGlobalTemperatureFact"<<(bool)(mClockStructFactor<mClockGlobalTemperatureFact)<<endl
1303  <<"mClockStructFactor<mClockGeomStructFact" <<(bool)(mClockStructFactor<mClockGeomStructFact)<<endl
1304  <<"mClockStructFactor<mClockScattFactorResonant" <<(bool)(mClockStructFactor<mClockScattFactorResonant)<<endl
1305  <<"mClockStructFactor<mClockThermicFact" <<(bool)(mClockStructFactor<mClockThermicFact)<<endl
1306  <<"mClockStructFactor<mClockLuzzatiFactor" <<(bool)(mClockStructFactor<mClockLuzzatiFactor)<<endl
1307  ,2)
1312  &&(mClockStructFactor>mClockFhklCalcVariance)
1313  &&(mClockStructFactor>mClockLuzzatiFactor)) return;
1314  VFN_DEBUG_ENTRY("ScatteringData::CalcStructFactor()",3)
1315  TAU_PROFILE("ScatteringData::CalcStructFactor()","void ()",TAU_DEFAULT);
1316  //TAU_PROFILE_START(timer4);
1317  //reset Fcalc
1318  mFhklCalcReal.resize(nbRefl);
1319  mFhklCalcImag.resize(nbRefl);
1320  mFhklCalcReal=0;
1321  mFhklCalcImag=0;
1322  //Add all contributions
1323  for(map<const ScatteringPower*,CrystVector_REAL>::const_iterator pos=mvRealGeomSF.begin();
1324  pos!=mvRealGeomSF.end();++pos)
1325  {
1326  const ScatteringPower* pScattPow=pos->first;
1327  VFN_DEBUG_MESSAGE("ScatteringData::CalcStructFactor():Fhkl Recalc, "<<pScattPow->GetName(),2)
1328  const REAL * RESTRICT pGeomR=mvRealGeomSF[pScattPow].data();
1329  const REAL * RESTRICT pGeomI=mvImagGeomSF[pScattPow].data();
1330  const REAL * RESTRICT pScatt=mvScatteringFactor[pScattPow].data();
1331  const REAL * RESTRICT pTemp=mvTemperatureFactor[pScattPow].data();
1332 
1333  REAL * RESTRICT pReal=mFhklCalcReal.data();
1334  REAL * RESTRICT pImag=mFhklCalcImag.data();
1335 
1336  VFN_DEBUG_MESSAGE("->mvRealGeomSF[i] "
1337  <<mvRealGeomSF[pScattPow].numElements()<<"elements",2)
1338  VFN_DEBUG_MESSAGE("->mvImagGeomSF[i] "
1339  <<mvImagGeomSF[pScattPow].numElements()<<"elements",2)
1340  VFN_DEBUG_MESSAGE("->mvScatteringFactor[i]"
1341  <<mvScatteringFactor[pScattPow].numElements()<<"elements",1)
1342  VFN_DEBUG_MESSAGE("->mvTemperatureFactor[i]"
1343  <<mvTemperatureFactor[pScattPow].numElements()<<"elements",1)
1344  VFN_DEBUG_MESSAGE("->mFhklCalcReal "<<mFhklCalcReal.numElements()<<"elements",2)
1345  VFN_DEBUG_MESSAGE("->mFhklCalcImag "<<mFhklCalcImag.numElements()<<"elements",2)
1346  VFN_DEBUG_MESSAGE("-> H K L sin(t/l) Re(F) Im(F) scatt Temp->"<<pScattPow->GetName(),1)
1347 
1348  VFN_DEBUG_MESSAGE(FormatVertVectorHKLFloats<REAL>(mH,mK,mL,mSinThetaLambda,
1349  mvRealGeomSF[pScattPow],
1350  mvImagGeomSF[pScattPow],
1351  mvScatteringFactor[pScattPow],
1352  mvTemperatureFactor[pScattPow],10,4,mNbReflUsed
1353  ),1);
1354  if(mvLuzzatiFactor[pScattPow].numElements()>0)
1355  {// using maximum likelihood
1356  const REAL* RESTRICT pLuzzati=mvLuzzatiFactor[pScattPow].data();
1357  if(false==mIgnoreImagScattFact)
1358  {
1359  const REAL fsecond=mvFsecond[pScattPow];
1360  VFN_DEBUG_MESSAGE("->fsecond= "<<fsecond,10)
1361  for(long j=mNbReflUsed;j>0;j--)
1362  {
1363  VFN_DEBUG_MESSAGE("-->"<<j<<" "<<*pReal<<" "<<*pImag<<" "<<*pGeomR<<" "<<*pGeomI<<" "<<*pScatt<<" "<<*pTemp,1)
1364  *pReal++ += (*pGeomR * *pScatt - *pGeomI * fsecond)* *pTemp * *pLuzzati;
1365  *pImag++ += (*pGeomI++ * *pScatt++ + *pGeomR++ * fsecond)* *pTemp++ * *pLuzzati++;
1366  }
1367  }
1368  else
1369  {
1370  for(long j=mNbReflUsed;j>0;j--)
1371  {
1372  *pReal++ += *pGeomR++ * *pTemp * *pScatt * *pLuzzati;
1373  *pImag++ += *pGeomI++ * *pTemp++ * *pScatt++ * *pLuzzati++;
1374  }
1375  }
1376  VFN_DEBUG_MESSAGE("ScatteringData::CalcStructFactor():"<<mIgnoreImagScattFact
1377  <<",f\"="<<mvFsecond[pScattPow]<<endl<<
1379  mvRealGeomSF[pScattPow],
1380  mvImagGeomSF[pScattPow],
1381  mvScatteringFactor[pScattPow],
1382  mvTemperatureFactor[pScattPow],
1383  mvLuzzatiFactor[pScattPow],
1384  mFhklCalcReal,
1385  mFhklCalcImag,10,4,mNbReflUsed
1386  ),2);
1387  }
1388  else
1389  {
1390  if(false==mIgnoreImagScattFact)
1391  {
1392  const REAL fsecond=mvFsecond[pScattPow];
1393  VFN_DEBUG_MESSAGE("->fsecond= "<<fsecond,2)
1394  for(long j=mNbReflUsed;j>0;j--)
1395  {
1396  *pReal += (*pGeomR * *pScatt - *pGeomI * fsecond)* *pTemp;
1397  *pImag += (*pGeomI * *pScatt + *pGeomR * fsecond)* *pTemp;
1398  VFN_DEBUG_MESSAGE("-->"<<j<<" "<<*pReal<<" "<<*pImag<<" "<<*pGeomR<<" "<<*pGeomI<<" "<<*pScatt<<" "<<*pTemp,1)
1399  pGeomR++;pGeomI++;pTemp++;pScatt++;pReal++;pImag++;
1400  }
1401  }
1402  else
1403  {
1404  for(long j=mNbReflUsed;j>0;j--)
1405  {
1406  *pReal++ += *pGeomR++ * *pTemp * *pScatt;
1407  *pImag++ += *pGeomI++ * *pTemp++ * *pScatt++;
1408  }
1409  }
1410  VFN_DEBUG_MESSAGE(FormatVertVectorHKLFloats<REAL>(mH,mK,mL,mSinThetaLambda,
1411  mvRealGeomSF[pScattPow],
1412  mvImagGeomSF[pScattPow],
1413  mvScatteringFactor[pScattPow],
1414  mvTemperatureFactor[pScattPow],
1415  mFhklCalcReal,
1416  mFhklCalcImag,10,4,mNbReflUsed
1417  ),2);
1418  }
1419  }
1420  //TAU_PROFILE_STOP(timer4);
1421  {
1422  //this->CalcGlobalTemperatureFactor();
1423  if(mGlobalTemperatureFactor.numElements()>0)
1424  {//else for some reason it's useless
1425  REAL *pReal=mFhklCalcReal.data();
1426  REAL *pImag=mFhklCalcImag.data();
1427  const REAL *pTemp=mGlobalTemperatureFactor.data();
1428  for(long j=0;j<mNbReflUsed;j++)
1429  {
1430  *pReal++ *= *pTemp;
1431  *pImag++ *= *pTemp++;
1432  }
1433  }
1434  }
1436 
1437  VFN_DEBUG_EXIT("ScatteringData::CalcStructFactor()",3)
1438 }
1439 
1440 void ScatteringData::CalcStructFactor_FullDeriv(std::set<RefinablePar *> &vPar)
1441 {
1442  TAU_PROFILE("ScatteringData::CalcStructFactor_FullDeriv()","void ()",TAU_DEFAULT);
1444  this->CalcSinThetaLambda();
1445  this->CalcGeomStructFactor_FullDeriv(vPar);
1446  this->CalcStructFactor();//called after CalcGeomStructFactor_FullDeriv, so that CalcGeomStructFactor is not redone
1447 
1448  mFhklCalcReal_FullDeriv.clear();//:TODO: avoid full clear
1449  mFhklCalcImag_FullDeriv.clear();
1450  mFhklCalcReal_FullDeriv[0]=mFhklCalcReal;
1451  mFhklCalcImag_FullDeriv[0]=mFhklCalcImag;
1452  for(std::set<RefinablePar*>::iterator par=vPar.begin();par!=vPar.end();++par)
1453  {
1454  if(*par==0) continue;
1455  if((*par)->GetType()->IsDescendantFromOrSameAs(gpRefParTypeScatt)==false)
1456  {//:TODO: allow derivatives from other parameters (ML, temperature factors, etc..)
1457  // No derivatives -> empty vectors
1458  mFhklCalcReal_FullDeriv[*par].resize(0);
1459  mFhklCalcImag_FullDeriv[*par].resize(0);
1460  continue;
1461  }
1462  for(map<const ScatteringPower*,CrystVector_REAL>::const_iterator pos=mvRealGeomSF.begin();
1463  pos!=mvRealGeomSF.end();++pos)
1464  {
1465  const ScatteringPower* pScattPow=pos->first;
1466  if(mvRealGeomSF_FullDeriv[*par][pScattPow].size()==0)
1467  {
1468  continue;//null derivative, so the array was empty
1469  }
1470  if(mFhklCalcReal_FullDeriv[*par].size()==0)
1471  {
1472  mFhklCalcReal_FullDeriv[*par].resize(mNbRefl);
1473  mFhklCalcImag_FullDeriv[*par].resize(mNbRefl);
1474  mFhklCalcReal_FullDeriv[*par]=0;
1475  mFhklCalcImag_FullDeriv[*par]=0;
1476  }
1477  const REAL * RESTRICT pGeomRd=mvRealGeomSF_FullDeriv[*par][pScattPow].data();
1478  const REAL * RESTRICT pGeomId=mvImagGeomSF_FullDeriv[*par][pScattPow].data();
1479  const REAL * RESTRICT pScatt=mvScatteringFactor[pScattPow].data();
1480  const REAL * RESTRICT pTemp=mvTemperatureFactor[pScattPow].data();
1481 
1482  REAL * RESTRICT pReal=mFhklCalcReal_FullDeriv[*par].data();
1483  REAL * RESTRICT pImag=mFhklCalcImag_FullDeriv[*par].data();
1484  if(mvLuzzatiFactor[pScattPow].numElements()>0)
1485  {// using maximum likelihood
1486  const REAL* RESTRICT pLuzzati=mvLuzzatiFactor[pScattPow].data();
1487  if(false==mIgnoreImagScattFact)
1488  {
1489  const REAL fsecond=mvFsecond[pScattPow];
1490  for(long j=mNbReflUsed;j>0;j--)
1491  {
1492  *pReal++ += (*pGeomRd * *pScatt - *pGeomId * fsecond)* *pTemp * *pLuzzati;
1493  *pImag++ += (*pGeomId++ * *pScatt++ + *pGeomRd++ * fsecond)* *pTemp++ * *pLuzzati++;
1494  }
1495  }
1496  else
1497  {
1498  for(long j=mNbReflUsed;j>0;j--)
1499  {
1500  *pReal++ += *pGeomRd++ * *pTemp * *pScatt * *pLuzzati;
1501  *pImag++ += *pGeomId++ * *pTemp++ * *pScatt++ * *pLuzzati++;
1502  }
1503  }
1504  }
1505  else
1506  {
1507  if(false==mIgnoreImagScattFact)
1508  {
1509  const REAL fsecond=mvFsecond[pScattPow];
1510  for(long j=mNbReflUsed;j>0;j--)
1511  {
1512  *pReal += (*pGeomRd * *pScatt - *pGeomId * fsecond)* *pTemp;
1513  *pImag += (*pGeomId * *pScatt + *pGeomRd * fsecond)* *pTemp;
1514  pGeomRd++;pGeomId++;pTemp++;pScatt++;pReal++;pImag++;
1515  }
1516  }
1517  else
1518  {
1519  for(long j=mNbReflUsed;j>0;j--)
1520  {
1521  *pReal++ += *pGeomRd++ * *pTemp * *pScatt;
1522  *pImag++ += *pGeomId++ * *pTemp++ * *pScatt++;
1523  }
1524  }
1525  }
1526  }
1527  //TAU_PROFILE_STOP(timer4);
1528  {
1529  //this->CalcGlobalTemperatureFactor();
1530  if( (mGlobalTemperatureFactor.numElements()>0)
1531  &&(mFhklCalcReal_FullDeriv[*par].size()>0)
1532  &&(mFhklCalcImag_FullDeriv[*par].size()>0))
1533  {//else for some reason it's useless
1534  REAL * RESTRICT pReal=mFhklCalcReal_FullDeriv[*par].data();
1535  REAL * RESTRICT pImag=mFhklCalcImag_FullDeriv[*par].data();
1536  const REAL *pTemp=mGlobalTemperatureFactor.data();
1537  for(long j=0;j<mNbReflUsed;j++)
1538  {
1539  *pReal++ *= *pTemp;
1540  *pImag++ *= *pTemp++;
1541  }
1542  }
1543  }
1544  }
1545  #if 0
1546  std::vector<const CrystVector_REAL*> v;
1547  v.push_back(&mH);
1548  v.push_back(&mK);
1549  v.push_back(&mL);
1550  std::map<RefinablePar*, CrystVector_REAL> oldDerivR,oldDerivI;
1551  for(std::set<RefinablePar*>::iterator par=vPar.begin();par!=vPar.end();++par)
1552  {
1553  const REAL step=(*par)->GetDerivStep();
1554  (*par)->Mutate(step);
1555  this->CalcStructFactor();
1556  oldDerivR[*par]=mFhklCalcReal;
1557  oldDerivI[*par]=mFhklCalcImag;
1558  (*par)->Mutate(-2*step);
1559  this->CalcStructFactor();
1560  oldDerivR[*par]-=mFhklCalcReal;
1561  oldDerivR[*par]/=2*step;
1562  oldDerivI[*par]-=mFhklCalcImag;
1563  oldDerivI[*par]/=2*step;
1564  (*par)->Mutate(step);
1565 
1566  v.push_back(&(mFhklCalcReal_FullDeriv[*par]));
1567  v.push_back(&(oldDerivR[*par]));
1568  v.push_back(&(mFhklCalcImag_FullDeriv[*par]));
1569  v.push_back(&(oldDerivI[*par]));
1570  if(v.size()>14) break;
1571  }
1572  cout<<"############################ Fhkl Deriv Real, Imag ##############################"
1573  <<endl<<FormatVertVectorHKLFloats<REAL>(v,14,4,20)
1574  <<"############################ END Fhkl Deriv Real, Imag ##############################"<<endl;
1575  //exit(0);
1576  #endif
1577 }
1578 
1580 {
1581  // This also updates the ScattCompList if necessary.
1582  const ScatteringComponentList *pScattCompList
1583  =&(this->GetCrystal().GetScatteringComponentList());
1588  TAU_PROFILE("ScatteringData::GeomStructFactor()","void (Vx,Vy,Vz,data,M,M,bool)",TAU_DEFAULT);
1589  VFN_DEBUG_ENTRY("ScatteringData::GeomStructFactor(Vx,Vy,Vz,...)",3)
1590  VFN_DEBUG_MESSAGE("-->Using fast functions:"<<mUseFastLessPreciseFunc,2)
1591  VFN_DEBUG_MESSAGE("-->Number of translation vectors:"
1592  <<this->GetCrystal().GetSpaceGroup().GetNbTranslationVectors()-1,2)
1593  VFN_DEBUG_MESSAGE("-->Has an inversion Center:"
1594  <<this->GetCrystal().GetSpaceGroup().HasInversionCenter(),2)
1595  VFN_DEBUG_MESSAGE("-->Number of symetry operations (w/o transl&inv cent.):"\
1596  <<this->GetCrystal().GetSpaceGroup().GetNbSymmetrics(true,true),2)
1597  VFN_DEBUG_MESSAGE("-->Number of Scattering Components :"
1598  <<this->GetCrystal().GetScatteringComponentList().GetNbComponent(),2)
1599  VFN_DEBUG_MESSAGE("-->Number of reflections:"
1600  <<this->GetNbRefl()<<" (actually used:"<<mNbReflUsed<<")",2)
1601  #ifdef __DEBUG__
1602  static long counter=0;
1603  VFN_DEBUG_MESSAGE("-->Number of GeomStructFactor calculations so far:"<<counter++,3)
1604  #endif
1605 
1606  //:TODO: implement for geometrical structure factor calculation
1607  //bool useGeomStructFactor=mUseGeomStructFactor;
1608 
1609  //if((mfpRealGeomStructFactor==0)||(mfpImagGeomStructFactor==0)) useGeomStructFactor=false ;
1610 
1611  //if(useGeomStructFactor==true)
1612  //{
1613  // (*mfpRealGeomStructFactor)(x,y,z,data.H2Pi(),data.K2Pi(),data.L2Pi(),rsf);
1614  // if(this->IsCentrosymmetric())return;
1615  // (*mfpImagGeomStructFactor)(x,y,z,data.H2Pi(),data.K2Pi(),data.L2Pi(),isf);
1616  // return;
1617  //}
1618  //else
1619  {
1620  const SpaceGroup *pSpg=&(this->GetCrystal().GetSpaceGroup());
1621 
1622  const int nbSymmetrics=pSpg->GetNbSymmetrics(true,true);
1623  const int nbTranslationVectors=pSpg->GetNbTranslationVectors();
1624  const long nbComp=pScattCompList->GetNbComponent();
1625  const std::vector<SpaceGroup::TRx> *pTransVect=&(pSpg->GetTranslationVectors());
1626  CrystMatrix_REAL allCoords(nbSymmetrics,3);
1627  CrystVector_REAL tmpVect(mNbReflUsed);
1628  #ifndef HAVE_SSE_MATHFUN
1629  const int nbRefl=this->GetNbRefl();
1630  CrystVector_long intVect(nbRefl);//not used if mUseFastLessPreciseFunc==false
1631  #endif
1632  // which scattering powers are actually used ?
1633  map<const ScatteringPower*,bool> vUsed;
1634  // Add existing previously used scattering power to the test;
1635  for(map<const ScatteringPower*,CrystVector_REAL>::const_iterator pos=mvRealGeomSF.begin();pos!=mvRealGeomSF.end();++pos)
1636  vUsed[pos->first]=false;// this will be changed to true later if they are actually used
1637 
1638  for(int i=mpCrystal->GetScatteringPowerRegistry().GetNb()-1;i>=0;i--)
1639  {// Here we make sure scattering power that only contribute ghost atoms are taken into account
1640  const ScatteringPower*pow=&(mpCrystal->GetScatteringPowerRegistry().GetObj(i));
1641  if(pow->GetMaximumLikelihoodNbGhostAtom()>0) vUsed[pow]=true;
1642  else vUsed[pow]=false;
1643  }
1644  for(long i=0;i<nbComp;i++)
1645  vUsed[(*pScattCompList)(i).mpScattPow]=true;
1646  //Resize all arrays and set them to 0
1647  for(map<const ScatteringPower*,bool>::const_iterator pos=vUsed.begin();pos!=vUsed.end();++pos)
1648  {
1649  if(pos->second)
1650  {// this will create the entry if it does not already exist
1651  mvRealGeomSF[pos->first].resize(mNbReflUsed);
1652  mvImagGeomSF[pos->first].resize(mNbReflUsed);
1653  mvRealGeomSF[pos->first]=0;
1654  mvImagGeomSF[pos->first]=0;
1655  }
1656  else
1657  {// erase entries that are not useful any more (e.g. ScatteringPower that were
1658  // used but are not any more).
1659  map<const ScatteringPower*,CrystVector_REAL>::iterator
1660  poubelle=mvRealGeomSF.find(pos->first);
1661  if(poubelle!=mvRealGeomSF.end()) mvRealGeomSF.erase(poubelle);
1662  poubelle=mvImagGeomSF.find(pos->first);
1663  if(poubelle!=mvImagGeomSF.end()) mvImagGeomSF.erase(poubelle);
1664  }
1665  }
1666 
1667  REAL centrMult=1.0;
1668  if(true==pSpg->HasInversionCenter()) centrMult=2.0;
1669  for(long i=0;i<nbComp;i++)
1670  {
1671  VFN_DEBUG_MESSAGE("ScatteringData::GeomStructFactor(),comp"<<i,3)
1672  const REAL x=(*pScattCompList)(i).mX;
1673  const REAL y=(*pScattCompList)(i).mY;
1674  const REAL z=(*pScattCompList)(i).mZ;
1675  const ScatteringPower *pScattPow=(*pScattCompList)(i).mpScattPow;
1676  const REAL popu= (*pScattCompList)(i).mOccupancy
1677  *(*pScattCompList)(i).mDynPopCorr
1678  *centrMult;
1679 
1680  allCoords=pSpg->GetAllSymmetrics(x,y,z,true,true);
1681  if((true==pSpg->HasInversionCenter()) && (false==pSpg->IsInversionCenterAtOrigin()))
1682  {
1683  const REAL STBF=2.*pSpg->GetCCTbxSpg().inv_t().den();
1684  for(int j=0;j<nbSymmetrics;j++)
1685  {
1686  //The phase of the structure factor will be wrong
1687  //This is fixed a bit further...
1688  allCoords(j,0) -= ((REAL)pSpg->GetCCTbxSpg().inv_t()[0])/STBF;
1689  allCoords(j,1) -= ((REAL)pSpg->GetCCTbxSpg().inv_t()[1])/STBF;
1690  allCoords(j,2) -= ((REAL)pSpg->GetCCTbxSpg().inv_t()[2])/STBF;
1691  }
1692  }
1693  for(int j=0;j<nbSymmetrics;j++)
1694  {
1695  VFN_DEBUG_MESSAGE("ScatteringData::GeomStructFactor(),comp #"<<i<<", sym #"<<j,3)
1696 
1697  #ifndef HAVE_SSE_MATHFUN
1698  if(mUseFastLessPreciseFunc==true)
1699  {
1700  REAL * RESTRICT rrsf=mvRealGeomSF[pScattPow].data();
1701  REAL * RESTRICT iisf=mvImagGeomSF[pScattPow].data();
1702 
1703  const long intX=(long)(allCoords(j,0)*sLibCrystNbTabulSine);
1704  const long intY=(long)(allCoords(j,1)*sLibCrystNbTabulSine);
1705  const long intZ=(long)(allCoords(j,2)*sLibCrystNbTabulSine);
1706 
1707  const long * RESTRICT intH=mIntH.data();
1708  const long * RESTRICT intK=mIntK.data();
1709  const long * RESTRICT intL=mIntL.data();
1710 
1711  register long * RESTRICT tmpInt=intVect.data();
1712  // :KLUDGE: using a AND to bring back within [0;sLibCrystNbTabulSine[ may
1713  // not be portable, depending on the model used to represent signed integers
1714  // a test should be added to throw up in that case.
1715  //
1716  // This work if we are using "2's complement" to represent negative numbers,
1717  // but not with a "sign magnitude" approach
1718  for(int jj=mNbReflUsed;jj>0;jj--)
1719  *tmpInt++ = (*intH++ * intX + *intK++ * intY + *intL++ *intZ)
1720  &sLibCrystNbTabulSineMASK;
1721  if(false==pSpg->HasInversionCenter())
1722  {
1723 
1724  tmpInt=intVect.data();
1725  for(int jj=mNbReflUsed;jj>0;jj--)
1726  {
1727  const REAL *pTmp=&spLibCrystTabulCosineSine[*tmpInt++ <<1];
1728  *rrsf++ += popu * *pTmp++;
1729  *iisf++ += popu * *pTmp;
1730  }
1731 
1732  }
1733  else
1734  {
1735  tmpInt=intVect.data();
1736  for(int jj=mNbReflUsed;jj>0;jj--)
1737  *rrsf++ += popu * spLibCrystTabulCosine[*tmpInt++];
1738  }
1739  }
1740  else
1741  #endif
1742  {
1743  const REAL x=allCoords(j,0);
1744  const REAL y=allCoords(j,1);
1745  const REAL z=allCoords(j,2);
1746  const register REAL *hh=mH2Pi.data();
1747  const register REAL *kk=mK2Pi.data();
1748  const register REAL *ll=mL2Pi.data();
1749 
1750  #ifdef HAVE_SSE_MATHFUN
1751  #if 0
1752  // This not much faster and is incorrect (does not take into account sign of h k l)
1753 
1754  //cout<<__FILE__<<":"<<__LINE__<<":"<<mMaxHKL<<","<<mMaxH<<","<<mMaxK<<","<<mMaxL<<":"<<mNbReflUsed<<endl;
1755  // cos&sin for 2pix 2piy 2piz
1756  static const float twopi=6.283185307179586f;
1757  sincos_ps(_mm_mul_ps(_mm_load1_ps(&twopi),_mm_set_ps(x,y,z,0)),cnxyz0,snxyz0);
1758  // harmonics: cos&sin for 2npix 2npiy 2npiz
1759  for(long k=1;k<mMaxHKL;k++)
1760  {
1761  cnxyz0[k]=_mm_sub_ps(_mm_mul_ps(cnxyz0[k-1],cnxyz0[0]),_mm_mul_ps(snxyz0[k-1],snxyz0[0]));//cos((n+1)x)=cos(nx)cos(x)-sin(nx)sinx
1762  snxyz0[k]=_mm_add_ps(_mm_mul_ps(snxyz0[k-1],cnxyz0[0]),_mm_mul_ps(cnxyz0[k-1],snxyz0[0]));//sin((n+1)x)=sin(nx)cos(x)+cos(nx)sinx
1763  }
1764  //
1765  for(long k=0;k<4;++k){*(pcnxyz0+k)=1.0f;*(psnxyz0+k)=0.0f;}
1766  for(long k=1;k<mMaxHKL;k++)
1767  {
1768  _mm_store_ps(pcnxyz0+4*k,cnxyz0[k-1]);
1769  _mm_store_ps(psnxyz0+4*k,snxyz0[k-1]);
1770  }
1771  // Actual structure factor calculations
1772  if(false==pSpg->HasInversionCenter())
1773  {// Slow ?
1774  REAL *rsf=mvRealGeomSF[pScattPow].data();
1775  REAL *isf=mvImagGeomSF[pScattPow].data();
1776  const long *h=mIntH.data();
1777  const long *k=mIntK.data();
1778  const long *l=mIntL.data();
1779  int jj;
1780  const v4sf v4popu=_mm_set1_ps(popu);
1781  for(jj=mNbReflUsed;jj>3;jj-=4)
1782  {
1783  //cout<<__FILE__<<":"<<__LINE__<<":"<<mNbReflUsed<<","<<jj<<"("<<*h<<','<<*k<<","<<*l<<")"<<endl;
1784  const v4sf ck=_mm_set_ps(pcnxyz0[(*(k))*4+1],pcnxyz0[(*(k+1))*4+1],pcnxyz0[(*(k+2))*4+1],pcnxyz0[(*(k+3))*4+1]);//cos 2pi kx =ck
1785  const v4sf cl=_mm_set_ps(pcnxyz0[(*(l))*4+1],pcnxyz0[(*(l+1))*4+1],pcnxyz0[(*(l+2))*4+1],pcnxyz0[(*(l+3))*4+1]);//cos 2pi lz =cl
1786  const v4sf sk=_mm_set_ps(psnxyz0[(*(k))*4+1],psnxyz0[(*(k+1))*4+1],psnxyz0[(*(k+2))*4+1],psnxyz0[(*(k+3))*4+1]);//sin 2pi kx =sk
1787  const v4sf sl=_mm_set_ps(psnxyz0[(*(l))*4+1],psnxyz0[(*(l+1))*4+1],psnxyz0[(*(l+2))*4+1],psnxyz0[(*(l+3))*4+1]);//sin 2pi lz =sl
1788  #define CH _mm_set_ps(pcnxyz0[*(h)*4],pcnxyz0[*(h+1)*4],pcnxyz0[*(h+2)*4],pcnxyz0[*(h+3)*4])
1789  #define SH _mm_set_ps(psnxyz0[*(h)*4],psnxyz0[*(h+1)*4],psnxyz0[*(h+2)*4],psnxyz0[*(h+3)*4])
1790  // popu *( ch*( ck*cl - sk*sl) - sh*( ck*sl + sk*cl))
1791  _mm_store_ps(rsf,_mm_mul_ps(v4popu,_mm_sub_ps(_mm_mul_ps(CH,_mm_sub_ps(_mm_mul_ps(ck,cl),_mm_mul_ps(sk,sl))),_mm_mul_ps(SH,_mm_add_ps(_mm_mul_ps(ck,sl),_mm_mul_ps(sk,cl))))));
1792  // popu *( sh*( ck*cl - sk*sl) + ch*( ck*sl + sk*cl))
1793  _mm_store_ps(isf,_mm_mul_ps(v4popu,_mm_add_ps(_mm_mul_ps(SH,_mm_sub_ps(_mm_mul_ps(ck,cl),_mm_mul_ps(sk,sl))),_mm_mul_ps(CH,_mm_add_ps(_mm_mul_ps(ck,sl),_mm_mul_ps(sk,cl))))));
1794  rsf+=4;isf+=4;h+=4;k+=4,l+=4;
1795  }
1796  for(;jj>0;jj--)
1797  {
1798  const float ch=pcnxyz0[*h *4];
1799  const float sh=psnxyz0[*h++ *4];
1800  const float ck=pcnxyz0[*k *4+1];
1801  const float sk=psnxyz0[*k++ *4+1];
1802  const float cl=pcnxyz0[*l *4+2];
1803  const float sl=psnxyz0[*l++ *4+2];
1804  *rsf++ += popu*(ch*(ck*cl-sk*sl)-sh*(sk*cl+ck*sl));
1805  *isf++ += popu*(sh*(ck*cl-sk*sl)+ch*(sk*cl+ck*sl));
1806  }
1807  }
1808  else
1809  {
1810  REAL *rsf=mvRealGeomSF[pScattPow].data();
1811  const long *h=mIntH.data();
1812  const long *k=mIntK.data();
1813  const long *l=mIntL.data();
1814  int jj;
1815  const v4sf v4popu=_mm_set1_ps(popu);
1816  for(jj=mNbReflUsed;jj>3;jj-=4)
1817  {
1818  //cout<<__FILE__<<":"<<__LINE__<<":"<<mNbReflUsed<<","<<jj<<"("<<*h<<','<<*k<<","<<*l<<")"<<endl;
1819  const v4sf ck=_mm_set_ps(pcnxyz0[(*(k))*4+1],pcnxyz0[(*(k+1))*4+1],pcnxyz0[(*(k+2))*4+1],pcnxyz0[(*(k+3))*4+1]);//cos 2pi kx =ck
1820  const v4sf cl=_mm_set_ps(pcnxyz0[(*(l))*4+1],pcnxyz0[(*(l+1))*4+1],pcnxyz0[(*(l+2))*4+1],pcnxyz0[(*(l+3))*4+1]);//cos 2pi lz =cl
1821  const v4sf sk=_mm_set_ps(psnxyz0[(*(k))*4+1],psnxyz0[(*(k+1))*4+1],psnxyz0[(*(k+2))*4+1],psnxyz0[(*(k+3))*4+1]);//sin 2pi kx =sk
1822  const v4sf sl=_mm_set_ps(psnxyz0[(*(l))*4+1],psnxyz0[(*(l+1))*4+1],psnxyz0[(*(l+2))*4+1],psnxyz0[(*(l+3))*4+1]);//sin 2pi lz =sl
1823  #define CH _mm_set_ps(pcnxyz0[*(h)*4],pcnxyz0[*(h+1)*4],pcnxyz0[*(h+2)*4],pcnxyz0[*(h+3)*4])
1824  #define SH _mm_set_ps(psnxyz0[*(h)*4],psnxyz0[*(h+1)*4],psnxyz0[*(h+2)*4],psnxyz0[*(h+3)*4])
1825  // popu *( ch*( ck*cl - sk*sl) - sh*( ck*sl + sk*cl))
1826  _mm_store_ps(rsf,_mm_mul_ps(v4popu,_mm_sub_ps(_mm_mul_ps(CH,_mm_sub_ps(_mm_mul_ps(ck,cl),_mm_mul_ps(sk,sl))),_mm_mul_ps(SH,_mm_add_ps(_mm_mul_ps(ck,sl),_mm_mul_ps(sk,cl))))));
1827  rsf+=4;h+=4;k+=4,l+=4;
1828  }
1829  for(;jj>0;jj--)
1830  {
1831  const float ch=pcnxyz0[*h *4];
1832  const float sh=psnxyz0[*h++ *4];
1833  const float ck=pcnxyz0[*k *4+1];
1834  const float sk=psnxyz0[*k++ *4+1];
1835  const float cl=pcnxyz0[*l *4+2];
1836  const float sl=psnxyz0[*l++ *4+2];
1837  *rsf++ += popu*(ch*(ck*cl-sk*sl)-sh*(sk*cl+ck*sl));
1838  }
1839  }
1840 
1841 
1842  #else
1843  const v4sf v4x=_mm_load1_ps(&x);
1844  const v4sf v4y=_mm_load1_ps(&y);
1845  const v4sf v4z=_mm_load1_ps(&z);
1846  const v4sf v4popu=_mm_load1_ps(&popu);// Can't multiply directly a vector by a scalar ?
1847  if(false==pSpg->HasInversionCenter())
1848  {
1849  REAL *rsf=mvRealGeomSF[pScattPow].data();
1850  REAL *isf=mvImagGeomSF[pScattPow].data();
1851  int jj=mNbReflUsed;
1852  for(;jj>3;jj-=4)
1853  {
1854  v4sf v4sin,v4cos;
1855 // sincos_ps(_mm_setr_ps(*(hh )*x+ *(kk )*y + *(ll )*z,
1856 // *(hh+1)*x+ *(kk+1)*y + *(ll+1)*z,
1857 // *(hh+2)*x+ *(kk+2)*y + *(ll+2)*z,
1858 // *(hh+3)*x+ *(kk+3)*y + *(ll+3)*z),&v4sin,&v4cos);
1859  sincos_ps(_mm_add_ps(_mm_add_ps(_mm_mul_ps(_mm_loadu_ps(hh),v4x),
1860  _mm_mul_ps(_mm_loadu_ps(kk),v4y)
1861  ),
1862  _mm_mul_ps(_mm_loadu_ps(ll),v4z)
1863  ),&v4sin,&v4cos);// A bit faster
1864  _mm_storeu_ps(rsf,_mm_add_ps(_mm_mul_ps(v4cos,v4popu),_mm_loadu_ps(rsf)));
1865  _mm_storeu_ps(isf,_mm_add_ps(_mm_mul_ps(v4sin,v4popu),_mm_loadu_ps(isf)));
1866 
1867  hh+=4;kk+=4;ll+=4;rsf+=4;isf+=4;
1868  }
1869  for(;jj>0;jj--)
1870  {
1871  const REAL tmp = *hh++ * x + *kk++ * y + *ll++ *z;
1872  *rsf++ += popu * cos(tmp);
1873  *isf++ += popu * sin(tmp);
1874  }
1875  }
1876  else
1877  {
1878  REAL *rsf=mvRealGeomSF[pScattPow].data();
1879  int jj=mNbReflUsed;
1880  for(;jj>3;jj-=4)
1881  {
1882 // const v4sf v4cos=cos_ps(_mm_setr_ps(*(hh )*x+ *(kk )*y + *(ll )*z,
1883 // *(hh+1)*x+ *(kk+1)*y + *(ll+1)*z,
1884 // *(hh+2)*x+ *(kk+2)*y + *(ll+2)*z,
1885 // *(hh+3)*x+ *(kk+3)*y + *(ll+3)*z));
1886  const v4sf v4cos=cos_ps(_mm_add_ps(_mm_add_ps(_mm_mul_ps(_mm_loadu_ps(hh),v4x),
1887  _mm_mul_ps(_mm_loadu_ps(kk),v4y)
1888  ),
1889  _mm_mul_ps(_mm_loadu_ps(ll),v4z)));
1890  _mm_storeu_ps(rsf,_mm_add_ps(_mm_loadu_ps(rsf),_mm_mul_ps(v4cos,v4popu)));
1891  hh+=4;kk+=4;ll+=4;rsf+=4;
1892  }
1893  for(;jj>0;jj--)
1894  {
1895  const REAL tmp = *hh++ * x + *kk++ * y + *ll++ *z;
1896  *rsf++ += popu * cos(tmp);
1897  }
1898  }
1899  #endif
1900  #else
1901  register REAL *tmp=tmpVect.data();
1902  for(int jj=0;jj<mNbReflUsed;jj++) *tmp++ = *hh++ * x + *kk++ * y + *ll++ *z;
1903 
1904  REAL *sf=mvRealGeomSF[pScattPow].data();
1905  tmp=tmpVect.data();
1906 
1907  for(int jj=0;jj<mNbReflUsed;jj++) *sf++ += popu * cos(*tmp++);
1908 
1909  if(false==pSpg->HasInversionCenter())
1910  {
1911  sf=mvImagGeomSF[pScattPow].data();
1912  tmp=tmpVect.data();
1913  for(int jj=0;jj<mNbReflUsed;jj++) *sf++ += popu * sin(*tmp++);
1914  }
1915  #endif
1916  }
1917  }
1918  }//for all components...
1919  if(nbTranslationVectors > 1)
1920  {
1921  tmpVect=1;
1922  if( (pSpg->GetSpaceGroupNumber()>= 143) && (pSpg->GetSpaceGroupNumber()<= 167))
1923  {//Special case for trigonal groups R3,...
1924  REAL * RESTRICT p1=tmpVect.data();
1925  const register REAL * RESTRICT hh=mH2Pi.data();
1926  const register REAL * RESTRICT kk=mK2Pi.data();
1927  const register REAL * RESTRICT ll=mL2Pi.data();
1928  for(long j=mNbReflUsed;j>0;j--) *p1++ += 2*cos((*hh++ - *kk++ - *ll++)/3.);
1929  }
1930  else
1931  {
1932  for(int j=1;j<nbTranslationVectors;j++)
1933  {
1934  const REAL x=(*pTransVect)[j].tr[0];
1935  const REAL y=(*pTransVect)[j].tr[1];
1936  const REAL z=(*pTransVect)[j].tr[2];
1937  REAL *p1=tmpVect.data();
1938  const register REAL *hh=mH2Pi.data();
1939  const register REAL *kk=mK2Pi.data();
1940  const register REAL *ll=mL2Pi.data();
1941  for(long j=mNbReflUsed;j>0;j--) *p1++ += cos(*hh++ *x + *kk++ *y + *ll++ *z );
1942  }
1943  }
1944  for(map<const ScatteringPower*,CrystVector_REAL>::iterator
1945  pos=mvRealGeomSF.begin();pos!=mvRealGeomSF.end();++pos)
1946  pos->second *= tmpVect;
1947 
1948  if(false==pSpg->HasInversionCenter())
1949  for(map<const ScatteringPower*,CrystVector_REAL>::iterator
1950  pos=mvImagGeomSF.begin();pos!=mvImagGeomSF.end();++pos)
1951  pos->second *= tmpVect;
1952  }
1953  if(true==pSpg->HasInversionCenter())
1954  {
1955  // we already multiplied real geom struct factor by 2
1956  if(false==pSpg->IsInversionCenterAtOrigin())
1957  {
1958  VFN_DEBUG_MESSAGE("ScatteringData::GeomStructFactor(Vx,Vy,Vz):\
1959  Inversion Center not at the origin...",2)
1960  //fix the phase of each reflection when the inversion center is not
1961  //at the origin, using :
1962  // Re(F) = RSF*cos(2pi(h*Xc+k*Yc+l*Zc))
1963  // Re(F) = RSF*sin(2pi(h*Xc+k*Yc+l*Zc))
1964  //cout << "Glop Glop"<<endl;
1965  const REAL STBF=2*pSpg->GetCCTbxSpg().inv_t().den();
1966  {
1967  const REAL xc=((REAL)pSpg->GetCCTbxSpg().inv_t()[0])/STBF;
1968  const REAL yc=((REAL)pSpg->GetCCTbxSpg().inv_t()[1])/STBF;
1969  const REAL zc=((REAL)pSpg->GetCCTbxSpg().inv_t()[2])/STBF;
1970  #ifdef __LIBCRYST_VECTOR_USE_BLITZ__
1971  tmpVect = mH2Pi() * xc + mK2PI() * yc + mL2PI() * zc;
1972  #else
1973  {
1974  const REAL * RESTRICT hh=mH2Pi.data();
1975  const REAL * RESTRICT kk=mK2Pi.data();
1976  const REAL * RESTRICT ll=mL2Pi.data();
1977  REAL * RESTRICT ttmpVect=tmpVect.data();
1978  for(long ii=mNbReflUsed;ii>0;ii--)
1979  *ttmpVect++ = *hh++ * xc + *kk++ * yc + *ll++ * zc;
1980  }
1981  #endif
1982  }
1983  CrystVector_REAL cosTmpVect;
1984  CrystVector_REAL sinTmpVect;
1985  cosTmpVect=cos(tmpVect);
1986  sinTmpVect=sin(tmpVect);
1987 
1988  map<const ScatteringPower*,CrystVector_REAL>::iterator posi=mvImagGeomSF.begin();
1989  map<const ScatteringPower*,CrystVector_REAL>::iterator posr=mvRealGeomSF.begin();
1990  for(;posi!=mvImagGeomSF.end();)
1991  {
1992  posi->second = posr->second;
1993  posi->second *= sinTmpVect;
1994  posr->second *= cosTmpVect;
1995  posi++;posr++;
1996  }
1997  }
1998  }
1999  }
2000  //cout << FormatVertVector<REAL>(*mvRealGeomSF,*mvImagGeomSF)<<endl;
2002  VFN_DEBUG_EXIT("ScatteringData::GeomStructFactor(Vx,Vy,Vz,...)",3)
2003 }
2004 void ScatteringData::CalcGeomStructFactor_FullDeriv(std::set<RefinablePar*> &vPar)
2005 {
2006  TAU_PROFILE("ScatteringData::CalcGeomStructFactor_FullDeriv()","void (..)",TAU_DEFAULT);
2007  TAU_PROFILE_TIMER(timer1,"ScatteringData::CalcGeomStructFactor_FullDeriv:1-ScattCompList deriv","", TAU_FIELD);
2008  TAU_PROFILE_TIMER(timer2,"ScatteringData::CalcGeomStructFactor_FullDeriv:2-Geom SF","", TAU_FIELD);
2009  this->CalcGeomStructFactor();//:TODO: avoid calling CalcGeomStructFactor()
2010  //:TODO: this->GetCrystal().GetScatteringComponentList_FullDeriv()
2011  const ScatteringComponentList *pScattCompList
2012  =&(this->GetCrystal().GetScatteringComponentList());
2013 
2014  const SpaceGroup *pSpg=&(this->GetCrystal().GetSpaceGroup());
2015 
2016  const int nbSymmetrics=pSpg->GetNbSymmetrics(true,true);
2017  const int nbTranslationVectors=pSpg->GetNbTranslationVectors();
2018  const unsigned long nbComp=pScattCompList->GetNbComponent();
2019  const std::vector<SpaceGroup::TRx> *pTransVect=&(pSpg->GetTranslationVectors());
2020  CrystMatrix_REAL allCoords(nbSymmetrics,3);
2021 
2022  const bool hasinv=pSpg->HasInversionCenter();
2023 
2024  TAU_PROFILE_START(timer1);
2025  // Calculate derivatives of the scattering component list vs all parameters
2026  std::map<RefinablePar*,CrystVector_REAL> vdx,vdy,vdz,vdocc;
2027  for(std::set<RefinablePar*>::iterator par=vPar.begin();par!=vPar.end();++par)
2028  {// :TODO: get this done in Crystal or Scatterers, and use analytical derivatives
2029  if(*par==0) continue;
2030  CrystVector_REAL *pdx =&(vdx[*par]);
2031  CrystVector_REAL *pdy =&(vdy[*par]);
2032  CrystVector_REAL *pdz =&(vdz[*par]);
2033  CrystVector_REAL *pdocc=&(vdocc[*par]);
2034  pdx->resize(nbComp);
2035  pdy->resize(nbComp);
2036  pdz->resize(nbComp);
2037  pdocc->resize(nbComp);
2038 
2039  const REAL p0=(*par)->GetValue();
2040  const REAL step=(*par)->GetDerivStep();
2041  (*par)->Mutate(step);
2042  pScattCompList=&(this->GetCrystal().GetScatteringComponentList());
2043  REAL *ppdx =pdx->data();
2044  REAL *ppdy =pdy->data();
2045  REAL *ppdz =pdz->data();
2046  REAL *ppdocc=pdocc->data();
2047  for(unsigned long i=0;i<nbComp;++i)
2048  {
2049  *ppdx++ =(*pScattCompList)(i).mX;
2050  *ppdy++ =(*pScattCompList)(i).mY;
2051  *ppdz++ =(*pScattCompList)(i).mZ;
2052  *ppdocc++=(*pScattCompList)(i).mOccupancy*(*pScattCompList)(i).mDynPopCorr;
2053  }
2054  (*par)->Mutate(-2*step);
2055  pScattCompList=&(this->GetCrystal().GetScatteringComponentList());
2056  ppdx =pdx->data();
2057  ppdy =pdy->data();
2058  ppdz =pdz->data();
2059  ppdocc=pdocc->data();
2060  for(unsigned long i=0;i<nbComp;++i)
2061  {
2062  *ppdx -=(*pScattCompList)(i).mX;
2063  *ppdx++/=2*step;
2064  *ppdy -=(*pScattCompList)(i).mY;
2065  *ppdy++/=2*step;
2066  *ppdz -=(*pScattCompList)(i).mZ;
2067  *ppdz++/=2*step;
2068  *ppdocc-=(*pScattCompList)(i).mOccupancy*(*pScattCompList)(i).mDynPopCorr;
2069  *ppdocc++/=2*step;
2070  }
2071  (*par)->SetValue(p0);
2072  if( (MaxAbs(vdx[*par])==0)&&(MaxAbs(vdy[*par])==0)&&(MaxAbs(vdz[*par])==0)&&(MaxAbs(vdocc[*par])==0))
2073  {
2074  pdx->resize(0);
2075  pdy->resize(0);
2076  pdz->resize(0);
2077  pdocc->resize(0);
2078  }
2079  }
2080  TAU_PROFILE_STOP(timer1);
2081  TAU_PROFILE_START(timer2);
2082  CrystVector_REAL transMult(mNbReflUsed);
2083  if(!hasinv) transMult=1;
2084  else transMult=2;
2085  if(nbTranslationVectors > 1)
2086  {
2087  if( (pSpg->GetSpaceGroupNumber()>= 143) && (pSpg->GetSpaceGroupNumber()<= 167))
2088  {//Special case for trigonal groups R3,...
2089  REAL * RESTRICT p1=transMult.data();
2090  const register REAL * RESTRICT hh=mH2Pi.data();
2091  const register REAL * RESTRICT kk=mK2Pi.data();
2092  const register REAL * RESTRICT ll=mL2Pi.data();
2093  for(long j=mNbReflUsed;j>0;j--) *p1++ += 2*cos((*hh++ - *kk++ - *ll++)/3.);
2094  }
2095  else
2096  {
2097  for(int j=1;j<nbTranslationVectors;j++)
2098  {
2099  const REAL x=(*pTransVect)[j].tr[0];
2100  const REAL y=(*pTransVect)[j].tr[1];
2101  const REAL z=(*pTransVect)[j].tr[2];
2102  REAL *p1=transMult.data();
2103  const register REAL * RESTRICT hh=mH2Pi.data();
2104  const register REAL * RESTRICT kk=mK2Pi.data();
2105  const register REAL * RESTRICT ll=mL2Pi.data();
2106  for(long j=mNbReflUsed;j>0;j--) *p1++ += cos(*hh++ *x + *kk++ *y + *ll++ *z );
2107  }
2108  }
2109  }
2110 
2111  pScattCompList=&(this->GetCrystal().GetScatteringComponentList());
2112 
2113  mvRealGeomSF_FullDeriv.clear();//:TODO: avoid clearing memory as much as possible
2114  mvImagGeomSF_FullDeriv.clear();
2115  CrystVector_REAL c(mNbReflUsed),s(mNbReflUsed);
2116  CrystMatrix_REAL allCoordsDeriv(nbSymmetrics,3);
2117  for(unsigned long i=0;i<nbComp;i++)
2118  {
2119  const REAL x0=(*pScattCompList)(i).mX;
2120  const REAL y0=(*pScattCompList)(i).mY;
2121  const REAL z0=(*pScattCompList)(i).mZ;
2122  const ScatteringPower *pScattPow=(*pScattCompList)(i).mpScattPow;
2123  const REAL popu= (*pScattCompList)(i).mOccupancy
2124  *(*pScattCompList)(i).mDynPopCorr;
2125  allCoords=pSpg->GetAllSymmetrics(x0,y0,z0,true,true);
2126  for(int j=0;j<nbSymmetrics;j++)
2127  {
2128  const REAL x=allCoords(j,0);
2129  const REAL y=allCoords(j,1);
2130  const REAL z=allCoords(j,2);
2131  {
2132  REAL *pc=c.data();
2133  REAL *ps=s.data();
2134  const register REAL *hh=mH2Pi.data();
2135  const register REAL *kk=mK2Pi.data();
2136  const register REAL *ll=mL2Pi.data();
2137  #ifdef HAVE_SSE_MATHFUN
2138  const v4sf v4x=_mm_load1_ps(&x);
2139  const v4sf v4y=_mm_load1_ps(&y);
2140  const v4sf v4z=_mm_load1_ps(&z);
2141  // Can't multiply directly a vector by a scalar ?
2142  const v4sf v4popu=_mm_load1_ps(&popu); POSSIBLY_UNUSED(v4popu)
2143  int jj=mNbReflUsed;
2144  for(;jj>3;jj-=4)
2145  {
2146  v4sf v4sin,v4cos;
2147 // sincos_ps(_mm_setr_ps(*(hh )*x+ *(kk )*y + *(ll )*z,
2148 // *(hh+1)*x+ *(kk+1)*y + *(ll+1)*z,
2149 // *(hh+2)*x+ *(kk+2)*y + *(ll+2)*z,
2150 // *(hh+3)*x+ *(kk+3)*y + *(ll+3)*z),&v4sin,&v4cos);
2151  sincos_ps(_mm_add_ps(_mm_add_ps(_mm_mul_ps(_mm_load_ps(hh),v4x),
2152  _mm_mul_ps(_mm_load_ps(kk),v4y)
2153  ),
2154  _mm_mul_ps(_mm_load_ps(ll),v4z)
2155  ),&v4sin,&v4cos);
2156  _mm_store_ps(pc,v4cos);
2157  _mm_store_ps(ps,v4sin);
2158 
2159  hh+=4;kk+=4;ll+=4;pc+=4;ps+=4;
2160  }
2161  for(;jj>0;jj--)
2162  {
2163  const REAL tmp = *hh++ * x + *kk++ * y + *ll++ *z;
2164  *pc++ =cos(tmp);
2165  *ps++ =sin(tmp);
2166  }
2167  #else
2168  for(int jj=0;jj<mNbReflUsed;jj++)
2169  {
2170  const REAL tmp = *hh++ * x + *kk++ * y + *ll++ *z;
2171  *pc++ =cos(tmp);
2172  *ps++ =sin(tmp);
2173  }
2174  #endif
2175  }
2176  for(std::set<RefinablePar*>::iterator par=vPar.begin();par!=vPar.end();++par)
2177  {
2178  if((*par)==0) continue;
2179  if(vdx[*par].size()==0) continue;
2180  REAL dx =vdx[*par](i);
2181  REAL dy =vdy[*par](i);
2182  REAL dz =vdz[*par](i);
2183  const REAL dpopu=vdocc[*par](i);
2184 
2185  if((abs(dx)+abs(dy)+abs(dz)+abs(dpopu))==0) continue;
2186  if(mvRealGeomSF_FullDeriv[*par][pScattPow].size()==0)
2187  {
2188  mvRealGeomSF_FullDeriv[*par][pScattPow].resize(mNbRefl);
2189  mvRealGeomSF_FullDeriv[*par][pScattPow]=0;
2190  }
2191  if(mvImagGeomSF_FullDeriv[*par][pScattPow].size()==0)
2192  {
2193  mvImagGeomSF_FullDeriv[*par][pScattPow].resize(mNbRefl);
2194  mvImagGeomSF_FullDeriv[*par][pScattPow]=0;
2195  }
2196  pSpg->GetSymmetric(j,dx,dy,dz,true,true,true);
2197  const register REAL *hh=mH2Pi.data();
2198  const register REAL *kk=mK2Pi.data();
2199  const register REAL *ll=mL2Pi.data();
2200  const register REAL *pmult=transMult.data();
2201  register REAL *rsf=mvRealGeomSF_FullDeriv[*par][pScattPow].data();
2202  register REAL *isf=mvImagGeomSF_FullDeriv[*par][pScattPow].data();
2203  VFN_DEBUG_MESSAGE("ScatteringData::CalcGeomStructFactor_FullDeriv()comp="<<i<<", par="<<(*par)->GetName()<<", rs="<<mvRealGeomSF_FullDeriv[*par][pScattPow].size(),1)
2204  const REAL *pc=c.data();
2205  const REAL *ps=s.data();
2206  //cout<<setw(12)<<(*par)->GetName()<<":"<<setw(12)<<pScattPow->GetName()<<":"<<i<<","<<j
2207  // <<":x="<<setw(12)<<x<<",y="<<setw(12)<<y<<",z="<<setw(12)<<z
2208  // <<":dx="<<setw(12)<<dx<<",dy="<<setw(12)<<dy<<",dz="<<setw(12)<<dz<<",dpopu="<<setw(12)<<dpopu<<",popu="<<setw(12)<<popu<<",c0="<<setw(12)<<*pc<<",s0="<<setw(12)<<*ps<<endl;
2209 
2210  for(int jj=0;jj<mNbReflUsed;jj++)
2211  {// :TODO: directly calculate corrected intensities, instead of 1) geom 2) F and 3) F^2 4) corrected intensity... Faster, less storage !
2212  *rsf += (dpopu * *pc - popu* *ps * (*hh * dx + *kk * dy + *ll * dz))* *pmult;
2213  if(!hasinv) *isf += (dpopu * *ps + popu* *pc * (*hh * dx + *kk * dy + *ll * dz))* *pmult;
2214  //if(jj<6) cout<<" rsf0+="<<setw(12)<<(dpopu * *pc - popu* *ps * (*hh * dx + *kk * dy + *ll * dz))* *pmult<<" ("<<setw(12)<<*rsf
2215  // <<"),isf0+="<<setw(12)<<(dpopu * *ps + popu* *pc * (*hh * dx + *kk * dy + *ll * dz))* *pmult
2216  // <<" ("<<setw(12)<<*isf<<")"<<endl;
2217  ps++;pc++;hh++;kk++;ll++;pmult++;rsf++;isf++;
2218  }
2219  }
2220  }
2221  }
2222  if(true==pSpg->HasInversionCenter())
2223  {
2224  if(false==pSpg->IsInversionCenterAtOrigin())
2225  {
2226  //:TODO: if there is an inversion center not in (0,0,0), apply a constant phase
2227  }
2228  }
2229 
2230  TAU_PROFILE_STOP(timer2);
2231  #if 0
2232  std::vector<const CrystVector_REAL*> v;
2233  v.push_back(&mH);
2234  v.push_back(&mK);
2235  v.push_back(&mL);
2236  std::map< std::pair<const ScatteringPower *,RefinablePar*>,CrystVector_REAL> mr,mi;
2238 
2239  for(std::set<RefinablePar*>::iterator par=vPar.begin();par!=vPar.end();++par)
2240  {
2241  for(std::map<const ScatteringPower*,CrystVector_REAL>::iterator pos=mvRealGeomSF.begin();pos!=mvRealGeomSF.end();++pos)
2242  {
2243  const ScatteringPower *pScattPow=pos->first;
2244  cout<<(*par)->GetName()<<","<<pScattPow->GetName();
2245  if(mvRealGeomSF_FullDeriv[*par][pScattPow].size()==0)
2246  {
2247  cout<<" => skipped (deriv==0)"<<endl;
2248  continue;
2249  }
2250  else cout <<endl;
2251  const REAL step=(*par)->GetDerivStep();
2252  (*par)->Mutate(step);
2253  this->CalcGeomStructFactor();
2254  mr[make_pair(pScattPow,*par)]=mvRealGeomSF[pScattPow];
2255  mi[make_pair(pScattPow,*par)]=mvImagGeomSF[pScattPow];
2256  (*par)->Mutate(-2*step);
2257  this->CalcGeomStructFactor();
2258  mr[make_pair(pScattPow,*par)]-=mvRealGeomSF[pScattPow];
2259  mr[make_pair(pScattPow,*par)]/=step*2;
2260  mi[make_pair(pScattPow,*par)]-=mvImagGeomSF[pScattPow];
2261  mi[make_pair(pScattPow,*par)]/=step*2;
2262  (*par)->Mutate(step);
2263 
2264  v.push_back(&(mvRealGeomSF_FullDeriv[*par][pScattPow]));
2265  v.push_back(&(mr[make_pair(pScattPow,*par)]));
2266  if(!hasinv)
2267  {
2268  v.push_back(&(mvImagGeomSF_FullDeriv[*par][pScattPow]));
2269  v.push_back(&(mi[make_pair(pScattPow,*par)]));
2270  }
2271  if(v.size()>20)break;
2272  }
2273  if(v.size()>20)break;
2274  }
2275  cout<<"############################ Geom Fhkl Deriv Real, Imag ##############################"
2276  <<endl<<FormatVertVectorHKLFloats<REAL>(v,11,4,20)
2277  <<"############################ END GeomF hkl Deriv Real, Imag ##############################"<<endl;
2278  cout<<"############################ X Y Z Occ ##############################"
2279  <<endl<<FormatVertVector<REAL>(vdx[*(vPar.begin())],vdy[*(vPar.begin())],vdz[*(vPar.begin())],vdocc[*(vPar.begin())],8,4,20)<<endl
2280  <<"############################ END X Y Z Occ ##############################"<<endl;
2281  //exit(0);
2282  #endif
2283 
2284  // We can use geom struct fact calculated at the beginning, since parameters are back to the same values.
2286 }
2287 
2289 {
2290  // Assume this is called by ScatteringData::CalcStructFactor()
2291  // and that we already have computed geometrical structure factors
2292  VFN_DEBUG_ENTRY("ScatteringData::CalcLuzzatiFactor",3)
2293  bool useLuzzati=false;
2294  for(map<const ScatteringPower*,CrystVector_REAL>::const_iterator
2295  pos=mvRealGeomSF.begin();pos!=mvRealGeomSF.end();++pos)
2296  {
2297  if(pos->first->GetMaximumLikelihoodPositionError()!=0)
2298  {
2299  useLuzzati=true;
2300  break;
2301  }
2302  }
2303  if(!useLuzzati)
2304  {
2305  mvLuzzatiFactor.clear();
2306  VFN_DEBUG_EXIT("ScatteringData::CalcLuzzatiFactor(): not needed, no positionnal errors",3)
2307  return;
2308  }
2309  bool recalc=false;
2310  if( (mClockTheta >mClockLuzzatiFactor)
2311  ||(mClockGeomStructFact>mClockLuzzatiFactor)//checks if occupancies changed
2312  ||(mClockNbReflUsed>mClockLuzzatiFactor))
2313  {
2314  //if(mClockTheta >mClockLuzzatiFactor)cout<<"1"<<endl;
2315  //if(mClockGeomStructFact>mClockLuzzatiFactor)cout<<"2"<<endl;
2316  //if(mClockNbReflUsed>mClockLuzzatiFactor)cout<<"3"<<endl;
2317  recalc=true;
2318  }
2319  else
2320  {
2321  for(int i=mpCrystal->GetScatteringPowerRegistry().GetNb()-1;i>=0;i--)
2322  {
2323  if(mpCrystal->GetScatteringPowerRegistry().GetObj(i)
2324  .GetMaximumLikelihoodParClock()>mClockLuzzatiFactor)
2325  {
2326  recalc=true;
2327  break;
2328  }
2329  }
2330  }
2331 
2332  if(false==recalc)
2333  {
2334  VFN_DEBUG_EXIT("ScatteringData::CalcLuzzatiFactor(): no recalc needed",3)
2335  return;
2336  }
2337  TAU_PROFILE("ScatteringData::CalcLuzzatiFactor()","void ()",TAU_DEFAULT);
2338 
2339  for(int i=mpCrystal->GetScatteringPowerRegistry().GetNb()-1;i>=0;i--)
2340  {
2341  const ScatteringPower* pScattPow=&(mpCrystal->GetScatteringPowerRegistry().GetObj(i));
2342  if(0 == pScattPow->GetMaximumLikelihoodPositionError())
2343  {
2344  mvLuzzatiFactor[pScattPow].resize(0);
2345  }
2346  else
2347  {
2348  mvLuzzatiFactor[pScattPow].resize(mNbRefl);
2349  const REAL b=-(8*M_PI*M_PI)* pScattPow->GetMaximumLikelihoodPositionError()
2350  * pScattPow->GetMaximumLikelihoodPositionError();
2351  const REAL *stol=this->GetSinThetaOverLambda().data();
2352  REAL *fact=mvLuzzatiFactor[pScattPow].data();
2353  for(long j=0;j<mNbReflUsed;j++) {*fact++ = exp(b * *stol * *stol);stol++;}
2354  VFN_DEBUG_MESSAGE("ScatteringData::CalcLuzzatiFactor():"<<pScattPow->GetName()<<endl<<
2356  mvRealGeomSF[pScattPow],mvImagGeomSF[pScattPow],
2357  mvScatteringFactor[pScattPow],mvLuzzatiFactor[pScattPow],10,4,mNbReflUsed
2358  ),2);
2359  }
2360  }
2361  mClockLuzzatiFactor.Click();
2362  VFN_DEBUG_EXIT("ScatteringData::CalcLuzzatiFactor(): no recalc needed",3)
2363 }
2364 
2366 {
2367  // this is called by CalcStructFactor(), after the calculation of the structure factors,
2368  // and the recomputation of Luzzati factors has already been asked
2369  // So we only recompute if these clocks have changed.
2370  //
2371  // The Crystal::mMasterClockScatteringPower will tell the last time the number of ghost
2372  // atoms has been changed in any of the scattpow.
2373 
2374  if( (mClockFhklCalcVariance>mClockLuzzatiFactor)
2375  &&(mClockFhklCalcVariance>mClockStructFactor)
2376  &&(mClockFhklCalcVariance>mpCrystal->GetMasterClockScatteringPower())) return;
2377 
2378  bool hasGhostAtoms=false;
2379  for(map<const ScatteringPower*,CrystVector_REAL>::const_iterator
2380  pos=mvRealGeomSF.begin();pos!=mvRealGeomSF.end();++pos)
2381  {
2382  if(pos->first->GetMaximumLikelihoodNbGhostAtom()!=0)
2383  {
2384  hasGhostAtoms=true;
2385  break;
2386  }
2387  }
2388 
2389  if( (0==mvLuzzatiFactor.size())&&(!hasGhostAtoms))
2390  {
2391  mFhklCalcVariance.resize(0);
2392  return;
2393  }
2394  VFN_DEBUG_ENTRY("ScatteringData::CalcStructFactVariance()",3)
2395  TAU_PROFILE("ScatteringData::CalcStructFactVariance()","void ()",TAU_DEFAULT);
2396  bool needVar=false;
2397 
2398  map<const ScatteringPower*,REAL> vComp;
2399  {
2401  const long nbComp=pList->GetNbComponent();
2402  const ScatteringComponent *pComp;
2403  for(long i=0;i<nbComp;i++)
2404  {
2405  pComp=&((*pList)(i));
2406  vComp[pComp->mpScattPow]=0;
2407  }
2408  for(long i=0;i<nbComp;i++)
2409  {
2410  pComp=&((*pList)(i));
2411  vComp[pComp->mpScattPow]+= pComp->mOccupancy * pComp->mDynPopCorr;
2412  }
2413  for(map<const ScatteringPower*,REAL>::iterator
2414  pos=vComp.begin();pos!=vComp.end();++pos)
2415  pos->second *= this->GetCrystal().GetSpaceGroup().GetNbSymmetrics();
2416  }
2417  // Ghost atoms
2418  map<const ScatteringPower*,REAL> vGhost;
2419  {
2420  const long nbScattPow=mpCrystal->GetScatteringPowerRegistry().GetNb();
2421  const long mult=this->GetCrystal().GetSpaceGroup().GetNbSymmetrics();
2422  for(int i=0;i<nbScattPow;i++)
2423  {
2424  const ScatteringPower* pow=&(mpCrystal->GetScatteringPowerRegistry().GetObj(i));
2425  const REAL nb=pow->GetMaximumLikelihoodNbGhostAtom();
2426  vGhost[pow]=nb*mult;
2427  }
2428  }
2429 
2430  if(mFhklCalcVariance.numElements() == mNbRefl)
2431  {
2432  REAL *pVar=mFhklCalcVariance.data();
2433  for(long j=0;j<mNbReflUsed;j++) *pVar++ = 0;
2434  }
2435 
2436  for(int i=mpCrystal->GetScatteringPowerRegistry().GetNb()-1;i>=0;i--)
2437  {
2438  const ScatteringPower* pScattPow=&(mpCrystal->GetScatteringPowerRegistry().GetObj(i));
2439  if( (mvLuzzatiFactor[pScattPow].numElements()==0)
2440  &&(vGhost[pScattPow]==0)) continue;
2441  needVar=true;
2442  if(mFhklCalcVariance.numElements() != mNbRefl)
2443  {
2444  mFhklCalcVariance.resize(mNbRefl);
2445  REAL *pVar=mFhklCalcVariance.data();
2446  for(long j=0;j<mNbReflUsed;j++) *pVar++ = 0;
2447  }
2448  // variance on real & imag parts of the structure factor
2449  const REAL *pScatt=mvScatteringFactor[pScattPow].data();
2450  const int *pExp=mExpectedIntensityFactor.data();
2451  REAL *pVar=mFhklCalcVariance.data();
2452  if(mvLuzzatiFactor[pScattPow].numElements()==0)
2453  {
2454  const REAL nbghost=vGhost[pScattPow];
2455  for(long j=0;j<mNbReflUsed;j++)
2456  {
2457  *pVar++ += *pExp++ * *pScatt * *pScatt * nbghost;
2458  pScatt++;
2459  }
2460  }
2461  else
2462  {
2463  const REAL *pLuz=mvLuzzatiFactor[pScattPow].data();
2464  const REAL occ=vComp[pScattPow];
2465  const REAL nbghost=vGhost[pScattPow];
2466  for(long j=0;j<mNbReflUsed;j++)
2467  {
2468  *pVar++ += *pExp++ * *pScatt * *pScatt * ( occ*(1 - *pLuz * *pLuz) + nbghost);
2469  pScatt++; pLuz++;
2470  }
2471  }
2472  }
2473  if(false == needVar) mFhklCalcVariance.resize(0);
2474 
2475  mClockFhklCalcVariance.Click();
2476  VFN_DEBUG_EXIT("ScatteringData::CalcStructFactVariance()",3)
2477 }
2478 
2479 }//namespace ObjCryst
CrystVector_REAL GetLatticePar() const
Lattice parameters (a,b,c,alpha,beta,gamma) as a 6-element vector in Angstroems and radians...
Definition: UnitCell.cpp:92
const RefinableObjClock & GetClockNbReflBelowMaxSinThetaOvLambda() const
Clock the last time the number of reflections used was changed.
const CrystVector_REAL & GetH() const
Return the 1D array of H coordinates for all reflections.
void Print() const
Print to screen/console the charcteristics of the radiation.
const RefinableObjClock & GetMasterClockScatteringPower() const
Get the clock which reports all changes in ScatteringPowers.
Definition: Crystal.cpp:280
const CrystVector_REAL & GetFhklCalcImag() const
Access to imaginary part of F(hkl)calc.
map< const ScatteringPower *, CrystVector_REAL > mvRealGeomSF
Geometrical Structure factor for each ScatteringPower, as vectors with NbRefl elements.
REAL GetMaximumLikelihoodNbGhostAtom() const
Maximum Likelihood: get the number of ghost elements per asymmetric unit.
virtual long GetNbReflBelowMaxSinThetaOvLambda() const
Recalc, and get the number of reflections which should be actually used, due to the maximuml sin(thet...
RefinableObjClock mClockGeomStructFact
Clock the last time the geometrical structure factors were computed.
map< const ScatteringPower *, REAL > mvFprime
Anomalous X-Ray scattering term f' and f" are stored here for each ScatteringPower We store here only...
bool mIgnoreImagScattFact
Ignore imaginary part of scattering factor.
REAL mXRayTubeDeltaLambda
Absolute difference between alpha1 and alpha2, in angstroems.
RefinableObjClock mClockHKL
Clock for the list of hkl.
void AddSubRefObj(RefinableObj &)
void AddPar(const RefinablePar &newRefPar)
Add a refinable parameter.
RadiationType GetRadiationType() const
Neutron or x-ray experiment ? Wavelength ?
int GetNbSymmetrics(const bool noCenter=false, const bool noTransl=false) const
Return the number of equivalent positions in the spacegroup, ie the multilicity of the general positi...
Definition: SpaceGroup.cpp:419
const CrystMatrix_REAL & GetBMatrix() const
Get the 'B' matrix (UnitCell::mBMatrix)for the UnitCell (orthogonalization matrix for the given latti...
Definition: UnitCell.cpp:235
const CrystVector_REAL & GetSinThetaOverLambda() const
Return an array with for all reflections.
RefObjOpt mWavelengthType
monochromatic ? Alpha1 & Alpha2 ? Multi-Wavelength ?
const RefinableObjClock & GetClockLatticePar() const
last time the Lattice parameters were changed
Definition: UnitCell.cpp:323
void CalcGlobalTemperatureFactor() const
Compute the overall temperature factor affecting all reflections.
virtual const ScatteringComponentList & GetScatteringComponentList() const
Get the list of all scattering components.
Definition: Crystal.cpp:283
const RefParType * gpRefParTypeScattDataCorrInt
Generic type for correction to calculated intensities.
const cctbx::sgtbx::space_group & GetCCTbxSpg() const
Get the underlying cctbx Spacegroup object.
Definition: SpaceGroup.cpp:465
RefinableObjClock mClockNbReflUsed
Clock recording the last time the number of reflections used has increased.
const CrystVector_REAL & GetTheta() const
Return an array with theta values for all reflections.
virtual CrystMatrix_REAL GetResonantScattFactReal(const ScatteringData &data, const int spgSymPosIndex=-1) const =0
Get the real part of the resonant scattering factor.
virtual void PrepareHKLarrays()
const RefParType * gpRefParTypeScattDataCorrPos
Parameter type for correction to peak positions.
const RefinableObjClock & GetClockSpaceGroup() const
Get the SpaceGroup Clock (corresponding to the time of the initialization of the SpaceGroup) ...
Definition: SpaceGroup.cpp:467
Output vectors as column arrays, with the first 3 columns printed as integers.
RefinableObjClock mClockScattFactor
Clock the last time scattering factors were computed.
CrystVector_REAL mGlobalTemperatureFactor
Global Biso factor.
A scattering position in a crystal, associated with the corresponding occupancy and a pointer to the ...
REAL mMaxSinThetaOvLambda
Maximum sin(theta)/lambda for all calculations (10 by default).
long GetNbComponent() const
Number of components.
We need to record exactly when refinable objects have been modified for the last time (to avoid re-co...
Definition: RefinableObj.h:138
bool IsBeingRefined() const
Is the object being refined ? (Can be refined by one algorithm at a time only.)
virtual void CalcSinThetaLambda() const
virtual void BeginOptimization(const bool allowApproximations=false, const bool enableRestraints=false)
This should be called by any optimization class at the begining of an optimization.
virtual void BeginOptimization(const bool allowApproximations=false, const bool enableRestraints=false)
This should be called by any optimization class at the begining of an optimization.
CrystVector_REAL GetWavelength() const
wavelength of the experiment (in Angstroems)
REAL mLinearPolarRate
Linear Polarization Rate (default:0, X-Ray tube unmonochromatized)
REAL GetMaxSinThetaOvLambda() const
Get the maximum value for sin(theta)/lambda.
CrystVector_REAL mFhklCalcVariance
The variance on all calculated structure factors, taking into account the positionnal errors and the ...
virtual void CalcResonantScattFactor() const
bool HasCrystal() const
Has a Crystal structure associated yet ?
REAL mXRayTubeAlpha2Alpha1Ratio
Ratio alpha2/alpha1 (should be 0.5)
REAL GetMaximumLikelihoodPositionError() const
Maximum Likelihood: get the estimated error (sigma) on the positions for this kind of element...
const RefinableObjClock & GetClockRadiation() const
Last time the nature (X-Rays/Neutron, number of wavelengths)radiation has been changed.
void Click()
Record an event for this clock (generally, the 'time' an object has been modified, or some computation has been made)
bool HasInversionCenter() const
Is centrosymmetric ?
Definition: SpaceGroup.cpp:463
void AddChild(const RefinableObjClock &)
Add a 'child' clock.
RefinableObjClock mClockStructFactor
Clock for the structure factor.
const RefParType * gpRefParTypeScattDataCorrInt_Ellipsoid
Parameter type for the ellipsoid coefficient.
string mXRayTubeName
Name of the X-Ray tube used, if relevant.
long mNbRefl
Number of H,K,L reflections.
virtual void PrintFhklCalcDetail(ostream &os=cout) const
Print H, K, L sin(theta)/lambda theta F^2 Re(F) Im(F) [Re(F) Im(F)]_i, where [Re(F) Im(F)]_i are the ...
virtual void SetCrystal(Crystal &crystal)
Set the crystal for this experiment.
CrystVector_REAL mWavelength
Wavelength of the Experiment, in Angstroems.
const RefParType * gpRefParTypeScattDataProfile
Type for reflection profile.
map< const ScatteringPower *, CrystVector_REAL > mvTemperatureFactor
Thermic factors for each ScatteringPower, as vectors with NbRefl elements.
const RefParType * gpRefParTypeScattDataProfileAsym
Type for reflection profile asymmetry.
RadiationType GetRadiationType() const
Get the radiation type (X-Rays, Neutron)
virtual void RegisterClient(RefinableObj &) const
Register a new object using this object.
CrystVector_REAL mH2Pi
H,K,L coordinates, multiplied by 2PI.
Class to compute structure factors for a set of reflections and a Crystal.
RefinableObjClock mClockMaster
Master clock, which is changed whenever the object has been altered.
void GetSymmetric(unsigned int i, REAL &x, REAL &y, REAL &z, const bool noCenter=false, const bool noTransl=false, const bool derivative=false) const
Get all equivalent positions of a (xyz) position.
Definition: SpaceGroup.cpp:366
virtual const Radiation & GetRadiation() const =0
Get the radiation object for this data.
virtual void SetApproximationFlag(const bool allow)
Enable or disable numerical approximations.
const RefinableObjClock & GetClockWavelength() const
Last time the wavelength has been changed.
Radiation()
Default constructor.
const CrystVector_REAL & GetL() const
Return the 1D array of L coordinates for all reflections.
REAL mDynPopCorr
Dynamical Population Correction.
const RefParType * gpRefParTypeScattDataScale
Type for scattering data scale factors.
CrystVector_int mMultiplicity
Multiplicity for each reflections (mostly for powder diffraction)
virtual CrystVector_REAL GetTemperatureFactor(const ScatteringData &data, const int spgSymPosIndex=-1) const =0
Get the temperature factor for all reflections of a given ScatteringData object.
virtual void GenHKLFullSpace2(const REAL maxsithsl, const bool unique=false)
Generate a list of h,k,l to describe a full reciprocal space, up to a given maximum theta value...
const SpaceGroup & GetSpaceGroup() const
Access to the SpaceGroup object.
Definition: UnitCell.cpp:320
RefinableObjClock mClockTheta
Clock the last time theta was computed.
void CalcGeomStructFactor() const
Compute the 'Geometrical Structure Factor' for each ScatteringPower of the Crystal.
map< const ScatteringPower *, CrystVector_REAL > mvLuzzatiFactor
The Luzzati 'D' factor for each scattering power and each reflection.
const CrystVector_REAL & GetReflZ() const
Return the 1D array of orthonormal z coordinates for all reflections (recipr. space) ...
virtual void GenHKLFullSpace(const REAL maxTheta, const bool unique=false)
Generate a list of h,k,l to describe a full reciprocal space, up to a given maximum theta value...
virtual CrystVector_long SortReflectionBySinThetaOverLambda(const REAL maxSTOL=-1.)
const RefParType * gpRefParTypeScattDataBackground
Parameter type for background intensity.
REAL GetXRayTubeDeltaLambda() const
Get the wavelength difference for Alpha1 and Alpha2.
int GetNbTranslationVectors() const
Number of translation vectors (1 for 'P' cells, 2 for 'I', 4 for 'F',etc..)
Definition: SpaceGroup.cpp:259
CrystVector_REAL mTheta
theta for the crystal and the HKL in ReciprSpace (in radians)
const RefParType * gpRefParTypeScattDataCorrIntPO_Amplitude
Parameter type for the amplitude of preferred orientation.
RefinableObjClock mClockThermicFact
Clock the last time temperature factors were computed.
The crystallographic space group, and the cell choice.
Definition: SpaceGroup.h:104
const ScatteringPower * mpScattPow
The ScatteringPower associated with this position.
void CalcStructFactor() const
Compute the overall structure factor (real and imaginary part).
CrystVector_REAL mH
H,K,L coordinates.
virtual const CrystMatrix_REAL & GetBMatrix() const
Get access to the B matrix used to compute reflection positions.
const map< const ScatteringPower *, CrystVector_REAL > & GetScatteringFactor() const
Scattering factors for each ScatteringPower, as vectors with NbRefl elements.
void CalcLuzzatiFactor() const
Calculate the Luzzati factor associated to each ScatteringPower and each reflection, for maximum likelihood optimization.
const RefinableObjClock & GetClockScattCompList() const
Get the list of all scattering components.
Definition: Crystal.cpp:319
WavelengthType GetWavelengthType() const
Get the Wavelength type (monochromatic, Alpha1+Alpha2, Time Of Flight...)
CrystVector_long mIntH
H,K,L integer coordinates.
CrystVector_REAL mSinThetaLambda
for the crystal and the reflections in ReciprSpace
virtual CrystMatrix_REAL GetResonantScattFactImag(const ScatteringData &data, const int spgSymPosIndex=-1) const =0
Get the imaginary part of the resonant scattering factor.
const CrystVector_REAL & GetK2Pi() const
Return the 1D array of K coordinates for all reflections, multiplied by 2*pi.
void SetIsIgnoringImagScattFact(const bool b)
If true, then the imaginary part of the scattering factor is ignored during Structure factor computat...
void SetWavelength(const REAL)
Set the (monochromatic) wavelength of the beam.
const CrystVector_REAL & GetReflY() const
Return the 1D array of orthonormal y coordinates for all reflections (recipr. space) ...
const CrystVector_REAL & GetFhklCalcSq() const
Returns the Array of calculated |F(hkl)|^2 for all reflections.
CrystMatrix_REAL GetAllSymmetrics(const REAL x, const REAL y, const REAL z, const bool noCenter=false, const bool noTransl=false, const bool noIdentical=false) const
Get all equivalent positions of a (xyz) position.
Definition: SpaceGroup.cpp:274
const CrystVector_REAL & GetReflX() const
Return the 1D array of orthonormal x coordinates for all reflections (recipr. space) ...
virtual void SetMaxSinThetaOvLambda(const REAL max)
Set the maximum value for sin(theta)/lambda.
WavelengthType
Incident beam characteristics : monochromatic, X-Ray tube with Alpha1 and alpha2, MAD (a few waveleng...
Definition: General.h:100
const CrystVector_REAL & GetFhklObsSq() const
Returns the Array of observed |F(hkl)|^2 for all reflections.
RefinableObjClock mClockScattFactorResonant
Clock the last time resonant scattering factors were computed.
RefObjOpt mRadiationType
Neutron ? X-Ray ? (Electron: unimplemented)
const RefParType * gpRefParTypeScattDataProfileWidth
Type for reflection profile width.
REAL mGlobalBiso
Global Biso, affecting the overall structure factor for all reflections (but not the structure factor...
const RefParType * gpRefParTypeScattDataCorrIntAbsorp
Parameter type for absorption correction.
virtual const string & GetClassName() const
Name for this class ("RefinableObj", "Crystal",...).
void SetRadiationType(const RadiationType)
Set the radiation type (X-Rays, Neutron)
CrystVector_REAL mFhklObsSq
Observed squared structure factors (zero-sized if none)
const CrystVector_REAL & GetWavelength() const
Get the wavelength(s) in Angstroems.
ObjRegistry< ScatteringPower > & GetScatteringPowerRegistry()
Get the registry of ScatteringPower included in this Crystal.
Definition: Crystal.cpp:222
const Crystal & GetCrystal() const
Const access to the data's crystal.
output a string with a fixed length (adding necessary space or removing excess characters) : ...
const RefParType * gpRefParTypeScattDataProfileType
Type for reflection profiles type (e.g. gaussian/lorentzian mix)
const CrystVector_REAL & GetFhklCalcReal() const
Access to real part of F(hkl)calc.
CrystVector_REAL mFhklCalcSq
F(HKL)^2 calc for each reflection.
const RefParType * gpRefParTypeScattDataCorrIntPolar
Parameter type for polarization correction.
long mNbReflUsed
Number of reflections which are below the max.
bool IsInversionCenterAtOrigin() const
Is the center of symmetry at the origin ?
Definition: SpaceGroup.cpp:464
CrystVector_long EliminateExtinctReflections()
const CrystVector_REAL & GetL2Pi() const
Return the 1D array of L coordinates for all reflections, multiplied by 2*pi.
unsigned int GetExpectedIntensityFactor(const REAL h, const REAL k, const REAL l) const
Get the "expected intensity factor" for a given reflection.
Definition: SpaceGroup.cpp:562
map< const ScatteringPower *, CrystVector_REAL > mvScatteringFactor
Scattering factors for each ScatteringPower, as vectors with NbRefl elements.
Crystal * mpCrystal
Pointer to the crystal corresponding to this experiment.
const RefParType * gpRefParTypeScattDataCorr
Generic type for scattering data correction parameter.
void Reset()
Reset a Clock to 0, to force an update.
CrystVector_int mExpectedIntensityFactor
Expected intensity factor for all reflections.
void CalcTemperatureFactor() const
Exception class for ObjCryst++ library.
Definition: General.h:119
virtual void EndOptimization()
This should be called by any optimization class at the end of an optimization.
long GetNbRefl() const
Return the number of reflections in this experiment.
const RefinableObjClock & GetClockTheta() const
Clock the last time the sin(theta)/lambda and theta arrays were re-computed.
The namespace which includes all objects (crystallographic and algorithmic) in ObjCryst++.
Definition: Atom.cpp:47
REAL GetXRayTubeAlpha2Alpha1Ratio() const
Get the Kalpha2/Kalpha1 ratio.
Generic class for parameters of refinable objects.
Definition: RefinableObj.h:223
const RefParType * gpRefParTypeScattDataCorrIntPO_Direction
Parameter type for preferred orientation direction.
RefinableObjClock mClockStructFactorSq
Clock for the square modulus of the structure factor.
RefinableObjClock mClockFhklObsSq
Last time observed squared structure factors were altered.
RefinableObjClock mClockGlobalBiso
last time the global Biso factor was modified
virtual const string & GetName() const
Name of the object.
virtual void PrintFhklCalc(ostream &os=cout) const
Print H, K, L F^2 Re(F) Im(F) theta sin(theta)/lambda for all reflections.
bool mUseFastLessPreciseFunc
Use faster, but less precise, approximations for functions? (integer approximations to compute sin an...
int GetSpaceGroupNumber() const
Id number of the spacegroup.
Definition: SpaceGroup.cpp:249
RefinableObjClock mClockGlobalTemperatureFact
last time the global temperature factor was computed
virtual void SetApproximationFlag(const bool allow)
Enable or disable numerical approximations.
virtual void SetHKL(const CrystVector_REAL &h, const CrystVector_REAL &k, const CrystVector_REAL &l)
input H,K,L
const RefParType * gpRefParTypeScattDataCorrIntExtinc
Parameter type for extinction correction.
RadiationType
Type of radiation used.
Definition: General.h:94
virtual void EndOptimization()
This should be called by any optimization class at the end of an optimization.
const RefParType * gpRefParTypeScattDataCorrIntPO_Fraction
Parameter type for fraction of preferred orientation.
Crystal class: Unit cell, spacegroup, scatterers.
Definition: Crystal.h:97
void SetWavelengthType(const WavelengthType &type)
Set the Wavelength type (monochromatic, Alpha1+Alpha2, Time Of Flight...)
virtual CrystVector_REAL GetScatteringFactor(const ScatteringData &data, const int spgSymPosIndex=-1) const =0
Get the Scattering factor for all reflections of a given ScatteringData object.
class of refinable parameter types.
Definition: RefinableObj.h:78
list of scattering positions in a crystal, associated with the corresponding occupancy and a pointer ...
void CalcStructFactVariance() const
Calculate the variance associated to the calculated structure factor.
const CrystVector_REAL & GetK() const
Return the 1D array of K coordinates for all reflections.
virtual void DeRegisterClient(RefinableObj &) const
Deregister an object (which not any more) using this object.
const RefParType * gpRefParTypeScattData
Generic type for scattering data.
CrystVector_REAL mFhklCalcReal
real &imaginary parts of F(HKL)calc
Class to define the radiation (type, monochromaticity, wavelength(s)) of an experiment.
const CrystVector_REAL & GetH2Pi() const
Return the 1D array of H coordinates for all reflections, multiplied by 2*pi.
int mOptimizationDepth
Is the object being refined or optimized ? if mOptimizationDepth=0, no optimization is taking place...
void AddOption(RefObjOpt *opt)
bool IsIgnoringImagScattFact() const
If true, then the imaginary part of the scattering factor is ignored during Structure factor computat...
CrystVector_REAL mX
reflection coordinates in an orthonormal base
const std::vector< SpaceGroup::TRx > & GetTranslationVectors() const
Return all Translation Vectors, as a 3 columns-array.
Definition: SpaceGroup.cpp:264
Abstract Base Class to describe the scattering power of any Scatterer component in a crystal...