OpenAtom  Version1.5a
CP_State_Plane.h
Go to the documentation of this file.
1 //////////////////////////////////////////////////////////////////////////////
2 //////////////////////////////////////////////////////////////////////////////
3 //////////////////////////////////////////////////////////////////////////////
4 /** \file CP_State_Plane.h
5  *
6  */
7 //////////////////////////////////////////////////////////////////////////////
8 
9 #ifndef _PLANE_H_
10 #define _PLANE_H_
11 
12 //////////////////////////////////////////////////////////////////////////////
13 
14 #include "charm++.h"
15 #include "debug_flags.h"
16 #include "ckmulticast.h"
17 #include "structure_factor/StructFactorCache.h"
18 #include "structure_factor/StructureFactor.h"
19 //#include "ckPairCalculator.h"
20 void getSplitDecomp(int *,int *,int *,int , int ,int );
21 
22 
23 
24 
25 //////////////////////////////////////////////////////////////////////////////
26 //////////////////////////////////////////////////////////////////////////////
27 //////////////////////////////////////////////////////////////////////////////
29  public:
30  int datalen, hops;
31  int subplane;
32  double *data;
33  int idx;
34 };
35 //////////////////////////////////////////////////////////////////////////////
36 //////////////////////////////////////////////////////////////////////////////
37 //////////////////////////////////////////////////////////////////////////////
38 //////////////////////////////////////////////////////////////////////////////
40  public:
41  int nZmat;
42  int iterNL;
43  double *zmat;
44  void init(int _nzmat, double *_zmat, int _iterNL)
45  {
46  nZmat = _nzmat;
47  iterNL = _iterNL;
48  CmiMemcpy(zmat, _zmat, sizeof(double)* nZmat);
49  }
50  friend class CMessage_CompAtmForcMsg;
51 };
52 //////////////////////////////////////////////////////////////////////////////
53 
54 
55 //////////////////////////////////////////////////////////////////////////////
56 //////////////////////////////////////////////////////////////////////////////
57 //////////////////////////////////////////////////////////////////////////////
59 };
60 //////////////////////////////////////////////////////////////////////////////
61 
62 
63 //////////////////////////////////////////////////////////////////////////////
64 //////////////////////////////////////////////////////////////////////////////
65 //////////////////////////////////////////////////////////////////////////////
67 };
68 //////////////////////////////////////////////////////////////////////////////
69 
70 
71 
72 //////////////////////////////////////////////////////////////////////////////
73 //////////////////////////////////////////////////////////////////////////////
74 //////////////////////////////////////////////////////////////////////////////
75 class TMsg: public CMessage_TMsg {
76 public:
77  int datalen;
78  complex *data;
79 };
80 //////////////////////////////////////////////////////////////////////////////
81 
82 
83 //////////////////////////////////////////////////////////////////////////////
84 //////////////////////////////////////////////////////////////////////////////
85 //////////////////////////////////////////////////////////////////////////////
87 public:
88  int size;
89  int offset;
90  int offsetGx;
91  int iopt;
92  int num;
93  complex *data;
94 };
95 //////////////////////////////////////////////////////////////////////////////
96 
97 //////////////////////////////////////////////////////////////////////////////
98 //////////////////////////////////////////////////////////////////////////////
99 //////////////////////////////////////////////////////////////////////////////
101 public:
102  int size;
103  int senderIndex;
104  int offset;
105  int offsetGx;
106  int iter;
107  int iopt;
108  int num;
109  complex *data;
110 };
111 //////////////////////////////////////////////////////////////////////////////
112 
113 //////////////////////////////////////////////////////////////////////////////
114 //////////////////////////////////////////////////////////////////////////////
115 //////////////////////////////////////////////////////////////////////////////
117 public:
118  int size;
119  int senderIndex;
120  int iopt;
121  int iter;
122  complex *data;
123 };
124 //////////////////////////////////////////////////////////////////////////////
125 
126 
127 
128 //////////////////////////////////////////////////////////////////////////////
129 //////////////////////////////////////////////////////////////////////////////
130 //////////////////////////////////////////////////////////////////////////////
132 public:
133  int size;
134  int senderIndex;
135  int iopt;
136  complex *data;
137 };
138 //////////////////////////////////////////////////////////////////////////////
139 
140 //////////////////////////////////////////////////////////////////////////////
141 //////////////////////////////////////////////////////////////////////////////
142 //////////////////////////////////////////////////////////////////////////////
144 public:
145  int size;
146  int index;
147  int senderBigIndex;
148  int senderStrtLine;
149  int iopt;
150  complex *data;
151 };
152 //////////////////////////////////////////////////////////////////////////////
153 
154 
155 //////////////////////////////////////////////////////////////////////////////
156 //////////////////////////////////////////////////////////////////////////////
157 //////////////////////////////////////////////////////////////////////////////
159 public:
160  int size;
161  int senderIndex;
162  int numPlanes;
163  complex *data;
164 };
165 //////////////////////////////////////////////////////////////////////////////
166 
167 /** @addtogroup RealSpaceState
168  * @{
169  *
170  * \brief Chare Array 2D chare array nplanex X nstates. Handles
171  * electronic structure in real space. The points of plane-wave
172  * pseudo-potential are cut along the x-dimension for finer
173  * parallelization.
174 */
176  public:
177  CP_State_RealSpacePlane_SDAG_CODE
178  CP_State_RealSpacePlane(int, int,int,int,int,int,int, UberCollection);
179  CP_State_RealSpacePlane(CkMigrateMessage *m) {};
180  ~CP_State_RealSpacePlane() { if(cookie!=NULL) delete [] cookie; };
181  void unpackFFT(RSFFTMsg *);
182  void doFFT();
183  void doVksFFT();
184  void unpackProduct(ProductMsg *);
185  void doProductThenFFT();
186  void sendFPsiToGSP();
187  void setNumPlanesToExpect(int num);
188  void printData();
189  void init(ProductMsg *);
190  void doReduction();
191  void ResumeFromSync();
192  void pup(PUP::er &);
193  private:
194  const UberCollection thisInstance;
195  int forwardTimeKeep;
196  int backwardTimeKeep;
197  int iplane_ind;
198  int ibead_ind,kpoint_ind, itemper_ind;
199  int iteration;
200  int rhoRsubplanes;
201  int ngrida;
202  int ngridb;
203  int ngridc;
204  int count;
205  int rsize;
206  int csize;
207  int countProduct;
208  int numCookies;
209  int istate;
210  UberCollection RhoReductionDest;
211  RealStateSlab rs;
212  CkSectionInfo *cookie;
213  CProxy_CP_State_GSpacePlane gproxy;
214 };
215 //////////////////////////////////////////////////////////////////////////////
216 /*@}*/
217 
218 /** @addtogroup Density
219  @{
220 */
221 //////////////////////////////////////////////////////////////////////////////
222 //////////////////////////////////////////////////////////////////////////////
223 //////////////////////////////////////////////////////////////////////////////
225  public:
226  int listSubFlag;
227  int cp_grad_corr_on;
228  int ees_eext_on;
229  int nplane_rho_x;
230  int ngridcEext;
231  int ngrida;
232  int ngridb;
233  int ngridc;
234  int iplane_ind;
235  int myNgridb; // I don't have all ngridb lines of x now
236  int myNplane_rho; // I don't have all nplane_rho_x lines of y now
237  int nptsExpndA;
238  int nptsExpndB;
239  int nptsA;
240  int nptsB;
241  int myBoff;
242  int myAoff;
243  int countDebug;
244  int countFFTRyToGy;
245  int rhoRsubplanes;
246  int myTime;
247  bool doneRhoReal;
248  bool doneRHart;
249  int recvCountFromGRho;
250  int recvCountFromGHartExt;
251  CP_Rho_RealSpacePlane(CkMigrateMessage *m){}
252  CP_Rho_RealSpacePlane(int, bool,int,int,int, UberCollection);
253  void init();
255  void pup(PUP::er &);
256  void acceptDensity(CkReductionMsg *);
257  void handleDensityReduction();
258  void launchEextRNlG();
259  void energyComputation();
260  void fftRhoRtoRhoG();
261  void launchNLRealFFT();
262  void sendPartlyFFTRyToGy(int iopt);
264  void fftRhoRyToGy(int iopt);
265  void sendPartlyFFTtoRhoG(int );
267  void sendPartlyFFTGxToRx(int );
269  void GradCorr();
270  void whiteByrdFFT();
271  void acceptWhiteByrd(RhoRSFFTMsg *msg);
272  void addWhiteByrdVks();
274  void addHartEextVks();
275  void RHartReport();
276  void doMulticastCheck();
277  void doMulticast();
278  void exitForDebugging();
279  void isAtSync(int iter){AtSync();};
280  void ResumeFromSync();
281  void sendPartlyFFTtoRhoGall();
282  void acceptGradRhoVksAll(RhoRSFFTMsg *msg);
283 
284  private:
285  const UberCollection thisInstance;
286  int rhoKeeperId;
287  int rhoGHelpers;
288  int countGradVks[5]; // number of collections that have arrived
289  int countIntRtoG[5]; // our internal transpose friends.
290  int countIntGtoR[5]; // our internal transpose friends.
291  int countWhiteByrd; // number of collections that have arrived
292  int countRHart;
293  int countRHartValue;
294  int doneGradRhoVks; // count 1,2,3 = x,y,z all done
295  bool doneWhiteByrd;
296  bool doneHartVks;
297  double FFTscale;
298  double volumeFactor;
299  double probScale;
300  RhoRealSlab rho_rs;
301  //Comlib multicast proxies
302  CProxySection_CP_State_RealSpacePlane *realSpaceSectionProxyA;
303  CProxySection_CP_State_RealSpacePlane *realSpaceSectionCProxyA;
304  CProxy_CP_Rho_GSpacePlane rhoGProxy_com;
305  CProxy_CP_Rho_GSpacePlane rhoGProxyIGX_com;
306  CProxy_CP_Rho_GSpacePlane rhoGProxyIGY_com;
307  CProxy_CP_Rho_GSpacePlane rhoGProxyIGZ_com;
308  int redCount;
309  CkReductionMsg *RedMsg;
310 };
311 //////////////////////////////////////////////////////////////////////////////
312 /*@}*/
313 
314 /** @addtogroup Density
315  @{
316 */
317 
318 //////////////////////////////////////////////////////////////////////////////
319 //////////////////////////////////////////////////////////////////////////////
320 //////////////////////////////////////////////////////////////////////////////
322  public:
323  CP_Rho_GSpacePlane(CkMigrateMessage *m) {}
324  CP_Rho_GSpacePlane(int, int, int, bool, UberCollection);
326  void run();
327  void init();
328  void acceptRhoData(RhoGSFFTMsg *msg);
329  void doRhoFFT(); // refine the name
330  void ResumeFromSync();
331  void divRhoVksGspace();
332  void RhoGSendRhoR(int );
334  void acceptWhiteByrd();
335  void pup(PUP::er &p);
336  void isAtSync(int iter){AtSync();};
337  int cp_grad_corr_on;
338  int ees_eext_on;
339  int ngridcEext;
340  int rhoRsubplanes;
341  int countDebug;
342  void exitForDebugging();
343  void acceptWhiteByrdAll(RhoGSFFTMsg *msg);
344  void RhoGSendRhoRall();
345  void launchNlG();
346  private:
347  const UberCollection thisInstance;
348  CProxySection_GSpaceDriver nlsectproxy;
349  int myTime;
350  int recvCountFromRRho;
351  int nPacked;
352  int count_stuff;
353  int count;
354  int countWhiteByrd[4];
355  int doneWhiteByrd;
356  int rhoGHelpers;
357  int iplane_ind;
358  int *numSplit;
359  int *istrtSplit;
360  int *iendSplit;
361  RhoGSlab rho_gs;
362  CProxy_CP_Rho_RealSpacePlane rhoRealProxy0_com;
363  CProxy_CP_Rho_RealSpacePlane rhoRealProxy1_com;
364  CProxy_CP_Rho_RealSpacePlane rhoRealProxy2_com;
365  CProxy_CP_Rho_RealSpacePlane rhoRealProxy3_com;
366  CProxy_CP_Rho_RealSpacePlane rhoRealProxyByrd_com;
367 };
368 //////////////////////////////////////////////////////////////////////////////
369 /*@}*/
370 
371 /** @addtogroup Density
372  @{
373 */
374 //////////////////////////////////////////////////////////////////////////////
375 //////////////////////////////////////////////////////////////////////////////
376 //////////////////////////////////////////////////////////////////////////////
378  public:
379  const UberCollection thisInstance;
380  int listSubFlag;
381  int nplane_rho_x;
382  int rhoRsubplanes;
383  int registrationFlag;
384  int launchFlag;
385  int ngrida;
386  int ngridb;
387  int ngridc;
388  int iplane_ind;
389  int ees_eext_on;
390  int natmTyp;
391  int countFFT[2];
392  int countIntRtoG;
393  int countIntGtoR[2];
394  int iteration;
395  int iterAtmTyp;
396  int csize;
397  int nAtmTypRecv;
398  int csizeInt;
399  int myNgridb;
400  int myBoff;
401  int nptsB;
402  int nptsExpndB;
403  int myNplane_rho;
404  int myAoff;
405  int nptsA;
406  int nptsExpndA;
407  int countDebug;
408  int recvCountFromGHartExt;
409  int nchareHartAtmT;
410  int natmTypTot;
411  int atmTypoff;
412 
413  complex *atmSFC;
414  double *atmSFR;
415  complex *atmForcC;
416  double *atmForcR;
417 
418  complex *atmEwdSFC;
419  double *atmEwdSFR;
420  complex *atmEwdForcC;
421  double *atmEwdForcR;
422 
423  complex *atmSFCint;
424  double *atmSFRint;
425  complex *atmForcCint;
426  double *atmForcRint;
427 
428  complex *atmEwdSFCint;
429  double *atmEwdSFRint;
430  complex *atmEwdForcCint;
431  double *atmEwdForcRint;
432 
433  CProxy_CP_Rho_GHartExt rhoGHartProxy_com;
434  CP_Rho_RHartExt(CkMigrateMessage *m) {}
435  CP_Rho_RHartExt(int , int , int , int , int, UberCollection );
436  ~CP_Rho_RHartExt();
437  void init();
438  void pup(PUP::er &p);
439  void startEextIter();
440  void computeAtmSF();
441  void registrationDone(CkReductionMsg *msg);
442  void fftAtmSfRtoG();
443  void sendAtmSfRyToGy();
444  void recvAtmSfRyToGy(RhoGHartMsg *msg);
445  void sendAtmSfRhoGHart();
447  void fftAtmForcGtoR(int flagEwd);
448  void sendAtmForcGxToRx(int iopt);
449  void recvAtmForcGxToRx(RhoGHartMsg *msg);
450  void computeAtmForc(int);
451  void exitForDebugging();
452 };
453 //////////////////////////////////////////////////////////////////////////////
454 /*@}*/
455 
456 /** @addtogroup Density
457  @{
458 */
459 //////////////////////////////////////////////////////////////////////////////
460 //////////////////////////////////////////////////////////////////////////////
461 //////////////////////////////////////////////////////////////////////////////
463  public:
464  CP_Rho_GHartExt(CkMigrateMessage *m) {}
465  CP_Rho_GHartExt(int , int , int , int ,int, UberCollection );
466  void init();
468  void pup(PUP::er &);
469  void acceptData(RhoGHartMsg *msg);
470  void HartExtVksG();
471  void FFTVks();
472  void sendVks();
473  void acceptVks(int size, complex * inVks);
474  void acceptAtmSFTot(int size, complex * inAtm);
476  void FFTEesBck();
477  void getHartEextEes();
478  void FFTEesFwd(int );
479  void sendAtmSF(int );
480  void isAtSync(int iter){AtSync();};
481  int rhoRsubplanes;
482  int ngridaEext;
483  int ngridbEext;
484  int ngridcEext;
485  int ees_eext_on;
486  int natmTyp;
487  int iterAtmTyp;
488  int nsendAtmTyp;
489  int numFullEext;
490  int registrationFlag;
491  int launchFlag;
492  int CountDebug;
493  int iperd;
494  complex *atmSF;
495  complex *atmSFtot;
496  double ehart_ret;
497  double eext_ret;
498  double ewd_ret;
499  void registrationDone(CkReductionMsg *msg);
500  void exitForDebugging();
501  private:
502  const UberCollection thisInstance;
503  complex *atmSFtotRecv;
504  complex *VksRecv;
505  int countAtmSFtot;
506  int countVksTot;
507  int nchareHartAtmT;
508  int natmTypTot;
509  int atmTypoff;
510  int recvCountFromRHartExt;
511  RhoGSlab rho_gs;
512  int atmSFHere;
513  int densityHere;
514  int countEextFFT;
515  int iopt;
516  int iteration;
517  int ind_x; // This chares index=thisIndex.x.
518  int ind_xdiv; // This chare is a subcollection of rhog(ind_xdiv).
519  int ind_xrem; // The subcollection index 0<= ind_rem < rhoGHelpers.
520  int rhoGHelpers; // The number of subcolletions of each rhog(ind_xdiv).
521  int istrt_lines; // start of my subdivion of lines in rhog()
522  int iend_lines; // end of my subdivion of lines in rhog()
523  int numLines; // Number of lines in my subdivision
524  CProxy_CP_Rho_RealSpacePlane rhoRealProxy_com;
525  CProxy_CP_Rho_RHartExt rhoRHartProxy_com0;
526  CProxy_CP_Rho_RHartExt rhoRHartProxy_com1;
527  int **index_pack_tran;
528 };
529 //////////////////////////////////////////////////////////////////////////////
530 /*@}*/
531 
532 /** @addtogroup Particle
533  @{
534 */
535 
536 //////////////////////////////////////////////////////////////////////////////
537 //////////////////////////////////////////////////////////////////////////////
538 //////////////////////////////////////////////////////////////////////////////
540  public:
541  // Variables
542  const UberCollection thisInstance;
543  int ees_nonlocal;
544  int nChareR; // Real Space chares=# C-planes
545  int nChareG; // G Space chares
546  int Rstates_per_pe; // Real Space topomap variable
547  int myPlane; // Real space plane number
548  int ibead_ind,kpoint_ind, itemper_ind;
549 
550  int rhoRTime;
551  int numIterNl; // # of non-local iterations per time step
552  int countZ;
553  int countEnl;
554  int count; // fft communication counter
555  int iterNL; // Nl iteration counter
556  int itime; // time step counter;
557  int recvBlock;
558 
559  int ngridA; // FFT grid size along a
560  int ngridB; // FFT grid size along b
561  int ngridC; // FFT grid size along c
562  int planeSize; // expanded plane size for FFTing
563  int planeSizeT; // true plane size
564  int csize; // complex variable size for FFT
565  int zmatSizeMax; // zmatrix size for projector
566  int reductionPlaneNum; // Reduction Plane number
567  int itimeRed;
568 
569  int registrationFlag;
570 
571  bool launchFFT;
572  bool fftDataDone;
573  bool planeRedSectionComplete;
574  bool enlSectionComplete;
575 
576  double cp_enl; // Non-local energy
577  double cp_enlTot; // Reduced Non-local energy
578  double *projPsiR; // real/final form of projector (after gx,gy FFTs)
579  complex *projPsiC; // complex/intermediate form of projector (before gx,gy FFTs)
580  double *zmat; // Non-local matrix for doublePack
581  double *zmatScr; // Non-local matrix for doublePack
582  complex *zmatC; // Non-local matrix complex for not doublePack
583  complex *zmatScrC; // Non-local matrix complex for not doublePack
584 
585  //-----------
586  // Proxies
587 
588  CProxySection_CP_State_RealParticlePlane rPlaneSectProxy; // Section Red proxy zmat
589  CProxySection_CP_State_RealParticlePlane rPlaneENLProxy; // Section Red proxy cp_enl
590  CkSectionInfo rPlaneRedCookie; // Section Red cookie for zmat
591  CkSectionInfo rEnlCookie; // Section Red cookie for cp_enl
592  CProxy_CP_State_ParticlePlane gPP_proxy;
593 #ifdef _CP_GS_DEBUG_COMPARE_VKS_
594  complex *savedprojpsiC;
595  complex *savedProjpsiCScr;
596  double *savedProjpsiRScr;
597  double *savedzmat;
598  double **savedmn;
599  double **saveddmn_x;
600  double **saveddmn_y;
601  double **saveddmn_z;
602  int **savedigrid;
603 
604 #endif
605  //-----------
606  // Functions
607  CP_State_RealParticlePlane(CkMigrateMessage *m) {}
608  CP_State_RealParticlePlane(int , int , int ,int , int ,int ,int,int, UberCollection);
609  void init();
611  void planeRedSectDone(CkReductionMsg *m);
612  void enlSectDone(CkReductionMsg *m);
613  void initComplete();
614  void launchFFTControl(int );
615  void pup(PUP::er &);
616  void printEnlR(CkReductionMsg *m);
617  void printEnlRSimp(double,int,int);
618  void recvFromEesGPP(NLFFTMsg *);
619  void FFTNLEesFwdR();
620  void computeZmatEes();
621  void recvZMatEes(CkReductionMsg *);
623  void FFTNLEesBckR();
624  void sendToEesGPP();
626  void setEnlCookie(EnlCookieMsg *);
627  int calcReductionPlaneNum(int );
628  void registrationDone();
629  void recvZMatEesSimp(int , double *,int,int,int);
630 };
631 //////////////////////////////////////////////////////////////////////////////
632 /*@}*/
633 
634 /** @addtogroup LargeSparse
635  @{
636 */
637 
638 //////////////////////////////////////////////////////////////////////////////
639 //////////////////////////////////////////////////////////////////////////////
640 //////////////////////////////////////////////////////////////////////////////
641 /**
642  * The Large Sparse RhoG is where we interface with NAMD in QMMM for
643  * the large grid
644 */
646  public:
647  const UberCollection thisInstance;
648  // Functions
649  CP_LargeSP_RhoGSpacePlane(CkMigrateMessage *m) {}
651  void init();
653  void pup(PUP::er &);
654  void acceptMDSg();
655  void acceptLSPRhoR();
656 
657 };
658 //////////////////////////////////////////////////////////////////////////////
659 /*@}*/
660 
661 /** @addtogroup LargeSparse
662  @{
663 */
664 
665 //////////////////////////////////////////////////////////////////////////////
666 //////////////////////////////////////////////////////////////////////////////
667 //////////////////////////////////////////////////////////////////////////////
668 /**
669  * The Large Sparse RhoR is where we interpolate dense RhoR onto the large
670  * grid for QMMM
671 */
673  public:
674  const UberCollection thisInstance;
675  // constructor Functions
676  CP_LargeSP_RhoRealSpacePlane(CkMigrateMessage *m) {}
678  // entry methods
679  void init();
680  void acceptLSPRhoG();
681  void acceptRhoR();
682 
683  // sequential stuff
685  void pup(PUP::er &);
686 
687 };
688 //////////////////////////////////////////////////////////////////////////////
689 
690 //////////////////////////////////////////////////////////////////////////////
691 
692 /*@}*/
693 
694 
695 //////////////////////////////////////////////////////////////////////////////
696 #endif // #ifndef _PLANE_H_
697 //////////////////////////////////////////////////////////////////////////////
698 
699 
void acceptAtmSFTot(int size, complex *inAtm)
Collect the SF from all the atm type chares on chare 0 //////////////////////////////////////////////...
void recvZMatEes(CkReductionMsg *)
Reduction client of zmat : Here everyone must have the same iteration ///////////////////////////////...
void acceptRhoGradVksGxToRx(RhoGSFFTMsg *msg)
Double Transpose Bck Recv : A(gx,y,z) on the way to A(x,y,z) Recv (gx,z) parallel -> (y...
void getHartEextEes()
compute HartreeEextEes
void acceptDensity(CkReductionMsg *)
Data comes from StateRspacePlane once an algorithm step.
void sendPartlyFFTtoRhoG(int)
The Tranpose to G-space : A(gx,gy,z) on the way to A(gx,gy,gz) Change parallel by gx...
void planeRedSectDone(CkReductionMsg *m)
All cookies initialized for zmat plane reduction ////////////////////////////////////////////////////...
holds the UberIndex and the offset for proxies
Definition: Uber.h:68
void startEextIter()
Invoke by Rspace-density : Density has arrived in r-space and will soon arrive in g-space...
void printEnlR(CkReductionMsg *m)
Energy reduction client for ees method! /////////////////////////////////////////////////////////////...
int calcReductionPlaneNum(int)
spread the reduction plane numbers around to minimize map collisions
void sendAtmSfRyToGy()
Double Transpose Fwd Send : A(gx,y,z) on the way to A(gx,gy,z) Send so that (y,z) parallelism is swit...
void addHartEextVks()
Add the VksHartEext to VksTot : Set the done flag.
void init()
post constructor initialization
void computeAtmForcEes(CompAtmForcMsg *msg)
Use Zmat and ProjPsi to get atmForces, Energy and psiforces /////////////////////////////////////////...
void exitForDebugging()
Glenn's Rhart exit.
void divRhoVksGspace()
invoke fft subroutine : doFwFFT() unpack do an FFT 1D along z rho(gx,gy,gz) parallelized in full z-li...
void init()
Post construction initialization.
void addWhiteByrdVks()
Add the VksWhiteByrd to VksTot : Set the done flag.
void registrationDone(CkReductionMsg *msg)
= Make sure everyone is registered on the 1st time step
void acceptLSPRhoG()
Data comes from LSP_RhoG once an algorithm step.
void FFTNLEesBckR()
Do the FFT of projPsi(x,y,z) to get projPsi(gx,gy,z) ////////////////////////////////////////////////...
void exitForDebugging()
Glenn's RhoReal exit.
void sendToEesGPP()
Send the PsiForce NL FFT back to GparticlePlane.
void initComplete()
Initialization and registration done for this chare /////////////////////////////////////////////////...
The Large Sparse RhoR is where we interpolate dense RhoR onto the large grid for QMMM.
void computeAtmForc(int)
Get forces on atoms from Ewald or electrons /////////////////////////////////////////////////////////...
void whiteByrdFFT()
The white-bird term : First fwfft redefined delrho(r) to delrho(g) then send to RhoGspacePlane.
void printEnlRSimp(double, int, int)
Energy reduction client for ees method! /////////////////////////////////////////////////////////////...
~CP_State_RealParticlePlane()
The destructor : never called directly but I guess migration needs it ///////////////////////////////...
void doMulticastCheck()
If all the parts of exc-eext-hart are done, invoke blast of vks to states.
void acceptWhiteByrd(RhoRSFFTMsg *msg)
The white bird vks correction is returned from RhoG : VksW(gx,gy,z) This routine recvs the transpose ...
void recvAtmSFFromRhoRHart(RhoGHartMsg *msg)
Recv Atm SF from RhoRhart : Euler Exponential spline based method ///////////////////////////////////...
void FFTEesFwd(int)
Statr FFting to R-space atmSF(gx,gy,gz) -> atmSF(gx,gy,z) ///////////////////////////////////////////...
void init(ProductMsg *)
Setting up the multicast trees for Gengbin's library.
void setPlaneRedCookie(EnlCookieMsg *)
Recv a dummy message, 1 integer, and set the cookie monster /////////////////////////////////////////...
void launchEextRNlG()
The density is here : Launch ees NL and ees Ext routines.
void fftRhoRyToGy(int iopt)
Double Transpose Fwd FFT : A(gx,y,z) -> A(gx,gy,z)
void GradCorr()
The gradient of the density is now completed.
void FFTNLEesFwdR()
Complete the Forward FFT of psi-projector sent from g-space chares //////////////////////////////////...
void fftAtmSfRtoG()
FFT of SF(x,y,z,iatmTyp) -> SF(gx,gy,z,iatmTyp)
void acceptVks(int size, complex *inVks)
Collect the VKS from all the atm type chares on chare 1 /////////////////////////////////////////////...
void recvFromEesGPP(NLFFTMsg *)
Recv FFT data/psi-projector from CP_StateParticlePlane kicking things off for this Chare...
The Large Sparse RhoG is where we interface with NAMD in QMMM for the large grid. ...
void acceptWhiteByrdAll(RhoGSFFTMsg *msg)
RhoReal sends rho(gx,gy,z) here such that it is now decomposed with lines of constant gx...
void doMulticast()
Send vks back to the states /////////////////////////////////////////////////////////////////////////...
void enlSectDone(CkReductionMsg *m)
All cookies initialized for enl section, only reduction root receives this //////////////////////////...
void pup(PUP::er &)
Pup my variables for migration.
void acceptHartVks(RhoHartRSFFTMsg *)
Accept hartExt transpose data : receive VksHartEext(gx,gy,z) gx,z is parallel.
void doReduction()
No one else can use tmpdataR until I perform doReduction because I will not be rolled out until I fin...
void acceptGradRhoVks(RhoRSFFTMsg *)
Accept transpose data from RhoG : receive grad_rho(gy,gx,z)
void launchNLRealFFT()
Launch ees-nonlocal real here.
void sendAtmSF(int)
Send the SF data to back to Rhart to get atm forces /////////////////////////////////////////////////...
void sendAtmSfRhoGHart()
Send SF(gx,gy,z,iatmTYP) to g-space whence the FFT will be completed.
void FFTVks()
Partly fft vks(gx,gy,gz) -> vks(gx,gy,z)
void acceptData(RhoGHartMsg *msg)
The density arrives from RhoGspace ONCE a time step /////////////////////////////////////////////////...
void computeAtmSF()
Start the real space part of the EES interpolation for atmSF(iatmTyp) ///////////////////////////////...
void acceptLSPRhoR()
Data comes from LargeSP_RhoReal once an algorithm step.
void handleDensityReduction()
Handle the memory cleanup and setting of flags when density has all arrived.
void init()
post constructor initialization
void sendPartlyFFTGxToRx(int)
Double Transpose Bck Send : A(gx,y,z) on the way to A(x,y,z) Send so (gx,z) parallel -> (y...
void recvZMatEesSimp(int, double *, int, int, int)
A chare can be behind by 1 iteration only.
void energyComputation()
Compute one part of the EXC energy using PINY CP_exc_calc.
~CP_Rho_GHartExt()
Destructor /////////////////////////////////////////////////////////////////////////// //////////////...
void init()
Post constructor initialization /////////////////////////////////////////////////////////////////////...
void acceptRhoGradVksRyToGy(RhoGSFFTMsg *msg)
Double Transpose Fwd Recv : A(gx,y,z) on the way to A(gx,gy,z) Recv so that (y,z) parallel switched t...
void doFFT()
After receiving from G-space, FFT to real space /////////////////////////////////////////////////////...
void getSplitDecomp(int *, int *, int *, int, int, int)
Initialization Function /////////////////////////////////////////////////////////////////////////// /...
Definition: util.C:2250
void fftAtmForcGtoR(int flagEwd)
Complete the FFT : sf(gx,gy,z) ->sf(x,y,z) //////////////////////////////////////////////////////////...
void launchNlG()
The density is here : Launch ees NL.
void launchFFTControl(int)
Control launch of FFT based on Rho having more of its act together //////////////////////////////////...
void registrationDone()
= Make sure everyone is registered on the 1st time step
void acceptRhoData(RhoGSFFTMsg *msg)
RhoReal sends rho(gx,gy,z) here such that it is now decomposed with lines of constant gx...
void fftRhoRtoRhoG()
1) Perform FFT of density: Single Transpose method rho(x,y,z) —> rho(gx,gy,z) Double Transpose method...
void sendVks()
Send vks_hart_ext back to rho_real where fft(gx,gy) will be performed ///////////////////////////////...
void recvAtmForcGxToRx(RhoGHartMsg *msg)
void acceptGradRhoVksAll(RhoRSFFTMsg *msg)
Accept transpose data from RhoG : receive grad_rho(gy,gx,z)
void sendPartlyFFTRyToGy(int iopt)
Double Transpose Fwd Send : A(gx,y,z) on the way to A(gx,gy,z) Send so that (y,z) parallelism is swit...
void FFTEesBck()
Finish FFting to G-space :: 2D) atmSF(gx,gy,z) -> atmSF(gx,gy,gz)
void acceptRhoR()
Data comes from Rho_Real once an algorithm step.
Useful debugging flags.
void registrationDone(CkReductionMsg *msg)
= Make sure everyone is registered on the 1st time step
void acceptMDSg()
Data comes from MD once an algorithm step.
void pup(PUP::er &)
Pup my variables for migration.
void HartExtVksG()
Compute hartree eext and vks using the N^2 method ///////////////////////////////////////////////////...
void doRhoFFT()
Complete the FFT to give use Rho(gx,gy,gz) //////////////////////////////////////////////////////////...
void setEnlCookie(EnlCookieMsg *)
Recv a dummy message, 1 integer, and set the cookie monster /////////////////////////////////////////...
void computeZmatEes()
Compute the Zmat elements of this iteration : Spawn the section reduction ///////////////////////////...
void recvAtmSfRyToGy(RhoGHartMsg *msg)
Double Transpose Fwd Recv : A(gx,y,z) on the way to A(gx,gy,z) Recv so that (y,z) parallel switched t...
void unpackProduct(ProductMsg *)
In this method, we receive vks from the density.
void recvAtmForcFromRhoGHart(RhoRHartMsg *msg)
Hartree sends back atom forces from e-atm interation Depending on the flag, it is Ewald or e-atm inte...
void exitForDebugging()
Glenn's special exit.
void init()
post constructor initialization /////////////////////////////////////////////////////////////////////...
void RHartReport()
Under ees-eext Rhart chare reports its completion : Set the done flag.
CP_State_RealSpacePlane_SDAG_CODE CP_State_RealSpacePlane(int, int, int, int, int, int, int, UberCollection)
void exitForDebugging()
Glenn's RhoReal exit.
void sendAtmForcGxToRx(int iopt)
void init()
post constructor initialization
void pup(PUP::er &)
Pup my variables for migration.