-
Notifications
You must be signed in to change notification settings - Fork 10
/
fsgrid.hpp
1160 lines (1024 loc) · 46.2 KB
/
fsgrid.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
#pragma once
/*
Copyright (C) 2016 Finnish Meteorological Institute
Copyright (C) 2016-2024 CSC -IT Center for Science
This file is part of fsgrid
fsgrid is free software: you can redistribute it and/or modify
it under the terms of the GNU General Public License as published by
the Free Software Foundation, either version 3 of the License, or
(at your option) any later version.
fsgrid is distributed in the hope that it will be useful,
but WITHOUT ANY WARRANTY;
without even the implied warranty of
MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
GNU General Public License for more details.
You should have received a copy of the GNU General Public License
along with fsgrid. If not, see <http://www.gnu.org/licenses/>.
*/
#include <array>
#include <vector>
#include <mpi.h>
#include <iostream>
#include <limits>
#include <stdint.h>
#include <cassert>
#include <stdio.h>
#include <algorithm>
#ifndef FS_MASTER_RANK
#define FS_MASTER_RANK 0
#endif
struct FsGridTools{
typedef uint32_t FsSize_t; // Size type for global array indices
typedef int32_t FsIndex_t; // Size type for global/local array indices, incl. possible negative values
typedef int64_t LocalID;
typedef int64_t GlobalID;
typedef int Task_t;
//! Helper function: calculate position of the local coordinate space for the given dimension
// \param globalCells Number of cells in the global Simulation, in this dimension
// \param ntasks Total number of tasks in this dimension
// \param my_n This task's position in this dimension
// \return Cell number at which this task's domains cells start (actual cells, not counting ghost cells)
static FsIndex_t calcLocalStart(FsSize_t globalCells, Task_t ntasks, Task_t my_n) {
FsIndex_t n_per_task = globalCells / ntasks;
FsIndex_t remainder = globalCells % ntasks;
if(my_n < remainder) {
return my_n * (n_per_task+1);
} else {
return my_n * n_per_task + remainder;
}
}
//! Helper function: given a global cellID, calculate the global cell coordinate from it.
// This is then used do determine the task responsible for this cell, and the
// local cell index in it.
static std::array<FsIndex_t, 3> globalIDtoCellCoord(GlobalID id, std::array<FsSize_t, 3> &globalSize) {
// Transform globalID to global cell coordinate
std::array<FsIndex_t, 3> cell;
assert(id >= 0);
assert(id < globalSize[0] * globalSize[1] * globalSize[2]);
int stride=1;
for(int i=0; i<3; i++) {
cell[i] = (id / stride) % globalSize[i];
stride *= globalSize[i];
}
return cell;
}
//! Helper function: calculate size of the local coordinate space for the given dimension
// \param globalCells Number of cells in the global Simulation, in this dimension
// \param ntasks Total number of tasks in this dimension
// \param my_n This task's position in this dimension
// \return Number of cells for this task's local domain (actual cells, not counting ghost cells)
static FsIndex_t calcLocalSize(FsSize_t globalCells, Task_t ntasks, Task_t my_n) {
FsIndex_t n_per_task = globalCells/ntasks;
FsIndex_t remainder = globalCells%ntasks;
if(my_n < remainder) {
return n_per_task+1;
} else {
return n_per_task;
}
}
//! Helper function to optimize decomposition of this grid over the given number of tasks
static void computeDomainDecomposition(const std::array<FsSize_t, 3>& GlobalSize, Task_t nProcs, std::array<Task_t,3>& processDomainDecomposition, int stencilSize=1, int verbose = 0) {
int myRank, MPI_flag;
MPI_Initialized(&MPI_flag);
if(MPI_flag){
MPI_Comm_rank(MPI_COMM_WORLD, &myRank); // TODO allow for separate communicator
} else {
myRank = FS_MASTER_RANK;
}
std::array<FsSize_t, 3> systemDim;
std::array<FsSize_t, 3> processBox;
std::array<FsSize_t, 3> minDomainSize;
int64_t optimValue = std::numeric_limits<int64_t>::max();
std::vector<std::pair<int64_t,std::array<Task_t,3>>> scored_decompositions;
scored_decompositions.push_back(std::pair<int64_t,std::array<Task_t,3>>(optimValue,{0,0,0}));
for(int i = 0; i < 3; i++) {
systemDim[i] = GlobalSize[i];
if(GlobalSize[i] == 1) {
// In 2D simulation domains, the "thin" dimension can be a single cell thick.
minDomainSize[i] = 1;
} else {
// Otherwise, it needs to be at least as large as our ghost
// stencil, so that ghost communication remains consistent.
minDomainSize[i] = stencilSize;
}
}
processDomainDecomposition = {1, 1, 1};
for (Task_t i = 1; i <= std::min(nProcs, (Task_t)(GlobalSize[0]/minDomainSize[0])); i++) {
for (Task_t j = 1; j <= std::min(nProcs, (Task_t)(GlobalSize[1]/minDomainSize[1])) ; j++) {
if( i * j > nProcs ){
break;
}
Task_t k = nProcs/(i*j);
// No need to optimize an incompatible DD, also checks for missing remainders
if( i * j * k != nProcs ) {
continue;
}
if (k > (Task_t)(GlobalSize[2]/minDomainSize[2])) {
continue;
}
processBox[0] = calcLocalSize(systemDim[0],i,0);
processBox[1] = calcLocalSize(systemDim[1],j,0);
processBox[2] = calcLocalSize(systemDim[2],k,0);
int64_t value =
(i > 1 ? processBox[1] * processBox[2]: 0) +
(j > 1 ? processBox[0] * processBox[2]: 0) +
(k > 1 ? processBox[0] * processBox[1]: 0);
// account for singular domains
if (i!=1 && j!= 1 && k!=1) {
value *= 13; // 26 neighbours to communicate to
}
if (i==1 && j!= 1 && k!=1) {
value *= 4; // 8 neighbours to communicate to
}
if (i!=1 && j== 1 && k!=1) {
value *= 4; // 8 neighbours to communicate to
}
if (i!=1 && j!= 1 && k==1) {
value *= 4; // 8 neighbours to communicate to
}
// else: 2 neighbours to communicate to, no need to adjust
if(value <= optimValue ){
optimValue = value;
if(value < scored_decompositions.back().first){
scored_decompositions.clear();
}
scored_decompositions.push_back(std::pair<int64_t, std::array<Task_t,3>>(value, {i,j,k}));
}
}
}
if(myRank == FS_MASTER_RANK && verbose){
std::cout << "(FSGRID) Number of equal minimal-surface decompositions found: " << scored_decompositions.size() << "\n";
for (auto kv : scored_decompositions){
std::cout << "(FSGRID) Decomposition " << kv.second[0] <<","<<kv.second[1]<<","<<kv.second[2]<< " "<< " for processBox size " <<
systemDim[0]/kv.second[0] << " " << systemDim[1]/kv.second[1] << " " << systemDim[2]/kv.second[2] <<"\n";
}
}
// Taking the first scored_decomposition (smallest X decomposition)
processDomainDecomposition[0] = scored_decompositions[0].second[0];
processDomainDecomposition[1] = scored_decompositions[0].second[1];
processDomainDecomposition[2] = scored_decompositions[0].second[2];
if(optimValue == std::numeric_limits<int64_t>::max() ||
(Task_t)(processDomainDecomposition[0] * processDomainDecomposition[1] * processDomainDecomposition[2]) != nProcs) {
if(myRank==0){
std::cerr << "(FSGRID) Domain decomposition failed, are you running on a prime number of tasks?" << std::endl;
}
throw std::runtime_error("FSGrid computeDomainDecomposition failed");
}
if(myRank==0 && verbose){
std::cout << "(FSGRID) decomposition chosen as "<< processDomainDecomposition[0] << " " << processDomainDecomposition[1] << " " << processDomainDecomposition[2] << ", for processBox sizes " <<
systemDim[0]/processDomainDecomposition[0] << " " << systemDim[1]/processDomainDecomposition[1] << " " << systemDim[2]/processDomainDecomposition[2] <<
" \n";
}
}
};
/*! Simple cartesian, non-loadbalancing MPI Grid for use with the fieldsolver
*
* \param T datastructure containing the field in each cell which this grid manages
* \param stencil ghost cell width of this grid
*/
template <typename T, int stencil> class FsGrid : public FsGridTools{
template<typename ArrayT> void swapArray(std::array<ArrayT, 3>& array) {
ArrayT a = array[0];
array[0] = array[2];
array[2] = a;
}
public:
/*! Constructor for this grid.
* \param globalSize Cell size of the global simulation domain.
* \param MPI_Comm The MPI communicator this grid should use.
* \param isPeriodic An array specifying, for each dimension, whether it is to be treated as periodic.
*/
FsGrid(std::array<FsSize_t,3> globalSize, MPI_Comm parent_comm, std::array<bool,3> isPeriodic,
const std::array<Task_t, 3>& decomposition = {0,0,0}, bool verbose = false)
: globalSize(globalSize) {
int status;
int size;
// Get parent_comm info
int parentRank;
MPI_Comm_rank(parent_comm, &parentRank);
int parentSize;
MPI_Comm_size(parent_comm, &parentSize);
// If environment variable FSGRID_PROCS is set,
// use that for determining the number of FS-processes
size = parentSize;
if(getenv("FSGRID_PROCS") != NULL) {
const int fsgridProcs = atoi(getenv("FSGRID_PROCS"));
if(fsgridProcs > 0 && fsgridProcs < size)
size = fsgridProcs;
}
std::array<Task_t,3> emptyarr = {0,0,0};
if (decomposition == emptyarr){
// If decomposition isn't pre-defined, heuristically choose a good domain decomposition for our field size
computeDomainDecomposition(globalSize, size, ntasksPerDim, stencil, verbose);
} else {
ntasksPerDim = decomposition;
if (ntasksPerDim[0]*ntasksPerDim[1]*ntasksPerDim[2] != size){
std::cerr << "Given decomposition ("<<ntasksPerDim[0] << " " << ntasksPerDim[1] << " " << ntasksPerDim[2] << ") does not distribute to the number of tasks given" << std::endl;
throw std::runtime_error("Given decomposition does not distribute to the number of tasks given");
}
ntasksPerDim[0] = decomposition[0];
ntasksPerDim[1] = decomposition[1];
ntasksPerDim[2] = decomposition[2];
}
//set private array
periodic = isPeriodic;
//set temporary int arrays for MPI_Cart_create
std::array<int, 3> isPeriodicInt, ntasksInt;
for(unsigned int i=0; i < isPeriodic.size(); i++) {
isPeriodicInt[i] = (int)isPeriodic[i];
ntasksInt[i] = (int)ntasksPerDim[i];
}
// Create a temporary FS subcommunicator for the MPI_Cart_create
int colorFs = (parentRank < size) ? 1 : MPI_UNDEFINED;
MPI_Comm_split(parent_comm, colorFs, parentRank, &comm1d);
if(colorFs != MPI_UNDEFINED){
// Create cartesian communicator. Note, that reorder is false so
// ranks match the ones in parent_comm
status = MPI_Cart_create(comm1d, 3, ntasksPerDim.data(), isPeriodicInt.data(), 0, &comm3d);
if(status != MPI_SUCCESS) {
std::cerr << "Creating cartesian communicatior failed when attempting to create FsGrid!" << std::endl;
throw std::runtime_error("FSGrid communicator setup failed");
}
status = MPI_Comm_rank(comm3d, &rank);
if(status != MPI_SUCCESS) {
std::cerr << "Getting rank failed when attempting to create FsGrid!" << std::endl;
// Without a rank, there's really not much we can do. Just return an uninitialized grid
// (I suppose we'll crash after this, anyway)
return;
}
// Determine our position in the resulting task-grid
status = MPI_Cart_coords(comm3d, rank, 3, taskPosition.data());
if(status != MPI_SUCCESS) {
std::cerr << "Rank " << rank
<< " unable to determine own position in cartesian communicator when attempting to create FsGrid!"
<< std::endl;
}
}
// Create a temporary aux subcommunicator for the (Aux) MPI_Cart_create
int colorAux = (parentRank > (parentSize - 1) % size) ? (parentRank - (parentSize % size)) / size : MPI_UNDEFINED;
MPI_Comm_split(parent_comm, colorAux, parentRank, &comm1d_aux);
int rankAux;
std::array<int, 3> taskPositionAux;
if(colorAux != MPI_UNDEFINED){
// Create an aux cartesian communicator corresponding to the comm3d (but shidted).
status = MPI_Cart_create(comm1d_aux, 3, ntasksPerDim.data(), isPeriodicInt.data(), 0, &comm3d_aux);
if(status != MPI_SUCCESS) {
std::cerr << "Creating cartesian communicatior failed when attempting to create FsGrid!" << std::endl;
throw std::runtime_error("FSGrid communicator setup failed");
}
status = MPI_Comm_rank(comm3d_aux, &rankAux);
if(status != MPI_SUCCESS) {
std::cerr << "Getting rank failed when attempting to create FsGrid!" << std::endl;
return;
}
// Determine our position in the resulting task-grid
status = MPI_Cart_coords(comm3d_aux, rankAux, 3, taskPositionAux.data());
if(status != MPI_SUCCESS) {
std::cerr << "Rank " << rankAux
<< " unable to determine own position in cartesian communicator when attempting to create FsGrid!"
<< std::endl;
}
}
#ifdef FSGRID_DEBUG
// All FS ranks send their true comm3d rank and taskPosition data to dest
MPI_Request *request = new MPI_Request[(parentSize - 1) / size * 2 + 2];
for(int i = 0; i < (parentSize - 1) / size; i++){
int dest = (colorFs != MPI_UNDEFINED) ? parentRank + i * size + (parentSize - 1)
% size + 1 : MPI_PROC_NULL;
if(dest >= parentSize)
dest = MPI_PROC_NULL;
MPI_Isend(&rank, 1, MPI_INT, dest, 9274, parent_comm, &request[2 * i]);
MPI_Isend(taskPosition.data(), 3, MPI_INT, dest, 9275, parent_comm, &request[2 * i + 1]);
}
// All Aux ranks receive the true comm3d rank and taskPosition data from
// source and then compare that it matches their aux data
std::array<int, 3> taskPositionRecv;
int rankRecv;
int source = (colorAux != MPI_UNDEFINED) ? parentRank - (parentRank -
(parentSize % size)) / size * size - parentSize % size : MPI_PROC_NULL;
MPI_Irecv(&rankRecv, 1, MPI_INT, source, 9274, parent_comm, &request[(parentSize - 1) / size * 2]);
MPI_Irecv(taskPositionRecv.data(), 3, MPI_INT, source, 9275, parent_comm,
&request[(parentSize - 1) / size * 2 + 1]);
MPI_Waitall((parentSize - 1) / size * 2 + 2, request, MPI_STATUS_IGNORE);
if(colorAux != MPI_UNDEFINED){
if(rankRecv != rankAux ||
taskPositionRecv[0] != taskPositionAux[0] ||
taskPositionRecv[1] != taskPositionAux[1] ||
taskPositionRecv[2] != taskPositionAux[2]){
std::cerr << "Rank: " << parentRank
<< ". Aux cartesian communicator 'comm3d_aux' does not match with 'comm3d' !" << std::endl;
throw std::runtime_error("FSGrid aux communicator setup failed.");
}
}
delete[] request;
#endif // FSGRID_DEBUG
// Set correct task position for non-FS ranks
if(colorFs == MPI_UNDEFINED){
for(int i=0; i<3; i++){
taskPosition[i] = taskPositionAux[i];
}
}
// Determine size of our local grid
for(int i=0; i<3; i++) {
localSize[i] = calcLocalSize(globalSize[i],ntasksPerDim[i], taskPosition[i]);
localStart[i] = calcLocalStart(globalSize[i],ntasksPerDim[i], taskPosition[i]);
}
if( localSize[0] == 0 || (globalSize[0] > stencil && localSize[0] < stencil)
|| localSize[1] == 0 || (globalSize[1] > stencil && localSize[1] < stencil)
|| localSize[2] == 0 || (globalSize[2] > stencil && localSize[2] < stencil)) {
std::cerr << "FSGrid space partitioning leads to a space that is too small on Rank " << rank << "." <<std::endl;
std::cerr << "Please run with a different number of Tasks, so that space is better divisible." <<std::endl;
throw std::runtime_error("FSGrid too small domains");
}
// NULL-initialize for all procs
neighbourSendType.fill(MPI_DATATYPE_NULL);
neighbourReceiveType.fill(MPI_DATATYPE_NULL);
// If non-FS process, set rank to -1 and localSize to zero and return
if(colorFs == MPI_UNDEFINED){
rank = -1;
localSize[0] = 0;
localSize[1] = 0;
localSize[2] = 0;
comm3d = comm3d_aux;
comm3d_aux = MPI_COMM_NULL;
return;
}
// Allocate the array of neighbours
for(int i=0; i<size; i++) {
neighbour_index.push_back(MPI_PROC_NULL);
}
for(int i=0; i<27; i++) {
neighbour[i]=MPI_PROC_NULL;
}
// Get the IDs of the 26 direct neighbours
for(int x=-1; x<=1;x++) {
for(int y=-1; y<=1;y++) {
for(int z=-1; z<=1; z++) {
std::array<Task_t,3> neighPosition;
/*
* Figure out the coordinates of the neighbours in all three
* directions
*/
neighPosition[0]=taskPosition[0]+x;
if(isPeriodic[0]) {
neighPosition[0] += ntasksPerDim[0];
neighPosition[0] %= ntasksPerDim[0];
}
neighPosition[1]=taskPosition[1]+y;
if(isPeriodic[1]) {
neighPosition[1] += ntasksPerDim[1];
neighPosition[1] %= ntasksPerDim[1];
}
neighPosition[2]=taskPosition[2]+z;
if(isPeriodic[2]) {
neighPosition[2] += ntasksPerDim[2];
neighPosition[2] %= ntasksPerDim[2];
}
/*
* If those coordinates exist, figure out the responsible CPU
* and store its rank
*/
if(neighPosition[0]>=0 && neighPosition[0]<ntasksPerDim[0] && neighPosition[1]>=0
&& neighPosition[1]<ntasksPerDim[1] && neighPosition[2]>=0 && neighPosition[2]<ntasksPerDim[2]) {
// Calculate the rank
int neighRank;
status = MPI_Cart_rank(comm3d, neighPosition.data(), &neighRank);
if(status != MPI_SUCCESS) {
std::cerr << "Rank " << rank << " can't determine neighbour rank at position [";
for(int i=0; i<3; i++) {
std::cerr << neighPosition[i] << ", ";
}
std::cerr << "]" << std::endl;
}
// Forward lookup table
neighbour[(x+1)*9+(y+1)*3+(z+1)]=neighRank;
// Reverse lookup table
if(neighRank >= 0 && neighRank < size) {
neighbour_index[neighRank]=(char) ((x+1)*9+(y+1)*3+(z+1));
}
} else {
neighbour[(x+1)*9+(y+1)*3+(z+1)]=MPI_PROC_NULL;
}
}
}
}
// Allocate local storage array
size_t totalStorageSize=1;
for(int i=0; i<3; i++) {
if(globalSize[i] <= 1) {
// Collapsed dimension => only one cell thick
storageSize[i] = 1;
} else {
// Size of the local domain + 2* size for the ghost cell stencil
storageSize[i] = (localSize[i] + stencil*2);
}
totalStorageSize *= storageSize[i];
}
data.resize(totalStorageSize);
MPI_Datatype mpiTypeT;
MPI_Type_contiguous(sizeof(T), MPI_BYTE, &mpiTypeT);
// Compute send and receive datatypes
//loop through the shifts in the different directions
for(int x=-1; x<=1;x++) {
for(int y=-1; y<=1;y++) {
for(int z=-1; z<=1; z++) {
std::array<int,3> subarraySize;
std::array<int,3> subarrayStart;
const int shiftId = (x+1) * 9 + (y + 1) * 3 + (z + 1);
if((storageSize[0] == 1 && x!= 0 ) ||
(storageSize[1] == 1 && y!= 0 ) ||
(storageSize[2] == 1 && z!= 0 ) ||
(x == 0 && y == 0 && z == 0)){
//skip flat dimension for 2 or 1D simulations, and self
neighbourSendType[shiftId] = MPI_DATATYPE_NULL;
neighbourReceiveType[shiftId] = MPI_DATATYPE_NULL;
continue;
}
subarraySize[0] = (x == 0) ? localSize[0] : stencil;
subarraySize[1] = (y == 0) ? localSize[1] : stencil;
subarraySize[2] = (z == 0) ? localSize[2] : stencil;
if( x == 0 || x == -1 )
subarrayStart[0] = stencil;
else if (x == 1)
subarrayStart[0] = storageSize[0] - 2 * stencil;
if( y == 0 || y == -1 )
subarrayStart[1] = stencil;
else if (y == 1)
subarrayStart[1] = storageSize[1] - 2 * stencil;
if( z == 0 || z == -1 )
subarrayStart[2] = stencil;
else if (z == 1)
subarrayStart[2] = storageSize[2] - 2 * stencil;
for(int i = 0;i < 3; i++)
if(storageSize[i] == 1)
subarrayStart[i] = 0;
std::array<int,3> swappedStorageSize = {(int)storageSize[0],(int)storageSize[1],(int)storageSize[2]};
swapArray(swappedStorageSize);
swapArray(subarraySize);
swapArray(subarrayStart);
MPI_Type_create_subarray(3,
swappedStorageSize.data(),
subarraySize.data(),
subarrayStart.data(),
MPI_ORDER_C,
mpiTypeT,
&(neighbourSendType[shiftId]) );
if(x == 1 )
subarrayStart[0] = 0;
else if (x == 0)
subarrayStart[0] = stencil;
else if (x == -1)
subarrayStart[0] = storageSize[0] - stencil;
if(y == 1 )
subarrayStart[1] = 0;
else if (y == 0)
subarrayStart[1] = stencil;
else if (y == -1)
subarrayStart[1] = storageSize[1] - stencil;
if(z == 1 )
subarrayStart[2] = 0;
else if (z == 0)
subarrayStart[2] = stencil;
else if (z == -1)
subarrayStart[2] = storageSize[2] - stencil;
for(int i = 0;i < 3; i++)
if(storageSize[i] == 1)
subarrayStart[i] = 0;
swapArray(subarrayStart);
MPI_Type_create_subarray(3,
swappedStorageSize.data(),
subarraySize.data(),
subarrayStart.data(),
MPI_ORDER_C,
mpiTypeT,
&(neighbourReceiveType[shiftId]));
}
}
}
for(int i=0;i<27;i++){
if(neighbourReceiveType[i] != MPI_DATATYPE_NULL)
MPI_Type_commit(&(neighbourReceiveType[i]));
if(neighbourSendType[i] != MPI_DATATYPE_NULL)
MPI_Type_commit(&(neighbourSendType[i]));
}
}
std::vector<T>& getData(){
return this->data;
}
void copyData(FsGrid &other){
this->data = other.getData(); // Copy assignment
}
/*!
* MPI calls fail after the main program called MPI_Finalize(),
* so this can be used instead of the destructor
* Cleans up the cartesian communicator and datatypes
*/
void finalize() noexcept {
if (comm3d != MPI_COMM_NULL) {
MPI_Comm_free(&comm3d);
comm3d = MPI_COMM_NULL;
}
if(comm3d_aux != MPI_COMM_NULL) {
MPI_Comm_free(&comm3d_aux);
comm3d_aux = MPI_COMM_NULL;
}
if(comm1d != MPI_COMM_NULL) {
MPI_Comm_free(&comm1d);
comm1d = MPI_COMM_NULL;
}
if(comm1d_aux != MPI_COMM_NULL) {
MPI_Comm_free(&comm1d_aux);
comm1d_aux = MPI_COMM_NULL;
}
for (auto& type : neighbourReceiveType) {
if (type != MPI_DATATYPE_NULL) {
MPI_Type_free(&type);
type = MPI_DATATYPE_NULL;
}
}
for (auto& type : neighbourSendType) {
if (type != MPI_DATATYPE_NULL) {
MPI_Type_free(&type);
type = MPI_DATATYPE_NULL;
}
}
}
/*!
* If finalize() isn't called manually, object should be destroyed before MPI finalization
*/
~FsGrid() {
finalize();
}
friend void swap (FsGrid& first, FsGrid& second) noexcept {
using std::swap;
swap(first.DX, second.DX);
swap(first.DY, second.DY);
swap(first.DZ, second.DZ);
swap(first.physicalGlobalStart, second.physicalGlobalStart);
swap(first.comm1d, second.comm1d);
swap(first.comm1d_aux, second.comm1d_aux);
swap(first.comm3d, second.comm3d);
swap(first.comm3d_aux, second.comm3d_aux);
swap(first.rank, second.rank);
swap(first.requests, second.requests);
swap(first.numRequests, second.numRequests);
swap(first.neighbour, second.neighbour);
swap(first.neighbour_index, second.neighbour_index);
swap(first.ntasksPerDim, second.ntasksPerDim);
swap(first.taskPosition, second.taskPosition);
swap(first.periodic, second.periodic);
swap(first.globalSize, second.globalSize);
swap(first.localSize, second.localSize);
swap(first.storageSize, second.storageSize);
swap(first.localStart, second.localStart);
swap(first.neighbourSendType, second.neighbourSendType);
swap(first.neighbourReceiveType, second.neighbourReceiveType);
swap(first.data, second.data);
}
// Copy constructor
FsGrid(const FsGrid& other) :
DX {other.DX},
DY {other.DY},
DZ {other.DZ},
physicalGlobalStart {other.physicalGlobalStart},
comm1d {MPI_COMM_NULL},
comm1d_aux {MPI_COMM_NULL},
comm3d {MPI_COMM_NULL},
comm3d_aux {MPI_COMM_NULL},
rank {other.rank},
requests {},
numRequests {0},
neighbour {other.neighbour},
neighbour_index {other.neighbour_index},
ntasksPerDim {other.ntasksPerDim},
taskPosition {other.taskPosition},
periodic {other.periodic},
globalSize {other.globalSize},
localSize {other.localSize},
storageSize {other.storageSize},
localStart {other.localStart},
neighbourSendType {},
neighbourReceiveType {},
data {other.data}
{
if (other.comm3d != MPI_COMM_NULL) {
MPI_Comm_dup(other.comm3d, &comm3d);
}
neighbourSendType.fill(MPI_DATATYPE_NULL);
neighbourReceiveType.fill(MPI_DATATYPE_NULL);
for (size_t i = 0; i < neighbourSendType.size(); ++i) {
if (other.neighbourSendType[i] != MPI_DATATYPE_NULL) {
MPI_Type_dup(other.neighbourSendType[i], neighbourSendType.data() + i);
}
if (other.neighbourReceiveType[i] != MPI_DATATYPE_NULL) {
MPI_Type_dup(other.neighbourReceiveType[i], neighbourReceiveType.data() + i);
}
}
}
// Move constructor
// We don't have a default constructor, so just set the MPI stuff NULL
FsGrid(FsGrid&& other) noexcept :
comm1d {MPI_COMM_NULL},
comm1d_aux {MPI_COMM_NULL},
comm3d {MPI_COMM_NULL},
comm3d_aux {MPI_COMM_NULL},
neighbourSendType {},
neighbourReceiveType {}
{
// NULL all the MPI stuff so they won't get freed if destroyed
neighbourSendType.fill(MPI_DATATYPE_NULL);
neighbourReceiveType.fill(MPI_DATATYPE_NULL);
using std::swap;
swap(*this, other);
}
// Copy assignment
// Canonical copy assign is construct to temp and swap, but this would result in an extra allocation of the entire grid
// Copy assignment is currently not used in Vlasiator, so the operator is deleted as any use is likely either erroneous or just a bad idea
FsGrid& operator=(FsGrid other) = delete;
// Move assignment
FsGrid& operator=(FsGrid&& other) noexcept {
using std::swap;
swap(*this, other);
return *this;
}
/*! Returns the task responsible, and its localID for handling the cell with the given GlobalID
* \param id GlobalID of the cell for which task is to be determined
* \return a task for the grid's cartesian communicator
*/
std::pair<int,LocalID> getTaskForGlobalID(GlobalID id) {
// Transform globalID to global cell coordinate
std::array<FsIndex_t, 3> cell = FsGridTools::globalIDtoCellCoord(id, globalSize);
// Find the index in the task grid this Cell belongs to
std::array<int, 3> taskIndex;
for(int i=0; i<3; i++) {
int n_per_task = globalSize[i] / ntasksPerDim[i];
int remainder = globalSize[i] % ntasksPerDim[i];
if(cell[i] < remainder * (n_per_task+1)) {
taskIndex[i] = cell[i] / (n_per_task + 1);
} else {
taskIndex[i] = remainder + (cell[i] - remainder*(n_per_task+1)) / n_per_task;
}
}
// Get the task number from the communicator
std::pair<int,LocalID> retVal;
int status = MPI_Cart_rank(comm3d, taskIndex.data(), &retVal.first);
if(status != MPI_SUCCESS) {
std::cerr << "Unable to find FsGrid rank for global ID " << id << " (coordinates [";
for(int i=0; i<3; i++) {
std::cerr << cell[i] << ", ";
}
std::cerr << "]" << std::endl;
return std::pair<int,LocalID>(MPI_PROC_NULL,0);
}
// Determine localID of that cell within the target task
std::array<FsIndex_t, 3> thatTasksStart;
std::array<FsIndex_t, 3> thatTaskStorageSize;
for(int i=0; i<3; i++) {
thatTasksStart[i] = calcLocalStart(globalSize[i], ntasksPerDim[i], taskIndex[i]);
thatTaskStorageSize[i] = calcLocalSize(globalSize[i], ntasksPerDim[i], taskIndex[i]) + 2 * stencil;
}
retVal.second = 0;
int stride = 1;
for(int i=0; i<3; i++) {
if(globalSize[i] <= 1) {
// Collapsed dimension, doesn't contribute.
retVal.second += 0;
} else {
retVal.second += stride*(cell[i] - thatTasksStart[i] + stencil);
stride *= thatTaskStorageSize[i];
}
}
return retVal;
}
/*! Transform global cell coordinates into the local domain.
* If the coordinates are out of bounds, (-1,-1,-1) is returned.
* \param x The cell's global x coordinate
* \param y The cell's global y coordinate
* \param z The cell's global z coordinate
*/
std::array<FsIndex_t, 3> globalToLocal(FsSize_t x, FsSize_t y, FsSize_t z) {
std::array<FsIndex_t, 3> retval;
retval[0] = (FsIndex_t)x - localStart[0];
retval[1] = (FsIndex_t)y - localStart[1];
retval[2] = (FsIndex_t)z - localStart[2];
if(retval[0] >= localSize[0] || retval[1] >= localSize[1] || retval[2] >= localSize[2]
|| retval[0] < 0 || retval[1] < 0 || retval[2] < 0) {
return {-1,-1,-1};
}
return retval;
}
/*! Determine the cell's GlobalID from its local x,y,z coordinates
* \param x The cell's task-local x coordinate
* \param y The cell's task-local y coordinate
* \param z The cell's task-local z coordinate
*/
GlobalID GlobalIDForCoords(int x, int y, int z) {
return x + localStart[0] + globalSize[0] * (y + localStart[1]) + globalSize[0] * globalSize[1] * (z + localStart[2]);
}
/*! Determine the cell's LocalID from its local x,y,z coordinates
* \param x The cell's task-local x coordinate
* \param y The cell's task-local y coordinate
* \param z The cell's task-local z coordinate
*/
LocalID LocalIDForCoords(int x, int y, int z) {
LocalID index=0;
if(globalSize[2] > 1) {
index += storageSize[0]*storageSize[1]*(stencil+z);
}
if(globalSize[1] > 1) {
index += storageSize[0]*(stencil+y);
}
if(globalSize[0] > 1 ) {
index += stencil+x;
}
return index;
}
/*! Perform ghost cell communication.
*/
void updateGhostCells() {
if(rank == -1) return;
//TODO, faster with simultaneous isends& ireceives?
std::array<MPI_Request, 27> receiveRequests;
std::array<MPI_Request, 27> sendRequests;
for(int i = 0; i < 27; i++){
receiveRequests[i] = MPI_REQUEST_NULL;
sendRequests[i] = MPI_REQUEST_NULL;
}
for(int x=-1; x<=1;x++) {
for(int y=-1; y<=1;y++) {
for(int z=-1; z<=1; z++) {
int shiftId = (x+1) * 9 + (y + 1) * 3 + (z + 1);
int receiveId = (1 - x) * 9 + ( 1 - y) * 3 + ( 1 - z);
if(neighbour[receiveId] != MPI_PROC_NULL &&
neighbourSendType[shiftId] != MPI_DATATYPE_NULL) {
MPI_Irecv(data.data(), 1, neighbourReceiveType[shiftId], neighbour[receiveId], shiftId, comm3d, &(receiveRequests[shiftId]));
}
}
}
}
for(int x=-1; x<=1;x++) {
for(int y=-1; y<=1;y++) {
for(int z=-1; z<=1; z++) {
int shiftId = (x+1) * 9 + (y + 1) * 3 + (z + 1);
int sendId = shiftId;
if(neighbour[sendId] != MPI_PROC_NULL &&
neighbourSendType[shiftId] != MPI_DATATYPE_NULL) {
MPI_Isend(data.data(), 1, neighbourSendType[shiftId], neighbour[sendId], shiftId, comm3d, &(sendRequests[shiftId]));
}
}
}
}
MPI_Waitall(27, receiveRequests.data(), MPI_STATUSES_IGNORE);
MPI_Waitall(27, sendRequests.data(), MPI_STATUSES_IGNORE);
}
/*! Get the size of the local domain handled by this grid.
*/
std::array<FsIndex_t, 3>& getLocalSize() {
return localSize;
}
/*! Get the start coordinates of the local domain handled by this grid.
*/
std::array<FsIndex_t, 3>& getLocalStart() {
return localStart;
}
/*! Get global size of the fsgrid domain
*/
std::array<FsSize_t, 3>& getGlobalSize() {
return globalSize;
}
/*! Calculate global cell position (XYZ in global cell space) from local cell coordinates.
*
* \param x x-Coordinate, in cells
* \param y y-Coordinate, in cells
* \param z z-Coordinate, in cells
*
* \return Global cell coordinates
*/
std::array<FsIndex_t, 3> getGlobalIndices(int64_t x, int64_t y, int64_t z) {
std::array<FsIndex_t, 3> retval;
retval[0] = localStart[0] + x;
retval[1] = localStart[1] + y;
retval[2] = localStart[2] + z;
return retval;
}
/*! Get a reference to the field data in a cell
* \param x x-Coordinate, in cells
* \param y y-Coordinate, in cells
* \param z z-Coordinate, in cells
* \return A reference to cell data in the given cell
*/
T* get(int x, int y, int z) {
// Keep track which neighbour this cell actually belongs to (13 = ourself)
int isInNeighbourDomain=13;
int coord_shift[3] = {0,0,0};
if(x < 0) {
isInNeighbourDomain -= 9;
coord_shift[0] = 1;
}
if(x >= localSize[0]) {
isInNeighbourDomain += 9;
coord_shift[0] = -1;
}
if(y < 0) {
isInNeighbourDomain -= 3;
coord_shift[1] = 1;
}
if(y >= localSize[1]) {
isInNeighbourDomain += 3;
coord_shift[1] = -1;
}
if(z < 0) {
isInNeighbourDomain -= 1;
coord_shift[2] = 1;
}
if(z >= localSize[2]) {
isInNeighbourDomain += 1;
coord_shift[2] = -1;
}
// Santiy-Check that the requested cell is actually inside our domain
// TODO: ugh, this is ugly.
#ifdef FSGRID_DEBUG
bool inside=true;
if(localSize[0] <= 1 && !periodic[0]) {
if(x != 0) {
std::cerr << "x != 0 despite non-periodic x-axis with only one cell." << std::endl;
inside = false;
}
} else {
if(x < -stencil || x >= localSize[0] + stencil) {
std::cerr << "x = " << x << " is outside of [ " << -stencil <<
", " << localSize[0] + stencil << "[!" << std::endl;
inside = false;
}
}
if(localSize[1] <= 1 && !periodic[1]) {
if(y != 0) {
std::cerr << "y != 0 despite non-periodic y-axis with only one cell." << std::endl;
inside = false;
}
} else {
if(y < -stencil || y >= localSize[1] + stencil) {
std::cerr << "y = " << y << " is outside of [ " << -stencil <<
", " << localSize[1] + stencil << "[!" << std::endl;
inside = false;
}
}
if(localSize[2] <= 1 && !periodic[2]) {
if(z != 0) {
std::cerr << "z != 0 despite non-periodic z-axis with only one cell." << std::endl;
inside = false;
}
} else {
if(z < -stencil || z >= localSize[2] + stencil) {
inside = false;
std::cerr << "z = " << z << " is outside of [ " << -stencil <<
", " << localSize[2] + stencil << "[!" << std::endl;
}
}
if(!inside) {
std::cerr << "Out-of bounds access in FsGrid::get! Expect weirdness." << std::endl;
return NULL;
}
#endif // FSGRID_DEBUG
if(isInNeighbourDomain != 13) {