forked from fmihpc/vlasiator
-
Notifications
You must be signed in to change notification settings - Fork 0
/
Copy pathioread.cpp
803 lines (709 loc) · 34.5 KB
/
ioread.cpp
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
/* This file is part of Vlasiator.
*
* Copyright 2010-2015 Finnish Meteorological Institute
*/
#include <cstdlib>
#include <iostream>
#include <iomanip> // for setprecision()
#include <cmath>
#include <vector>
#include <sstream>
#include <ctime>
#include <sys/types.h>
#include <sys/stat.h>
#include "ioread.h"
#include "phiprof.hpp"
#include "parameters.h"
#include "logger.h"
#include "vlsv_reader_parallel.h"
#include "vlasovmover.h"
#include "object_wrapper.h"
using namespace std;
using namespace phiprof;
extern Logger logFile, diagnostic;
typedef Parameters P;
/*!
* \brief Checks for command files written to the local directory.
* If a file STOP was written and is readable, then a bailout with restart writing is initiated.
* If a file KILL was written and is readable, then a bailout without a restart is initiated.
* To avoid bailing out upfront on a new run the files are renamed with the date to keep a trace.
* The function should only be called by MASTER_RANK. This ensures that resetting P::bailout_write_restart works.
*/
void checkExternalCommands() {
struct stat tempStat;
if (stat("STOP", &tempStat) == 0) {
bailout(true, "Received an external STOP command. Setting bailout.write_restart to true.");
P::bailout_write_restart = true;
char newName[80];
// Get the current time.
const time_t rawTime = time(NULL);
const struct tm * timeInfo = localtime(&rawTime);
strftime(newName, 80, "STOP_%F_%H-%M-%S", timeInfo);
rename("STOP", newName);
return;
}
if(stat("KILL", &tempStat) == 0) {
bailout(true, "Received an external KILL command. Setting bailout.write_restart to false.");
P::bailout_write_restart = false;
char newName[80];
// Get the current time.
const time_t rawTime = time(NULL);
const struct tm * timeInfo = localtime(&rawTime);
strftime(newName, 80, "KILL_%F_%H-%M-%S", timeInfo);
rename("KILL", newName);
return;
}
}
/*!
\brief Collective exit on error functions
If any process(es) have a false success values then the program will
abort and write out the message to the logfile
\param success This process' value for successful execution
\param message Error message to be shown upon failure
\param comm MPI comm
\return Returns true if the operation was successful
*/
bool exitOnError(bool success,string message,MPI_Comm comm) {
int successInt;
int globalSuccessInt;
if(success)
successInt=1;
else
successInt=0;
MPI_Allreduce(&successInt,&globalSuccessInt,1,MPI_INT,MPI_MIN,comm);
if(globalSuccessInt==1) {
return true;
}
else{
logFile << message << endl<<write ;
exit(1);
}
}
/*!
\brief Read cell ID's
Read in cell ID's from file. Note: Uses the newer version of vlsv parallel reader
\param file Some vlsv reader with a file open
\param fileCells Vector in whic to store the cell ids
\param masterRank The simulation's master rank id (Vlasiator uses 0, which should be the default)
\param comm MPI comm (MPI_COMM_WORLD should be the default)
*/
bool readCellIds(vlsv::ParallelReader & file,
vector<CellID>& fileCells, const int masterRank,MPI_Comm comm){
// Get info on array containing cell Ids:
uint64_t arraySize = 0;
uint64_t vectorSize;
vlsv::datatype::type dataType;
uint64_t byteSize;
list<pair<string,string> > attribs;
bool success=true;
int rank;
MPI_Comm_rank(comm,&rank);
if (rank==masterRank) {
const short int readFromFirstIndex = 0;
//let's let master read cellId's, we anyway have at max ~1e6 cells
attribs.push_back(make_pair("name","CellID"));
attribs.push_back(make_pair("mesh","SpatialGrid"));
if (file.getArrayInfoMaster("VARIABLE",attribs,arraySize,vectorSize,dataType,byteSize) == false) {
logFile << "(RESTART) ERROR: Failed to read cell ID array info!" << endl << write;
return false;
}
//Make a routine error check:
if( vectorSize != 1 ) {
logFile << "(RESTART) ERROR: Bad vectorsize at " << __FILE__ << " " << __LINE__ << endl << write;
return false;
}
// Read cell Ids:
char* IDbuffer = new char[arraySize*vectorSize*byteSize];
if (file.readArrayMaster("VARIABLE",attribs,readFromFirstIndex,arraySize,IDbuffer) == false) {
logFile << "(RESTART) ERROR: Failed to read cell Ids!" << endl << write;
success = false;
}
// Convert global Ids into our local DCCRG 64 bit uints
const uint64_t& numberOfCells = arraySize;
fileCells.resize(numberOfCells);
if (dataType == vlsv::datatype::type::UINT && byteSize == 4) {
uint32_t* ptr = reinterpret_cast<uint32_t*>(IDbuffer);
//Input cell ids
for (uint64_t i=0; i<numberOfCells; ++i) {
const CellID cellID = ptr[i];
fileCells[i] = cellID;
}
} else if (dataType == vlsv::datatype::type::UINT && byteSize == 8) {
uint64_t* ptr = reinterpret_cast<uint64_t*>(IDbuffer);
for (uint64_t i=0; i<numberOfCells; ++i) {
const CellID cellID = ptr[i];
fileCells[i] = cellID;
}
} else {
logFile << "(RESTART) ERROR: ParallelReader returned an unsupported datatype for cell Ids!" << endl << write;
success = false;
}
delete[] IDbuffer;
}
//broadcast cellId's to everybody
MPI_Bcast(&arraySize,1,MPI_UINT64_T,masterRank,comm);
fileCells.resize(arraySize);
MPI_Bcast(&(fileCells[0]),arraySize,MPI_UINT64_T,masterRank,comm);
return success;
}
/* Read the total number of velocity blocks per spatial cell in the spatial mesh.
* The returned value for each cell is a sum of the velocity blocks associated in each particle species.
* The value is used to calculate an initial load balance after restart.
* @param file Some vlsv reader with a file open (can be old or new vlsv reader)
* @param nBlocks Vector for holding information on cells and the number of blocks in them -- this function saves data here
* @param masterRank The master rank of this process (Vlasiator uses masterRank = 0 and so it should be the default)
* @param comm MPI comm
* @return Returns true if the operation was successful
@ @see exec_readGrid
*/
bool readNBlocks(vlsv::ParallelReader& file,const std::string& meshName,
std::vector<size_t>& nBlocks,int masterRank,MPI_Comm comm) {
bool success = true;
// Get info on array containing cell IDs:
uint64_t arraySize;
uint64_t vectorSize;
vlsv::datatype::type dataType;
uint64_t byteSize;
// Read mesh bounding box to all processes, the info in bbox contains
// the number of spatial cells in the mesh.
uint64_t bbox[6];
uint64_t* bbox_ptr = bbox;
list<pair<string,string> > attribs;
attribs.push_back(make_pair("mesh",meshName));
if (file.read("MESH_BBOX",attribs,0,6,bbox_ptr,false) == false) return false;
// Resize the output vector and init to zero values
const uint64_t N_spatialCells = bbox[0]*bbox[1]*bbox[2];
nBlocks.resize(N_spatialCells);
#pragma omp parallel for
for (size_t i=0; i<nBlocks.size(); ++i) nBlocks[i] = 0;
// Note: the input file contains N particle species, which are also
// defined in the configuration file. We need to read the BLOCKSPERCELL
// array for each species, and sum the values for each spatial cell.
set<string> speciesNames;
if (file.getUniqueAttributeValues("BLOCKSPERCELL","name",speciesNames) == false) return false;
// Iterate over all particle species and read in BLOCKSPERCELL array
// to all processes, and add the values to nBlocks
uint64_t* buffer = new uint64_t[N_spatialCells];
for (set<string>::const_iterator s=speciesNames.begin(); s!=speciesNames.end(); ++s) {
attribs.clear();
attribs.push_back(make_pair("mesh",meshName));
attribs.push_back(make_pair("name",*s));
if (file.getArrayInfo("BLOCKSPERCELL",attribs,arraySize,vectorSize,dataType,byteSize) == false) return false;
if (file.read("BLOCKSPERCELL",attribs,0,arraySize,buffer) == false) {
delete [] buffer; buffer = NULL;
return false;
}
#pragma omp parallel for
for (size_t i=0; i<N_spatialCells; ++i) {
nBlocks[i] += buffer[i];
}
}
delete [] buffer; buffer = NULL;
return success;
}
//Outputs the velocity block indices of some given block into indices
//Input:
//[0] cellStruct -- some cell structure that has been constructed properly
//[1] block -- some velocity block id
//Output:
//[0] indices -- the array where to store the indices
/*! Outputs given block's velocity min coordinates (the corner of the block) in blockCoordinates
\param block The block's id
\param blockCoordinates An empty array where to store the block coordinates
\sa readBlockData
*/
void getVelocityBlockCoordinates(const vmesh::GlobalID& block, boost::array<Real, 3>& blockCoordinates ) {
#warning FIXME restart
cerr << "restart disabled" << endl;
exit(1);
/*
//Get indices:
boost::array<vmesh::LocalID, 3> blockIndices;
blockIndices[0] = block % P::vxblocks_ini;
blockIndices[1] = (block / P::vxblocks_ini) % P::vyblocks_ini;
blockIndices[2] = block / (P::vxblocks_ini * P::vyblocks_ini);
//Store the coordinates:
blockCoordinates[0] = P::vxmin + ((P::vxmax - P::vxmin) / P::vxblocks_ini) * blockIndices[0];
blockCoordinates[1] = P::vymin + ((P::vymax - P::vymin) / P::vyblocks_ini) * blockIndices[1];
blockCoordinates[2] = P::vzmin + ((P::vzmax - P::vzmin) / P::vzblocks_ini) * blockIndices[2];
return;*/
}
/** Read velocity block mesh data and distribution function data belonging to this process
* for the given particle species. This function must be called simultaneously by all processes.
* @param file VLSV reader with input file open.
* @param spatMeshName Name of the spatial mesh.
* @param fileCells List of all spatial cell IDs.
* @param localCellStartOffset The offset from which to start reading cells.
* @param localCells How many spatial cells after the offset to read.
* @param blocksPerCell Number of velocity blocks for this particle species in each spatial cell belonging to this process.
* @param localBlockStartOffset Offset into velocity block data arrays from which to start reading data.
* @param localBlocks Number of velocity blocks for this species assigned to this process.
* @param mpiGrid Parallel grid library.
* @param popID ID of the particle species who's data is to be read.
* @return If true, velocity block data was read successfully.*/
bool _readBlockData(
vlsv::ParallelReader& file,
const std::string& spatMeshName,
const std::vector<uint64_t>& fileCells,
const uint64_t localCellStartOffset,
const uint64_t localCells,
const vmesh::LocalID* blocksPerCell,
const uint64_t localBlockStartOffset,
const uint64_t localBlocks,
dccrg::Dccrg<SpatialCell,dccrg::Cartesian_Geometry>& mpiGrid,
const int& popID
) {
bool success = true;
const string popName = getObjectWrapper().particleSpecies[popID].name;
const string tagName = "BLOCKIDS";
uint64_t arraySize,vectorSize,dataSize;
vlsv::datatype::type datatype;
list<pair<string,string> > attribs;
attribs.push_back(make_pair("mesh",spatMeshName));
attribs.push_back(make_pair("name",popName));
if (file.multiReadStart(tagName,attribs) == false) success = false;
// Read velocity block global IDs
if (mpiGrid.get_rank() == 0) {
if (file.getArrayInfoMaster(tagName,attribs,arraySize,vectorSize,datatype,dataSize) == false) success = false;
if (sizeof(vmesh::GlobalID) != dataSize && datatype != vlsv::datatype::UINT) {
cerr << "(RESTART) ERROR: Datatype conversion not implemented in " << __FILE__ << ":" << __LINE__ << endl;
exit(1);
}
}
vector<vector<vmesh::GlobalID> > blockGIDs(localCells);
for (uint64_t c=0; c<localCells; ++c) {
const vmesh::LocalID N_blocks = blocksPerCell[c];
blockGIDs[c].resize(N_blocks);
char* ptr = reinterpret_cast<char*>(&(blockGIDs[c][0]));
if (file.multiReadAddUnit(N_blocks,ptr) == false) success = false;
}
if (file.multiReadEnd(localBlockStartOffset) == false) success = false;
// Add velocity blocks and calculate block parameters
for (uint64_t c=0; c<localCells; ++c) {
const CellID cellID = fileCells[localCellStartOffset+c];
SpatialCell* cell = mpiGrid[cellID];
cell->add_velocity_blocks(blockGIDs[c],popID);
}
// Deallocate memory
{
vector<vector<vmesh::GlobalID> >().swap(blockGIDs);
}
// Read distribution function data
if (mpiGrid.get_rank() == 0) {
if (file.getArrayInfoMaster("BLOCKVARIABLE",attribs,arraySize,vectorSize,datatype,dataSize) == false) success = false;
if (sizeof(Realf) != dataSize && datatype != vlsv::datatype::FLOAT) {
cerr << "(RESTART) ERROR: Datatype conversion not implemented in " << __FILE__ << ":" << __LINE__ << endl;
exit(1);
}
}
if (file.multiReadStart("BLOCKVARIABLE",attribs) == false) success = false;
for (uint64_t c=0; c<localCells; ++c) {
const CellID cellID = fileCells[localCellStartOffset+c];
SpatialCell* cell = mpiGrid[cellID];
const vmesh::LocalID N_blocks = cell->get_number_of_velocity_blocks(popID);
char* ptr = reinterpret_cast<char*>(cell->get_data(popID));
if (file.multiReadAddUnit(N_blocks,ptr) == false) success = false;
}
if (file.multiReadEnd(localBlockStartOffset) == false) success = false;
return success;
}
/** Read velocity block data of all existing particle species.
* @param file VLSV reader.
* @param meshName Name of the spatial mesh.
* @param fileCells Vector containing spatial cell IDs.
* @param localCellStartOffset Offset into fileCells, determines where the cells belonging
* to this process start.
* @param localCells Number of spatial cells assigned to this process.
* @param mpiGrid Parallel grid library.
* @return If true, velocity block data was read successfully.*/
bool readBlockData(
vlsv::ParallelReader& file,
const string& meshName,
const vector<CellID>& fileCells,
const uint64_t localCellStartOffset,
const uint64_t localCells,
dccrg::Dccrg<SpatialCell,dccrg::Cartesian_Geometry>& mpiGrid
) {
bool success = true;
const uint64_t bytesReadStart = file.getBytesRead();
int N_processes;
MPI_Comm_size(MPI_COMM_WORLD,&N_processes);
uint64_t* offsetArray = new uint64_t[N_processes];
for (int popID=0; popID<getObjectWrapper().particleSpecies.size(); ++popID) {
const string& popName = getObjectWrapper().particleSpecies[popID].name;
list<pair<string,string> > attribs;
attribs.push_back(make_pair("mesh",meshName));
attribs.push_back(make_pair("name",popName));
// In restart files each spatial cell has an entry in CELLSWITHBLOCKS.
// Each process calculates how many velocity blocks it has for this species.
vmesh::LocalID* blocksPerCell = NULL;
if (file.read("BLOCKSPERCELL",attribs,localCellStartOffset,localCells,blocksPerCell,true) == false) {
logFile << "(RESTART) ERROR: Failed to read BLOCKSPERCELL at " << __FILE__ << ":" << __LINE__ << endl << write;
success = false;
}
// Count how many velocity blocks this process gets
uint64_t blockSum = 0;
for (uint64_t i=0; i<localCells; ++i) blockSum += blocksPerCell[i];
// Gather all block sums to master process who will them broadcast
// the values to everyone
MPI_Gather(&blockSum,1,MPI_Type<uint64_t>(),offsetArray,1,MPI_Type<uint64_t>(),0,MPI_COMM_WORLD);
MPI_Bcast(offsetArray,N_processes,MPI_Type<uint64_t>(),0,MPI_COMM_WORLD);
// Calculate the offset from which this process starts reading block data
uint64_t myOffset = 0;
for (uint64_t i=0; i<mpiGrid.get_rank(); ++i) myOffset += offsetArray[i];
if (_readBlockData(file,meshName,fileCells,localCellStartOffset,localCells,blocksPerCell,
myOffset,blockSum,mpiGrid,popID) == false) success = false;
delete [] blocksPerCell; blocksPerCell = NULL;
} // for-loop over particle species
delete [] offsetArray; offsetArray = NULL;
const uint64_t bytesReadEnd = file.getBytesRead() - bytesReadStart;
logFile << "Velocity meshes and data read, approximate data rate is ";
logFile << vlsv::printDataRate(bytesReadEnd,file.getReadTime()) << endl << write;
return success;
}
/*! Reads cell parameters from the file and saves them in the right place in mpiGrid
\param file Some parallel vlsv reader with a file open
\param fileCells List of all cell ids
\param localCellStartOffset Offset in the fileCells list for this process ( calculated so that the amount of blocks is distributed somewhat evenly between processes)
\param localCells The amount of cells to read in this process after localCellStartOffset
\param cellParamsIndex The parameter of the cell index e.g. CellParams::RHO
\param expectedVectorSize The amount of elements in the parameter (parameter can be a scalar or a vector of size N)
\param mpiGrid Vlasiator's grid (the parameters are saved here)
\return Returns true if the operation is successful
*/
template <typename fileReal>
static bool _readCellParamsVariable(
vlsv::ParallelReader& file,
const vector<uint64_t>& fileCells,
const uint64_t localCellStartOffset,
const uint64_t localCells,
const string& variableName,
const size_t cellParamsIndex,
const size_t expectedVectorSize,
dccrg::Dccrg<SpatialCell,dccrg::Cartesian_Geometry>& mpiGrid
) {
uint64_t arraySize;
uint64_t vectorSize;
vlsv::datatype::type dataType;
uint64_t byteSize;
list<pair<string,string> > attribs;
fileReal *buffer;
bool success=true;
attribs.push_back(make_pair("name",variableName));
attribs.push_back(make_pair("mesh","SpatialGrid"));
if (file.getArrayInfo("VARIABLE",attribs,arraySize,vectorSize,dataType,byteSize) == false) {
logFile << "(RESTART) ERROR: Failed to read " << endl << write;
return false;
}
if(vectorSize!=expectedVectorSize){
logFile << "(RESTART) vectorsize wrong " << endl << write;
return false;
}
buffer=new fileReal[vectorSize*localCells];
if(file.readArray("VARIABLE",attribs,localCellStartOffset,localCells,(char *)buffer) == false ) {
logFile << "(RESTART) ERROR: Failed to read " << variableName << endl << write;
return false;
}
for(uint i=0;i<localCells;i++){
uint cell=fileCells[localCellStartOffset+i];
for(uint j=0;j<vectorSize;j++){
mpiGrid[cell]->parameters[cellParamsIndex+j]=buffer[i*vectorSize+j];
}
}
delete[] buffer;
return success;
}
/*! Reads cell parameters from the file and saves them in the right place in mpiGrid
\param file Some parallel vlsv reader with a file open
\param fileCells List of all cell ids
\param localCellStartOffset Offset in the fileCells list for this process ( calculated so that the amount of blocks is distributed somewhat evenly between processes)
\param localCells The amount of cells to read in this process after localCellStartOffset
\param cellParamsIndex The parameter of the cell index e.g. CellParams::RHO
\param expectedVectorSize The amount of elements in the parameter (parameter can be a scalar or a vector of size N)
\param mpiGrid Vlasiator's grid (the parameters are saved here)
\return Returns true if the operation is successful
*/
bool readCellParamsVariable(
vlsv::ParallelReader& file,
const vector<CellID>& fileCells,
const uint64_t localCellStartOffset,
const uint64_t localCells,
const string& variableName,
const size_t cellParamsIndex,
const size_t expectedVectorSize,
dccrg::Dccrg<SpatialCell,dccrg::Cartesian_Geometry>& mpiGrid
) {
uint64_t arraySize;
uint64_t vectorSize;
vlsv::datatype::type dataType;
uint64_t byteSize;
list<pair<string,string> > attribs;
attribs.push_back(make_pair("name",variableName));
attribs.push_back(make_pair("mesh","SpatialGrid"));
if (file.getArrayInfo("VARIABLE",attribs,arraySize,vectorSize,dataType,byteSize) == false) {
logFile << "(RESTART) ERROR: Failed to read " << endl << write;
return false;
}
// Call _readCellParamsVariable
if( dataType == vlsv::datatype::type::FLOAT ) {
switch (byteSize) {
case sizeof(double):
return _readCellParamsVariable<double>( file, fileCells, localCellStartOffset, localCells, variableName, cellParamsIndex, expectedVectorSize, mpiGrid );
break;
case sizeof(float):
return _readCellParamsVariable<float>( file, fileCells, localCellStartOffset, localCells, variableName, cellParamsIndex, expectedVectorSize, mpiGrid );
break;
}
} else if( dataType == vlsv::datatype::type::UINT ) {
switch (byteSize) {
case sizeof(uint32_t):
return _readCellParamsVariable<uint32_t>( file, fileCells, localCellStartOffset, localCells, variableName, cellParamsIndex, expectedVectorSize, mpiGrid );
break;
case sizeof(uint64_t):
return _readCellParamsVariable<uint64_t>( file, fileCells, localCellStartOffset, localCells, variableName, cellParamsIndex, expectedVectorSize, mpiGrid );
break;
}
} else if( dataType == vlsv::datatype::type::INT ) {
switch (byteSize) {
case sizeof(int32_t):
return _readCellParamsVariable<int32_t>( file, fileCells, localCellStartOffset, localCells, variableName, cellParamsIndex, expectedVectorSize, mpiGrid );
break;
case sizeof(int64_t):
return _readCellParamsVariable<int64_t>( file, fileCells, localCellStartOffset, localCells, variableName, cellParamsIndex, expectedVectorSize, mpiGrid );
break;
}
} else {
logFile << "(RESTART) ERROR: Failed to read data type at readCellParamsVariable" << endl << write;
return false;
}
}
/*! A function for reading parameters, e.g., 'timestep'.
\param file VLSV parallel reader with a file open.
\param name Name of the parameter.
\param value Variable in which to store the scalar variable (double, float, int .. ).
\param masterRank The master process' id (Vlasiator uses 0 so this should equal 0 by default).
\param comm MPI comm (MPI_COMM_WORLD should be the default).
\return Returns true if the operation is successful. */
template <typename T>
bool readScalarParameter(vlsv::ParallelReader& file,string name,T& value,int masterRank,MPI_Comm comm) {
if (file.readParameter(name,value) == false) {
logFile << "(RESTART) ERROR: Failed to read parameter '" << name << "' value in ";
logFile << __FILE__ << ":" << __LINE__ << endl << write;
return false;
}
return true;
}
/*! A function for checking the scalar parameter
\param file Some parallel vlsv reader with a file open
\param name Name of the parameter
\param correctValue The correct value of the parameter to compare to
\param masterRank The master process' id (Vlasiator uses 0 so this should be 0 by default)
\param comm MPI comm (Default should be MPI_COMM_WORLD)
\return Returns true if the operation is successful
*/
template <typename T>
bool checkScalarParameter(vlsv::ParallelReader& file,const string& name,T correctValue,int masterRank,MPI_Comm comm) {
T value;
readScalarParameter(file,name,value,masterRank,comm);
if (value != correctValue){
std::ostringstream s;
s << "(RESTART) Parameter " << name << " has mismatching value.";
s << " CFG value = " << correctValue;
s << " Restart file value = " << value;
exitOnError(false,s.str(),MPI_COMM_WORLD);
return false;
}
else{
exitOnError(true,"",MPI_COMM_WORLD);
return true;
}
}
/*!
\brief Read in state from a vlsv file in order to restart simulations
\param mpiGrid Vlasiator's grid
\param name Name of the restart file e.g. "restart.00052.vlsv"
\return Returns true if the operation was successful
\sa readGrid
*/
bool exec_readGrid(dccrg::Dccrg<SpatialCell,dccrg::Cartesian_Geometry>& mpiGrid,
const std::string& name) {
vector<CellID> fileCells; /*< CellIds for all cells in file*/
vector<size_t> nBlocks;/*< Number of blocks for all cells in file*/
bool success=true;
int myRank,processes;
#warning Spatial grid name hard-coded here
const string meshName = "SpatialGrid";
// Attempt to open VLSV file for reading:
MPI_Comm_rank(MPI_COMM_WORLD,&myRank);
MPI_Comm_size(MPI_COMM_WORLD,&processes);
phiprof::start("readGrid");
vlsv::ParallelReader file;
MPI_Info mpiInfo = MPI_INFO_NULL;
if (file.open(name,MPI_COMM_WORLD,MASTER_RANK,mpiInfo) == false) {
success=false;
}
exitOnError(success,"(RESTART) Could not open file",MPI_COMM_WORLD);
// Around May 2015 time was renamed from "t" to "time", we try to read both,
// new way is read first
if (readScalarParameter(file,"time",P::t,MASTER_RANK,MPI_COMM_WORLD) == false)
if (readScalarParameter(file,"t", P::t,MASTER_RANK,MPI_COMM_WORLD) == false)
success=false;
P::t_min=P::t;
// Around May 2015 timestep was renamed from "tstep" to "timestep", we to read
// both, new way is read first
if (readScalarParameter(file,"timestep",P::tstep,MASTER_RANK,MPI_COMM_WORLD) == false)
if (readScalarParameter(file,"tstep", P::tstep,MASTER_RANK,MPI_COMM_WORLD) ==false)
success = false;
P::tstep_min=P::tstep;
if(readScalarParameter(file,"dt",P::dt,MASTER_RANK,MPI_COMM_WORLD) ==false) success=false;
if(readScalarParameter(file,"fieldSolverSubcycles",P::fieldSolverSubcycles,MASTER_RANK,MPI_COMM_WORLD) ==false) {
// Legacy restarts do not have this field, it "should" be safe for one or two steps...
P::fieldSolverSubcycles = 1.0;
std::cout << " No P::fieldSolverSubcycles found in restart, setting 1." << std::endl;
}
MPI_Bcast(&(P::fieldSolverSubcycles),1,MPI_Type<Real>(),MASTER_RANK,MPI_COMM_WORLD);
checkScalarParameter(file,"xmin",P::xmin,MASTER_RANK,MPI_COMM_WORLD);
checkScalarParameter(file,"ymin",P::ymin,MASTER_RANK,MPI_COMM_WORLD);
checkScalarParameter(file,"zmin",P::zmin,MASTER_RANK,MPI_COMM_WORLD);
checkScalarParameter(file,"xmax",P::xmax,MASTER_RANK,MPI_COMM_WORLD);
checkScalarParameter(file,"ymax",P::ymax,MASTER_RANK,MPI_COMM_WORLD);
checkScalarParameter(file,"zmax",P::zmax,MASTER_RANK,MPI_COMM_WORLD);
checkScalarParameter(file,"xcells_ini",P::xcells_ini,MASTER_RANK,MPI_COMM_WORLD);
checkScalarParameter(file,"ycells_ini",P::ycells_ini,MASTER_RANK,MPI_COMM_WORLD);
checkScalarParameter(file,"zcells_ini",P::zcells_ini,MASTER_RANK,MPI_COMM_WORLD);
#warning Vel Mesh Parameters not checked anymore
//checkScalarParameter(file,"vxmin",P::vxmin,MASTER_RANK,MPI_COMM_WORLD);
//checkScalarParameter(file,"vymin",P::vymin,MASTER_RANK,MPI_COMM_WORLD);
//checkScalarParameter(file,"vzmin",P::vzmin,MASTER_RANK,MPI_COMM_WORLD);
//checkScalarParameter(file,"vxmax",P::vxmax,MASTER_RANK,MPI_COMM_WORLD);
//checkScalarParameter(file,"vymax",P::vymax,MASTER_RANK,MPI_COMM_WORLD);
//checkScalarParameter(file,"vzmax",P::vzmax,MASTER_RANK,MPI_COMM_WORLD);
//checkScalarParameter(file,"vxblocks_ini",P::vxblocks_ini,MASTER_RANK,MPI_COMM_WORLD);
//checkScalarParameter(file,"vyblocks_ini",P::vyblocks_ini,MASTER_RANK,MPI_COMM_WORLD);
//checkScalarParameter(file,"vzblocks_ini",P::vzblocks_ini,MASTER_RANK,MPI_COMM_WORLD);
phiprof::start("readDatalayout");
if (success == true) success = readCellIds(file,fileCells,MASTER_RANK,MPI_COMM_WORLD);
// Check that the cellID lists are identical in file and grid
if (myRank==0){
vector<CellID> allGridCells=mpiGrid.get_all_cells();
if (fileCells.size() != allGridCells.size()){
success=false;
}
}
exitOnError(success,"(RESTART) Wrong number of cells in restart file",MPI_COMM_WORLD);
if (success == true) {
success = readNBlocks(file,meshName,nBlocks,MASTER_RANK,MPI_COMM_WORLD);
}
//make sure all cells are empty, we will anyway overwrite everything and
// in that case moving cells is easier...
{
const vector<CellID>& gridCells = getLocalCells();
for (size_t i=0; i<gridCells.size(); i++) {
for (int popID=0; popID<getObjectWrapper().particleSpecies.size(); ++popID)
mpiGrid[gridCells[i]]->clear(popID);
}
}
uint64_t totalNumberOfBlocks=0;
unsigned int numberOfBlocksPerProcess;
for(uint i=0; i<nBlocks.size(); ++i){
totalNumberOfBlocks += nBlocks[i];
}
numberOfBlocksPerProcess= 1 + totalNumberOfBlocks/processes;
uint64_t localCellStartOffset=0; // This is where local cells start in file-list after migration.
uint64_t localCells=0;
uint64_t numberOfBlocksCount=0;
// Pin local cells to remote processes, we try to balance number of blocks so that
// each process has the same amount of blocks, more or less.
for (size_t i=0; i<fileCells.size(); ++i) {
numberOfBlocksCount += nBlocks[i];
int newCellProcess = numberOfBlocksCount/numberOfBlocksPerProcess;
if (newCellProcess == myRank) {
if (localCells == 0)
localCellStartOffset=i; //here local cells start
++localCells;
}
if (mpiGrid.is_local(fileCells[i])) {
mpiGrid.pin(fileCells[i],newCellProcess);
}
}
SpatialCell::set_mpi_transfer_type(Transfer::ALL_DATA);
mpiGrid.balance_load(false);
//update list of local gridcells
recalculateLocalCellsCache();
getObjectWrapper().meshData.reallocate();
//get new list of local gridcells
const vector<CellID>& gridCells = getLocalCells();
// Unpin cells, otherwise we will never change this initial bad balance
for (size_t i=0; i<gridCells.size(); ++i) {
mpiGrid.unpin(gridCells[i]);
}
// Check for errors, has migration succeeded
if (localCells != gridCells.size() ) {
success=false;
}
if (success == true) {
for (uint64_t i=localCellStartOffset; i<localCellStartOffset+localCells; ++i) {
if(mpiGrid.is_local(fileCells[i]) == false) success = false;
}
}
exitOnError(success,"(RESTART) Cell migration failed",MPI_COMM_WORLD);
// Set cell coordinates based on cfg (mpigrid) information
for (size_t i=0; i<gridCells.size(); ++i) {
std::array<double, 3> cell_min = mpiGrid.geometry.get_min(gridCells[i]);
std::array<double, 3> cell_length = mpiGrid.geometry.get_length(gridCells[i]);
mpiGrid[gridCells[i]]->parameters[CellParams::XCRD] = cell_min[0];
mpiGrid[gridCells[i]]->parameters[CellParams::YCRD] = cell_min[1];
mpiGrid[gridCells[i]]->parameters[CellParams::ZCRD] = cell_min[2];
mpiGrid[gridCells[i]]->parameters[CellParams::DX ] = cell_length[0];
mpiGrid[gridCells[i]]->parameters[CellParams::DY ] = cell_length[1];
mpiGrid[gridCells[i]]->parameters[CellParams::DZ ] = cell_length[2];
}
// Where local data start in the blocklists
uint64_t localBlocks=0;
for(uint64_t i=localCellStartOffset; i<localCellStartOffset+localCells; ++i) {
localBlocks += nBlocks[i];
}
phiprof::stop("readDatalayout");
//todo, check file datatype, and do not just use double
phiprof::start("readCellParameters");
if(success) { success=readCellParamsVariable(file,fileCells,localCellStartOffset,localCells,"perturbed_B",CellParams::PERBX,3,mpiGrid); }
// This has to be set anyway, there are also the derivatives that should be written/read if we want to only read in background field
// if(success)
// success=readCellParamsVariable(file,fileCells,localCellStartOffset,localCells,"background_B",CellParams::BGBX,3,mpiGrid);
if(success) { success=readCellParamsVariable(file,fileCells,localCellStartOffset,localCells,"moments",CellParams::RHO,4,mpiGrid); }
if(success) { success=readCellParamsVariable(file,fileCells,localCellStartOffset,localCells,"moments_dt2",CellParams::RHO_DT2,4,mpiGrid); }
if(success) { success=readCellParamsVariable(file,fileCells,localCellStartOffset,localCells,"moments_r",CellParams::RHO_R,4,mpiGrid); }
if(success) { success=readCellParamsVariable(file,fileCells,localCellStartOffset,localCells,"moments_v",CellParams::RHO_V,4,mpiGrid); }
if(success) { success=readCellParamsVariable(file,fileCells,localCellStartOffset,localCells,"pressure",CellParams::P_11,3,mpiGrid); }
if(success) { success=readCellParamsVariable(file,fileCells,localCellStartOffset,localCells,"pressure_dt2",CellParams::P_11_DT2,3,mpiGrid); }
if(success) { success=readCellParamsVariable(file,fileCells,localCellStartOffset,localCells,"pressure_r",CellParams::P_11_R,3,mpiGrid); }
if(success) { success=readCellParamsVariable(file,fileCells,localCellStartOffset,localCells,"pressure_v",CellParams::P_11_V,3,mpiGrid); }
if(success) { success=readCellParamsVariable(file,fileCells,localCellStartOffset,localCells,"LB_weight",CellParams::LBWEIGHTCOUNTER,1,mpiGrid); }
if(success) { success=readCellParamsVariable(file,fileCells,localCellStartOffset,localCells,"max_v_dt",CellParams::MAXVDT,1,mpiGrid); }
if(success) { success=readCellParamsVariable(file,fileCells,localCellStartOffset,localCells,"max_r_dt",CellParams::MAXRDT,1,mpiGrid); }
if(success) { success=readCellParamsVariable(file,fileCells,localCellStartOffset,localCells,"max_fields_dt",CellParams::MAXFDT,1,mpiGrid); }
// Read rho losses Note: vector size = 1 (In the older versions the rho loss wasn't recorded)
if(success) { success=readCellParamsVariable(file,fileCells,localCellStartOffset,localCells,"rho_loss_adjust",CellParams::RHOLOSSADJUST,1,mpiGrid); }
if(success) { success=readCellParamsVariable(file,fileCells,localCellStartOffset,localCells,"rho_loss_velocity_boundary",CellParams::RHOLOSSVELBOUNDARY,1,mpiGrid); }
phiprof::stop("readCellParameters");
phiprof::start("readBlockData");
if (success == true) {
success = readBlockData(file,meshName,fileCells,localCellStartOffset,localCells,mpiGrid);
}
phiprof::stop("readBlockData");
success = file.close();
phiprof::stop("readGrid");
exitOnError(success,"(RESTART) Other failure",MPI_COMM_WORLD);
return success;
}
//FIXME, readGrid has no support for checking or converting endianness
/*!
\brief Read in state from a vlsv file in order to restart simulations
\param mpiGrid Vlasiator's grid
\param name Name of the restart file e.g. "restart.00052.vlsv"
*/
bool readGrid(dccrg::Dccrg<SpatialCell,dccrg::Cartesian_Geometry>& mpiGrid,
const std::string& name){
//Check the vlsv version from the file:
return exec_readGrid(mpiGrid,name);
}