-
Notifications
You must be signed in to change notification settings - Fork 1
/
Copy pathLATfield2_Field.hpp
executable file
·1634 lines (1271 loc) · 48.2 KB
/
LATfield2_Field.hpp
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
747
748
749
750
751
752
753
754
755
756
757
758
759
760
761
762
763
764
765
766
767
768
769
770
771
772
773
774
775
776
777
778
779
780
781
782
783
784
785
786
787
788
789
790
791
792
793
794
795
796
797
798
799
800
801
802
803
804
805
806
807
808
809
810
811
812
813
814
815
816
817
818
819
820
821
822
823
824
825
826
827
828
829
830
831
832
833
834
835
836
837
838
839
840
841
842
843
844
845
846
847
848
849
850
851
852
853
854
855
856
857
858
859
860
861
862
863
864
865
866
867
868
869
870
871
872
873
874
875
876
877
878
879
880
881
882
883
884
885
886
887
888
889
890
891
892
893
894
895
896
897
898
899
900
901
902
903
904
905
906
907
908
909
910
911
912
913
914
915
916
917
918
919
920
921
922
923
924
925
926
927
928
929
930
931
932
933
934
935
936
937
938
939
940
941
942
943
944
945
946
947
948
949
950
951
952
953
954
955
956
957
958
959
960
961
962
963
964
965
966
967
968
969
970
971
972
973
974
975
976
977
978
979
980
981
982
983
984
985
986
987
988
989
990
991
992
993
994
995
996
997
998
999
1000
#ifndef LATFIELD2_FIELD_HPP
#define LATFIELD2_FIELD_HPP
/*! \file LATfield2_Field.hpp
\brief LATfield2_Field.hpp contain the class Field definition.
\author David Daverio, Neil Bevis
*/
//Matrix-style component symmetry
int symmetric = 1;
int unsymmetric = 0;
//PROTO-TYPEs============================================
template<class FieldType>
void defaultFieldSave(fstream& file, FieldType* siteData, int components);
template<class FieldType>
void defaultFieldLoad(fstream& file, FieldType* siteData, int components);
/*! \class Field
\brief The Field class represent a field on a given lattice.
It stores the description of the field i.e. the datatype, the number of components, and the pointer to the field array in the memory.
The datatype is under user control; a field of structur or class can be also used (to be able to use the update halo method, the = operator must be defined or overloaded). However, the I/O support only native datatype and 1d array of them.
A field can be a single element, a vector of element or a 2d matrix of elements. In the case of a matrix it is possible to define the symmetry of the matrix:
LATfield2d::unsymmetric : no symmetry.
LATfield2d::symmetric : symmetric matrix (Tij = Tji)
Antisymmetric field are not yet implemented.
*/
//FIELD CLASS DEFINITION==================================
template<class FieldType>
class Field
{
public:
//! Constructor.
Field();
/*!
Constructor of a "vector" field with initialization and allocation.
\sa initialize(Lattice& lattice, int components = 1)
\sa alloc()
\param lattice : lattice on which the field is defined
\param components : number of components. The default is 1.
*/
Field(Lattice& lattice, int components = 1);
/*!
Constructor of a "matrix" field with initialization and allocation.
\sa initialize(Lattice& lattice, int rows, int cols, int symmetry=unsymmetric)
\sa alloc()
\param lattice : lattice on which the field is defined
\param matrixRows : matrix number of row .
\param matrixCols : matrix number of colomn.
\param symmetry : symmetry of the matrix, default is unsymmetric. LATfield2d::symmetric can be passed to specify the symmetry, reducing memory usage.
*/
Field(Lattice& lattice, int matrixRows, int matrixCols, int symmetry=unsymmetric);
//!Destructor.
~Field();
//INITIALIZATION-TYPE FUNCTIONS
/*!
Initialization of a "vector" field. Without allocation.
\param lattice : lattice on which the field is defined
\param components : number of components. Default is 1.
*/
void initialize(Lattice& lattice, int components = 1);
/*!
Initialization of a "matrix" field. Without allocation.
\param lattice : lattice on which the field is defined
\param matrixRows : matrix number of row .
\param matrixCols : matrix number of colomn.
\param symmetry : symmetry of the matrix, default is unsymmetric. LATfield2d::symmetric can be pass to specify the symmetry.
*/
void initialize(Lattice& lattice, int rows, int cols, int symmetry=unsymmetric);
/*!
Memory allocation. Allocate the data_ array of this field. It allocated "components_*lattice_->sitesLocalGross()*sizeof(FieldType)" bytes. This method use malloc() to allocate the memory, in case the pointer is not allocated it will return a error message but not exiting the executable.
*/
void alloc();
/*!
Memory allocation. Allocate the data_ array of this field. It allocated "size" bytes if "size" > "components_*lattice_->sitesLocalGross()*sizeof(FieldType)", if not it call this->alloc(). This method use malloc() to allocate the memory, in case the pointer is not allocated it will return a error message but not exiting the executable.
*/
void alloc(long size);
/*!
Free the data_ array.
*/
void dealloc();
//INDEXING
/*!
Returns the value of the field stored in data_[index]. User should used operator()(const Site& site) to refer and access to the value of the field.
\param index: displacment on the data_ array.
*/
FieldType& operator()(long index);
/*!
Returns the value of the field stored in data_[component + index*components_]. User should used operator()(const Site& site, int component) to refer and access to the value of the field.
\param index : number of site to skip.
\param component: index of the desired component.
*/
FieldType& operator()(long index, int component);
/*!
Returns the value of the field stored in data_[j*rows_ + i + index*components_]. In the symmetric case, it returns data_[abs(i-j) + min(i,j)*(rows_+0.5-0.5*min(i,j)) + index*components_]. User should used operator()(const Site& site, int i, int j) to refer and access to the value of the field.
\param index : number of site to skip.
\param i : index of the row
\param j : index of the column
*/
FieldType& operator()(long index, int i, int j);
/*!
Returns the value of the field at the position pointed by the Site object (data_[site.index()]). Can be used only for field with one component!
\param site: a site instance which points to the desired lattice site.
\sa To have more description see the Site class documentation.
*/
FieldType& operator()(const Site& site);
/*!
Returns the value of a (vector) field's component at the position pointed by the Site object (data_[component + site.index()*components_]).
\param site: a site instance which points to the desired lattice site.
\param component: index of the desired component.
\sa To have more description see the Site class documentation.
*/
FieldType& operator()(const Site& site, int component);
/*!
Returns the value of the (i,j) matrix component of the field at the position pointed by the Site object (data_[j*rows_ + i + site.index*components_]). In the symmetric case, it returns data_[abs(i-j) + min(i,j)*(rows_+0.5-0.5*min(i,j)) + site.index()*components_].
\param site: a site instance which points to the desired lattice site.
\param i : index of the row
\param j : index of the column
\sa To have more description see the Site class documentation.
*/
FieldType& operator()(const Site& site, int i, int j);
#ifdef FFT3D
/*!
Equivalent to FieldType& operator()(const Site& site) for cKsite
*/
FieldType& operator()(const cKSite& site);
/*!
Equivalent to FieldType& operator()(const Site& site, int component) for cKsite
*/
FieldType& operator()(const cKSite& site, int component);
/*!
Equivalent to FieldType& operator()(const Site& site, int i, int j) for cKsite
*/
FieldType& operator()(const cKSite& site, int i, int j);
/*!
Equivalent to FieldType& operator()(const Site& site) for rKsite
*/
FieldType& operator()(const rKSite& site);
/*!
Equivalent to FieldType& operator()(const Site& site, int component) for rKsite
*/
FieldType& operator()(const rKSite& site, int component);
/*!
Equivalent to FieldType& operator()(const Site& site, int i, int j) for rKsite
*/
FieldType& operator()(const rKSite& site, int i, int j);
#endif
//BOUNDARY UPDATE
/*!
Update the halo sites (ghost cells) of the field. This method us the operator = to assign values, thefor be sure that this operator is defined or correctly overloaded when using field of class or struct.
*/
void updateHalo();
//FILE I/O FUNCTIONS
/*!
Method to write a field in binary. This method use serial I/O so can be very slow. Should never be used during production, but can be usefull during development.
\param filename: path to the file, from the executable folder.
*/
void write(const string filename);
/*!
Method to read a field in binary which have been writen by the void write(const string filename) method.
\param filename: path to the file, from the executable folder.
*/
void read(const string filename);
/*!
Method to write a field in Binary. This method uses serial I/O so can be very slow, but is faster than void save(const string filename). Should never be used! but it can be usefull on some architectures, where HDF5 is not installed and/or crashes the filesystem. There is no method to read back such a file. The file structure is dependent of the local geometry. This function dumps serially (in the paralle.lat_world_rank order) the data stored in each MPI process.
\param filename: path to the file, from the executable folder.
*/
void fastwrite(const string filename);
/*!
Method to write a field in ASCII. This method use serial I/O so can be very slow. Should never be used during production, but can be useful during development.
\param filename : path to the file, from the executable folder.
\param FormatFunction: format used for the writting procedure.
*/
void save(const string filename,
void (*FormatFunction)(fstream&,FieldType*,int) = defaultFieldSave<FieldType>);
/*!
Method to read a field in ASCII which have been written by the void write(const string filename) method.
\param filename : path to the file, from the executable folder.
\param FormatFunction: format used for the writting procedure.
*/
void load(const string filename,
void (*FormatFunction)(fstream&,FieldType*,int) = defaultFieldLoad<FieldType>);
/*!
Method to write a field in ASCII. This method use serial I/O so can be very slow, but is faster than void write(const string filename). Should never be used! but it can be usefull on some architectures, where HDF5 is not installed and/or crashes the filesystem. There is no method to read back such a file. The file structur is dependent of the local geometry. This function dumps serially (in the paralle.lat_world_rank order) the data stored in each MPI process.
*/
void fastsave(const string filename,
void (*FormatFunction)(fstream&,FieldType*,int) = defaultFieldSave<FieldType>);
/*!
Method to write a field with HDF5. To be able to use this method the flag HDF5 need to be set at compilation (-DHDF5). This method use serial HDF5 by default. For parallel HDF5 the flag -DH5_HAVE_PARALLEL must be used at compilation.
This methods will write 1 dataset (named "/field") which contain all components of the field. If one want to use a dataset per components (named "/comp0" to "/compN") the flag -DH5_HAVE_PIXIE need to be set at compilation.
\param filename : path to the file, from the executable folder.
*/
void saveHDF5(string filename);
/*!
Method to load a field with HDF5. To be able to use this method the flag HDF5 need to be set at compilation (-DHDF5). This method use serial HDF5 by default. For parallel HDF5 the flag -DH5_HAVE_PARALLEL must be set at compilation.
This methods will expect 1 dataset named field which contain all component of the field. If one want to use a dataset per components (named comp0 to compN) the flag -DH5_HAVE_PIXIE need to be set at compilation.
\param filename : path to the file, from the executable folder.
*/
void loadHDF5(string filename);
/*!
A way to save coarse grained version of the fields. To be able to use this method the flag HDF5 need to be set at compilation (-DHDF5). Work only for 3D lattice!!!
This methods will write 1 dataset (named "/field") which contain all component of the field. If one want to use a dataset per components (named "/comp0" to "/compN") the flag -DH5_HAVE_PIXIE need to be set at compilation.
\param filename : path to the file, from the executable folder.
\param ration : ration of the coarse graining. Must be an integer divider of the size of each dimension of this->lattice()
*/
void saveHDF5_coarseGrain3D(string filename,int ratio);
/*!
Save a slice perpendicular to the first coordinate, at xcoord. To be able to use this method the flag HDF5 need to be set at compilation (-DHDF5).
This methods will write 1 dataset (named "/field") which contain all component of the field. If one want to use a dataset per components (named "/comp0" to "/compN") the flag -DH5_HAVE_PIXIE need to be set at compilation.
\param filename : path to the file, from the executable folder.
\param xcoord : coordinate of the slice on the first dimension of the lattice.
\param thickness: thickness of the slice, the default is 1.
*/
void saveSliceHDF5(string filename, int xcoord, int thickness = 1);
//MISCELLANEOUS
/*!
Returns a pointer to the lattice on which the field is defined.
*/
Lattice& lattice();
/*!
Returns the number of components of the field at each sites.
*/
int components();
/*!
Returns the number of rows of the component matrix at each sites.
*/
int rows();
/*!
Returns the number of columns of the component matrix at each sites.
*/
int cols();
/*!
returns the symmetry of the component matrix at each sites.
*/
int symmetry();
/*!
Returns the pointer to the data_ array of the field.
*/
FieldType*& data();
private:
//PRIVATE FUNCTIONS
void updateHaloComms();
#ifdef HDF5
void get_h5type();
#endif
public:
FieldType* data_;
protected:
//MEMBER DATA
Lattice* lattice_;
int components_;
int rows_;
int cols_;
int symmetry_;
unsigned int sizeof_fieldType_;
int status_;
static int initialized;
static int allocated;
#ifdef HDF5
hid_t type_id;
int array_size;
#endif
};
template <class FieldType>
int Field<FieldType>::initialized = 1; //Status flag for initialized
template <class FieldType>
int Field<FieldType>::allocated = 2; //Status flag for allocated memory
//CONSTRUCTORS=================
template <class FieldType>
Field<FieldType>::Field() {
status_=0;
#ifdef HDF5
this->get_h5type();
#endif
}
template <class FieldType>
Field<FieldType>::Field(Lattice& lattice, int components)
{
status_=0;
this->initialize(lattice, components);
this->alloc();
#ifdef HDF5
this->get_h5type();
#endif
}
template <class FieldType>
Field<FieldType>::Field(Lattice& lattice, int rows, int cols, int symmetry)
{
status_=0;
this->initialize(lattice, rows, cols, symmetry);
this->alloc();
#ifdef HDF5
this->get_h5type();
#endif
}
#ifdef HDF5
template <class FieldType>
void Field<FieldType>::get_h5type()
{
string type_name;
type_name = typeid(FieldType).name();
char str_array_size[20];
//string str_array_size;
array_size = 1;
int nt;
nt = 0;
//get the array size
if(type_name[0]=='A')
{
while(type_name[nt+1]!='_')
{
str_array_size[nt]=type_name[nt+1];
nt++;
}
str_array_size[nt]='\0';
nt = nt + 2;
if(type_name[nt]=='A')
{
COUT<<"LATField2d::Field::save_hdf5 : Cannot recognize type of field: "<< type_name <<endl;
COUT<<"LATField2d::Field::save_hdf5 : ---------------------------------------------------------------------------------" <<endl;
COUT<<"LATField2d::Field::save_hdf5 : List of accepted type : short, unsigned short, int, unsigned int, long," <<endl;
COUT<<"LATField2d::Field::save_hdf5 : unsigned long, long long, unsigned long long, float, double, long double, Imag ." <<endl;
COUT<<"LATField2d::Field::save_hdf5 : 1 dimensional array of those type are also accepted" <<endl;
COUT<<"LATField2d::Field::save_hdf5 : ---------------------------------------------------------------------------------" <<endl;
return;
}
array_size=atoi(str_array_size);
}
//get the type
if(type_name[nt]=='s')
{
//COUT << " type : short ; size : "<<array_size<<endl;
type_id = H5T_NATIVE_SHORT;
}
else if(type_name[nt]=='t')
{
//COUT << " type : unsigned short ; size : "<<array_size<<endl;
type_id = H5T_NATIVE_USHORT;
}
else if(type_name[nt]=='i')
{
//COUT << " type : int ; size : "<<array_size<<endl;
type_id = H5T_NATIVE_INT;
}
else if(type_name[nt]=='j')
{
//COUT << " type : unsigned int ; size : "<<array_size<<endl;
type_id = H5T_NATIVE_UINT;
}
else if(type_name[nt]=='l')
{
//COUT << " type : long ; size : "<<array_size<<endl;
type_id = H5T_NATIVE_LONG;
}
else if(type_name[nt]=='m')
{
//COUT << " type : unsigned long ; size : "<<array_size<<endl;
type_id = H5T_NATIVE_ULONG;
}
else if(type_name[nt]=='x')
{
//COUT << " type : long long ; size : "<<array_size<<endl;
type_id = H5T_NATIVE_LLONG;
}
else if(type_name[nt]=='y')
{
//COUT << " type : unsigned long long ; size : "<<array_size<<endl;
type_id = H5T_NATIVE_ULLONG;
}
else if(type_name[nt]=='f')
{
//COUT << " type : float ; size : "<<array_size<<endl;
type_id = H5T_NATIVE_FLOAT;
}
else if(type_name[nt]=='d')
{
//COUT << " type : double ; size : "<<array_size<<endl;
type_id = H5T_NATIVE_DOUBLE;
}
else if(type_name[nt]=='e')
{
//COUT << " type : long double ; size : "<<array_size<<endl;
type_id = H5T_NATIVE_LDOUBLE;
}
else if (type_name[nt]=='N'&&type_name[nt+12]=='I'&&type_name[nt+13]=='m'&&type_name[nt+14]=='a'&&type_name[nt+15]=='g')
{
type_id = H5Tcreate (H5T_COMPOUND, sizeof (Imag));
#ifdef SINGLE
H5Tinsert (type_id, "real", 0,H5T_NATIVE_FLOAT);
H5Tinsert (type_id, "imaginary", sizeof(Real),H5T_NATIVE_FLOAT);
#else
H5Tinsert (type_id, "real", 0,H5T_NATIVE_DOUBLE);
H5Tinsert (type_id, "imaginary", sizeof(Real),H5T_NATIVE_DOUBLE);
#endif
}
else
{
COUT<<"LATField2d::Field::save_hdf5 : Cannot recognize type of field: " << type_name <<endl;
COUT<<"LATField2d::Field::save_hdf5 : ---------------------------------------------------------------------------------" <<endl;
COUT<<"LATField2d::Field::save_hdf5 : List of accepted type : short, unsigned short, int, unsigned int, long," <<endl;
COUT<<"LATField2d::Field::save_hdf5 : unsigned long, long long, unsigned long long, float, double, long double, Imag ." <<endl;
COUT<<"LATField2d::Field::save_hdf5 : 1 dimensional array of those type are also accepted" <<endl;
COUT<<"LATField2d::Field::save_hdf5 : ---------------------------------------------------------------------------------" <<endl;
}
}
#endif
//DESTRUCTOR=======================
template <class FieldType>
Field<FieldType>::~Field()
{
if(status_ & allocated) this->dealloc();
}
template <class FieldType>
void Field<FieldType>::initialize(Lattice& lattice, int components)
{
if(status_ == allocated) { this->dealloc(); }
sizeof_fieldType_ = sizeof(FieldType);
status_= initialized;
lattice_=&lattice;
components_=components;
rows_=components_;
cols_=1;
symmetry_=LATfield2::unsymmetric;
}
template <class FieldType>
void Field<FieldType>::initialize(Lattice& lattice, int rows, int cols, int symmetry)
{
int components;
if(status_ == allocated) { this->dealloc(); }
sizeof_fieldType_ = sizeof(FieldType);
status_= initialized;
lattice_=&lattice;
rows_=rows;
cols_=cols;
symmetry_=symmetry;
if(symmetry_==LATfield2::symmetric) { components = ( rows_ * (rows_+1) ) / 2; }
else { components = rows*cols; }
components_=components;
}
template <class FieldType>
void Field<FieldType>::alloc()
{
if((status_ & allocated) == 0)
{
data_= new FieldType[components_*lattice_->sitesLocalGross()];
status_ = status_ | allocated;
if(data_==NULL)
{
cout<<"LATField2d::Field::alloc(long size) :process "<< parallel.rank() <<" cannot allocate the field data array."<<endl;
}
}
}
template <class FieldType>
void Field<FieldType>::alloc(long size)
{
long alloc_number;
if(size < lattice_->sitesLocalGross()) alloc_number = lattice_->sitesLocalGross();
else alloc_number= size;
if((status_ & allocated) == 0)
{
data_= new FieldType[components_* size];
status_ = status_ | allocated;
if(data_==NULL)
{
cout<<"LATField2d::Field::alloc(long size) :process "<< parallel.rank() <<" cannot allocate the field data array."<<endl;
}
}
else if(size > lattice_->sitesLocalGross())
{
this->dealloc();
this->alloc(size);
}
}
template <class FieldType>
void Field<FieldType>::dealloc()
{
if((status_ & allocated) > 0)
{
delete[] data_;
status_= status_ ^ allocated;
}
}
//FIELD INDEXING===============
template <class FieldType>
inline FieldType& Field<FieldType>::operator()(long index)
{
return data_[index];
}
template <class FieldType>
inline FieldType& Field<FieldType>::operator()(long index, int component)
{
return data_[index*components_ + component];
}
template <class FieldType>
inline FieldType& Field<FieldType>::operator()(long index, int i, int j)
{
int component;
if(symmetry_==LATfield2::symmetric)
{
int min=i;
if(i>j) min=j;
component = int( abs(i-j) + min*(rows_+0.5-0.5*min) ); //Will always be an int anyway
}
else { component = j*rows_ + i; }
return data_[index*components_ + component];
}
template <class FieldType>
inline FieldType& Field<FieldType>::operator()(const Site& site)
{
return this->operator()(site.index());
}
template <class FieldType>
inline FieldType& Field<FieldType>::operator()(const Site& site, int component)
{
return this->operator()(site.index(),component);
}
template <class FieldType>
inline FieldType& Field<FieldType>::operator()(const Site& site, int i, int j)
{
return this->operator()(site.index(),i,j);
}
#ifdef FFT3D
template <class FieldType>
inline FieldType& Field<FieldType>::operator()(const cKSite& site)
{
return this->operator()(site.index());
}
template <class FieldType>
inline FieldType& Field<FieldType>::operator()(const cKSite& site, int component)
{
return this->operator()(site.index(),component);
}
template <class FieldType>
inline FieldType& Field<FieldType>::operator()(const cKSite& site, int i, int j)
{
return this->operator()(site.index(),i,j);
}
template <class FieldType>
inline FieldType& Field<FieldType>::operator()(const rKSite& site)
{
return this->operator()(site.index());
}
template <class FieldType>
inline FieldType& Field<FieldType>::operator()(const rKSite& site, int component)
{
return this->operator()(site.index(),component);
}
template <class FieldType>
inline FieldType& Field<FieldType>::operator()(const rKSite& site, int i, int j)
{
return this->operator()(site.index(),i,j);
}
#endif
//FIELD BOUNDARY UPDATE=========
////////////////////////////////
template <class FieldType>
void Field<FieldType>::updateHalo()
{
Site site(*lattice_);
int copyfrom;
int dim=lattice_->dim();
int* jump=new int[dim];
int* size=new int[dim];
for(int i=0; i<dim; i++)
{
jump[i]=lattice_->jump(i);
size[i]=lattice_->sizeLocal(i);
}
for(site.haloFirst(); site.haloTest(); site.haloNext())
{
//Work out where to copy from
copyfrom=site.index();
for(int i=0; i<dim; i++)
{
if( site.coordLocal(i)<0 )
{
copyfrom += jump[i] * size[i];
}
else if( site.coordLocal(i) >= size[i] )
{
copyfrom -= jump[i] * size[i];
}
}
//Copy data
for(int i=0; i<components_; i++)
{
data_[site.index()*components_+i] = data_[copyfrom*components_+i];
//memcpy(data_[site.index()*components_+i],data_[copyfrom*components_+i],sizeof_fieldType_);
}
}
delete[] jump;
delete[] size;
if( parallel.size()>1 ) { updateHaloComms(); }
}
template <class FieldType>
void Field<FieldType>::updateHaloComms()
{
int buffer_size0, buffer_size1,temp;
int i,j;
//Size of buffer : max size between 2 scatered dimension;
buffer_size0 = lattice_->halo() * components_*lattice_->jump(lattice_->dim()-1);
buffer_size1 = lattice_->halo() * components_*lattice_->jump(lattice_->dim()-2) *lattice_->sizeLocal(lattice_->dim()-1);
if(buffer_size0>buffer_size1)temp=buffer_size0;
else temp=buffer_size1;
FieldType* buffer_send = new FieldType[ temp ];
FieldType* buffer_rec = new FieldType[ temp ];
FieldType* pointer_send_up;
FieldType* pointer_send_down;
FieldType* pointer_rec_up;
FieldType* pointer_rec_down;
pointer_send_up = data_ + ((lattice_->halo()+1)*lattice_->jump(lattice_->dim()-1) - 2*lattice_->halo()*lattice_->jump(lattice_->dim()-2))*components_;
pointer_rec_up = data_ + ((lattice_->halo()+1)*lattice_->jump(lattice_->dim()-1) - lattice_->halo()*lattice_->jump(lattice_->dim()-2))*components_;
pointer_send_down = data_ + buffer_size0 + lattice_->jump(lattice_->dim()-2)*lattice_->halo()*components_ ;
pointer_rec_down = data_ + buffer_size0;
if(parallel.grid_rank()[1]%2==0)
{
for(j=0;j<lattice_->sizeLocal(lattice_->dim()-1);j++)
{
for(i=0;i<buffer_size1/lattice_->sizeLocal(lattice_->dim()-1);i++)
{
buffer_send[i+j*buffer_size1/lattice_->sizeLocal(lattice_->dim()-1)]= pointer_send_up[i+j*lattice_->jump(lattice_->dim()-1)*components_];
}
}
if(parallel.grid_rank()[1]!=parallel.grid_size()[1]-1)
{
parallel.send_dim1( buffer_send, buffer_size1, parallel.grid_rank()[1]+1);
parallel.receive_dim1( buffer_rec, buffer_size1, parallel.grid_rank()[1]+1);
}
for(j=0;j<lattice_->sizeLocal(lattice_->dim()-1);j++)
{
for(i=0;i<buffer_size1/lattice_->sizeLocal(lattice_->dim()-1);i++)
{
pointer_rec_up[i+j*lattice_->jump(lattice_->dim()-1)*components_] = buffer_rec[i+j*buffer_size1/lattice_->sizeLocal(lattice_->dim()-1)];
}
}
for(j=0;j<lattice_->sizeLocal(lattice_->dim()-1);j++)
{
for(i=0;i<buffer_size1/lattice_->sizeLocal(lattice_->dim()-1);i++)
{
buffer_send[i+j*buffer_size1/lattice_->sizeLocal(lattice_->dim()-1)]= pointer_send_down[i+j*lattice_->jump(lattice_->dim()-1)*components_];
}
}
if(parallel.grid_rank()[1] != 0)
{
parallel.send_dim1( buffer_send, buffer_size1, parallel.grid_rank()[1]-1);
parallel.receive_dim1( buffer_rec, buffer_size1, parallel.grid_rank()[1]-1);
}
else if(parallel.grid_size()[1]%2==0)
{
parallel.send_dim1( buffer_send, buffer_size1, parallel.grid_size()[1]-1);
parallel.receive_dim1( buffer_rec, buffer_size1, parallel.grid_size()[1]-1);
}
for(j=0;j<lattice_->sizeLocal(lattice_->dim()-1);j++)
{
for(i=0;i<buffer_size1/lattice_->sizeLocal(lattice_->dim()-1);i++)
{
pointer_rec_down[i+j*lattice_->jump(lattice_->dim()-1)*components_] = buffer_rec[i+j*buffer_size1/lattice_->sizeLocal(lattice_->dim()-1)];
}
}
}
else
{
for(j=0;j<lattice_->sizeLocal(lattice_->dim()-1);j++)
{
for(i=0;i<buffer_size1/lattice_->sizeLocal(lattice_->dim()-1);i++)
{
buffer_send[i+j*buffer_size1/lattice_->sizeLocal(lattice_->dim()-1)]= pointer_send_down[i+j*lattice_->jump(lattice_->dim()-1)*components_];
}
}
parallel.receive_dim1( buffer_rec, buffer_size1, parallel.grid_rank()[1]-1);
parallel.send_dim1( buffer_send, buffer_size1, parallel.grid_rank()[1]-1);
for(j=0;j<lattice_->sizeLocal(lattice_->dim()-1);j++)
{
for(i=0;i<buffer_size1/lattice_->sizeLocal(lattice_->dim()-1);i++)
{
pointer_rec_down[i+j*lattice_->jump(lattice_->dim()-1)*components_] = buffer_rec[i+j*buffer_size1/lattice_->sizeLocal(lattice_->dim()-1)];
}
}
for(j=0;j<lattice_->sizeLocal(lattice_->dim()-1);j++)
{
for(i=0;i<buffer_size1/lattice_->sizeLocal(lattice_->dim()-1);i++)
{
buffer_send[i+j*buffer_size1/lattice_->sizeLocal(lattice_->dim()-1)]= pointer_send_up[i+j*lattice_->jump(lattice_->dim()-1)*components_];
}
}
if(parallel.grid_rank()[1]!=parallel.grid_size()[1]-1)
{
parallel.receive_dim1( buffer_rec, buffer_size1, parallel.grid_rank()[1]+1);
parallel.send_dim1( buffer_send, buffer_size1, parallel.grid_rank()[1]+1);
}
else
{
parallel.receive_dim1( buffer_rec, buffer_size1,0);
parallel.send_dim1( buffer_send, buffer_size1, 0);
}
for(j=0;j<lattice_->sizeLocal(lattice_->dim()-1);j++)
{
for(i=0;i<buffer_size1/lattice_->sizeLocal(lattice_->dim()-1);i++)
{
pointer_rec_up[i+j*lattice_->jump(lattice_->dim()-1)*components_] = buffer_rec[i+j*buffer_size1/lattice_->sizeLocal(lattice_->dim()-1)];
}
}
}
if(parallel.grid_size()[1]%2!=0)
{
if(parallel.grid_rank()[1]==0)
{
for(j=0;j<lattice_->sizeLocal(lattice_->dim()-1);j++)
{
for(i=0;i<buffer_size1/lattice_->sizeLocal(lattice_->dim()-1);i++)
{
buffer_send[i+j*buffer_size1/lattice_->sizeLocal(lattice_->dim()-1)]= pointer_send_down[i+j*lattice_->jump(lattice_->dim()-1)*components_];
}
}
parallel.send_dim1( buffer_send, buffer_size1, parallel.grid_size()[1]-1);
parallel.receive_dim1( buffer_rec, buffer_size1, parallel.grid_size()[1]-1);
for(j=0;j<lattice_->sizeLocal(lattice_->dim()-1);j++)
{
for(i=0;i<buffer_size1/lattice_->sizeLocal(lattice_->dim()-1);i++)
{
pointer_rec_down[i+j*lattice_->jump(lattice_->dim()-1)*components_] = buffer_rec[i+j*buffer_size1/lattice_->sizeLocal(lattice_->dim()-1)];
}
}
}
if(parallel.grid_rank()[1]==parallel.grid_size()[1]-1)
{
for(j=0;j<lattice_->sizeLocal(lattice_->dim()-1);j++)
{
for(i=0;i<buffer_size1/lattice_->sizeLocal(lattice_->dim()-1);i++)
{
buffer_send[i+j*buffer_size1/lattice_->sizeLocal(lattice_->dim()-1)]= pointer_send_up[i+j*lattice_->jump(lattice_->dim()-1)*components_];
}
}
parallel.receive_dim1( buffer_rec, buffer_size1,0);
parallel.send_dim1( buffer_send, buffer_size1, 0);
for(j=0;j<lattice_->sizeLocal(lattice_->dim()-1);j++)
{
for(i=0;i<buffer_size1/lattice_->sizeLocal(lattice_->dim()-1);i++)
{
pointer_rec_up[i+j*lattice_->jump(lattice_->dim()-1)*components_] = buffer_rec[i+j*buffer_size1/lattice_->sizeLocal(lattice_->dim()-1)];
}
}
}
}
pointer_send_up = data_ + lattice_->sitesLocalGross() * components_ - 2*buffer_size0;
pointer_rec_up = data_ + lattice_->sitesLocalGross() * components_ - buffer_size0;
pointer_send_down = data_ + buffer_size0;
pointer_rec_down = data_;
if(parallel.grid_rank()[0]%2==0)
{
if(parallel.grid_rank()[0]!=parallel.grid_size()[0]-1)
{
parallel.send_dim0( pointer_send_up, buffer_size0, parallel.grid_rank()[0]+1);
parallel.receive_dim0( pointer_rec_up, buffer_size0, parallel.grid_rank()[0]+1);
}
if(parallel.grid_rank()[0] != 0)
{
parallel.send_dim0( pointer_send_down, buffer_size0, parallel.grid_rank()[0]-1);
parallel.receive_dim0( pointer_rec_down, buffer_size0, parallel.grid_rank()[0]-1);
}
else if(parallel.grid_size()[0]%2==0)
{
parallel.send_dim0( pointer_send_down, buffer_size0, parallel.grid_size()[0]-1);
parallel.receive_dim0( pointer_rec_down, buffer_size0, parallel.grid_size()[0]-1);
}
}
else
{
parallel.receive_dim0( pointer_rec_down, buffer_size0, parallel.grid_rank()[0]-1);
parallel.send_dim0( pointer_send_down, buffer_size0, parallel.grid_rank()[0]-1);
if(parallel.grid_rank()[0]!=parallel.grid_size()[0]-1)
{
parallel.receive_dim0( pointer_rec_up, buffer_size0, parallel.grid_rank()[0]+1);
parallel.send_dim0( pointer_send_up, buffer_size0, parallel.grid_rank()[0]+1);
}
else
{
parallel.receive_dim0( pointer_rec_up, buffer_size0,0);
parallel.send_dim0( pointer_send_up, buffer_size0, 0);
}
}
if(parallel.grid_size()[0]%2!=0)
{
if(parallel.grid_rank()[0]==0)