00001 #ifdef CH_LANG_CC
00002
00003
00004
00005
00006
00007
00008
00009 #endif
00010
00011 #ifndef _CH_HDF5_H_
00012 #define _CH_HDF5_H_
00013
00014 #define CHOFFSET(object, member) (int)((char*)&(object.member) - (char*)&object)
00015
00016 #ifdef CH_USE_HDF5 // if you don't have CH_USE_HDF5, then this file is useless
00017
00018 #include <iostream>
00019 using std::cout;
00020 using std::endl;
00021
00022 #ifdef CH_MPI
00023 #include "mpi.h"
00024 #endif
00025
00026 #include "LevelData.H"
00027 #include "HDF5Portable.H"
00028
00029 #undef inline
00030 #include <string>
00031 #include <map>
00032 #include "RealVect.H"
00033 #include "CH_Timer.H"
00034 #include "LoadBalance.H"
00035 #include "LayoutIterator.H"
00036 #include "Vector.H"
00037 #include "memtrack.H"
00038 #include "FluxBox.H"
00039
00040 #ifdef CH_MULTIDIM
00041 #include "BoxTools_ExternC_Mangler.H"
00042 #endif
00043 #include "NamespaceHeader.H"
00044 #ifdef H5_USE_16_API
00045 #define H516
00046 #endif
00047
00048 template <class T>
00049 void read(T& item, Vector<Vector<char> >& a_allocatedBuffers, const Box& box,
00050 const Interval& comps)
00051 {
00052 Vector<void*> b(a_allocatedBuffers.size());
00053 for (int i=0; i<b.size(); ++i) b[i]= &(a_allocatedBuffers[i][0]);
00054 read(item, b, box, comps);
00055 }
00056
00057 namespace CH_HDF5
00058 {
00059 enum IOPolicy
00060 {
00061 IOPolicyDefault = 0,
00062 IOPolicyMultiDimHyperslab = (1<<0),
00063
00064 IOPolicyCollectiveWrite = (1<<1)
00065 };
00066 }
00067
00068 using std::map;
00069
00070 class HDF5Handle;
00071
00072
00073
00074
00075
00076
00077
00078
00079
00080
00081
00082
00083
00084
00085
00086
00087 template <class T>
00088 int writeLevel(HDF5Handle& a_handle,
00089 const int& a_level,
00090 const T& a_data,
00091 const Real& a_dx,
00092 const Real& a_dt,
00093 const Real& a_time,
00094 const Box& a_domain,
00095 const int& a_refRatio,
00096 const IntVect& outputGhost = IntVect::Zero,
00097 const Interval& comps = Interval());
00098
00099 template <class T>
00100 int readLevel(HDF5Handle& a_handle,
00101 const int& a_level,
00102 LevelData<T>& a_data,
00103 Real& a_dx,
00104 Real& a_dt,
00105 Real& a_time,
00106 Box& a_domain,
00107 int& a_refRatio,
00108 const Interval& a_comps = Interval(),
00109 bool setGhost = false);
00110
00111 template <class T>
00112 int readLevel(HDF5Handle& a_handle,
00113 const int& a_level,
00114 LevelData<T>& a_data,
00115 RealVect& a_dx,
00116 Real& a_dt,
00117 Real& a_time,
00118 Box& a_domain,
00119 IntVect& a_refRatio,
00120 const Interval& a_comps = Interval(),
00121 bool setGhost = false);
00122
00123
00124
00125
00126
00127
00128
00129
00130
00131
00132
00133
00134
00135 int write(HDF5Handle& a_handle,
00136 const BoxLayout& a_layout,
00137 const std::string& name = "boxes");
00138
00139
00140
00141
00142
00143
00144
00145 template <class T>
00146 int write(HDF5Handle& a_handle,
00147 const BoxLayoutData<T>& a_data,
00148 const std::string& a_name,
00149 IntVect outputGhost = IntVect::Zero,
00150 const Interval& comps = Interval(),
00151 bool newForm = false);
00152
00153
00154
00155
00156
00157
00158
00159
00160
00161 template <class T>
00162 int write(HDF5Handle& a_handle,
00163 const LevelData<T>& a_data,
00164 const std::string& a_name,
00165 const IntVect& outputGhost = IntVect::Zero,
00166 const Interval& comps = Interval());
00167
00168
00169
00170
00171
00172
00173
00174
00175
00176
00177
00178
00179
00180
00181
00182
00183
00184
00185
00186
00187
00188
00189
00190
00191
00192
00193
00194
00195
00196
00197
00198
00199
00200
00201
00202
00203
00204
00205 int read(HDF5Handle& a_handle,
00206 Vector<Box>& boxes,
00207 const std::string& name = "boxes");
00208
00209
00210
00211
00212
00213
00214
00215
00216
00217 int readBoxes(HDF5Handle& a_handle,
00218 Vector<Vector<Box> >& boxes);
00219
00220
00221
00222
00223
00224
00225
00226
00227
00228
00229 int readFArrayBox(HDF5Handle& a_handle,
00230 FArrayBox& a_fab,
00231 int a_level,
00232 int a_boxNumber,
00233 const Interval& a_components,
00234 const std::string& a_dataName = "data" );
00235
00236
00237
00238
00239
00240
00241
00242
00243 template <class T>
00244 int read(HDF5Handle& a_handle,
00245 BoxLayoutData<T>& a_data,
00246 const std::string& a_name,
00247 const BoxLayout& a_layout,
00248 const Interval& a_comps = Interval(),
00249 bool redefineData = true);
00250
00251
00252
00253
00254
00255
00256
00257
00258
00259
00260
00261
00262
00263
00264
00265
00266
00267
00268
00269
00270
00271
00272
00273 template <class T>
00274 int read(HDF5Handle& a_handle,
00275 LevelData<T>& a_data,
00276 const std::string& a_name,
00277 const DisjointBoxLayout& a_layout,
00278 const Interval& a_comps = Interval(),
00279 bool redefineData = true);
00280
00281
00282
00283
00284
00285
00286
00287
00288
00289
00290
00291
00292 class HDF5Handle
00293 {
00294 public:
00295
00296
00297
00298
00299
00300
00301
00302
00303
00304
00305
00306
00307
00308
00309
00310
00311
00312
00313
00314 enum mode
00315 {
00316 CREATE,
00317 CREATE_SERIAL,
00318 OPEN_RDONLY,
00319 OPEN_RDWR
00320 };
00321
00322
00323
00324
00325
00326
00327
00328
00329 HDF5Handle();
00330
00331
00332
00333
00334
00335
00336
00337
00338
00339
00340
00341
00342
00343
00344 HDF5Handle(
00345 const std::string& a_filename,
00346 mode a_mode,
00347 const char *a_globalGroupName="Chombo_global");
00348
00349
00350
00351 ~HDF5Handle();
00352
00353
00354
00355
00356
00357
00358
00359
00360
00361
00362
00363
00364
00365
00366
00367
00368
00369
00370
00371
00372
00373
00374
00375
00376
00377
00378
00379
00380 int open(
00381 const std::string& a_filename,
00382 mode a_mode,
00383 const char *a_globalGroupName="Chombo_global");
00384
00385
00386
00387
00388
00389
00390
00391 bool isOpen() const;
00392
00393
00394
00395
00396
00397
00398
00399 void close();
00400
00401
00402
00403
00404
00405
00406
00407 void setGroupToLevel(int a_level);
00408
00409
00410
00411
00412
00413
00414
00415
00416
00417 int setGroup(const std::string& groupAbsPath);
00418
00419
00420
00421
00422
00423
00424
00425
00426 int pushGroup( const std::string& grp );
00427
00428
00429
00430
00431
00432
00433
00434
00435
00436 int popGroup();
00437
00438
00439
00440
00441
00442
00443
00444
00445 const std::string& getGroup() const;
00446
00447 HDF5Handle::mode openMode() const {return m_mode;}
00448 const hid_t& fileID() const;
00449 const hid_t& groupID() const;
00450 static hid_t box_id;
00451 static hid_t intvect_id;
00452 static hid_t realvect_id;
00453 static map<std::string, std::string> groups;
00454
00455 private:
00456
00457 HDF5Handle(const HDF5Handle&);
00458 HDF5Handle& operator=(const HDF5Handle&);
00459
00460 hid_t m_fileID;
00461 hid_t m_currentGroupID;
00462 mode m_mode;
00463 bool m_isOpen;
00464 std::string m_filename;
00465 std::string m_group;
00466 int m_level;
00467
00468
00469 static bool initialized;
00470 static void initialize();
00471
00472 };
00473
00474
00475
00476
00477
00478
00479
00480
00481
00482
00483
00484
00485
00486
00487
00488
00489
00490
00491
00492
00493
00494
00495
00496
00497
00498
00499
00500
00501
00502
00503
00504
00505
00506
00507
00508
00509
00510
00511
00512
00513
00514
00515
00516
00517 class HDF5HeaderData
00518 {
00519 public:
00520
00521
00522
00523
00524
00525
00526
00527
00528 int writeToFile(HDF5Handle& file) const;
00529
00530
00531
00532
00533
00534
00535
00536
00537
00538
00539 int readFromFile(HDF5Handle& file);
00540
00541
00542 void clear();
00543
00544
00545 map<std::string, Real> m_real;
00546
00547
00548 map<std::string, int> m_int;
00549
00550
00551 map<std::string, std::string> m_string;
00552
00553
00554 map<std::string, IntVect> m_intvect;
00555
00556
00557 map<std::string, Box> m_box;
00558
00559
00560 map<std::string, RealVect> m_realvect;
00561
00562
00563
00564 int writeToLocation(hid_t loc_id) const;
00565 int readFromLocation(hid_t loc_id);
00566
00567
00568 void dump() const;
00569
00570 private:
00571 static herr_t attributeScan(hid_t loc_id, const char *name, void *opdata);
00572 };
00573
00574 extern "C"
00575 {
00576 #ifdef H516
00577 herr_t HDF5HeaderDataattributeScan(hid_t loc_id, const char *name, void *opdata);
00578 #else
00579 herr_t HDF5HeaderDataattributeScan(hid_t loc_id, const char *name, const H5A_info_t* info, void *opdata);
00580 #endif
00581 }
00582
00583 std::ostream& operator<<(std::ostream& os, const HDF5HeaderData& data);
00584
00585
00586
00587
00588
00589
00590
00591
00592
00593
00594
00595
00596
00597
00598
00599
00600
00601
00602
00603
00604
00605
00606
00607
00608
00609
00610
00611
00612
00613
00614
00615
00616
00617
00618
00619
00620
00621
00622
00623
00624
00625
00626
00627
00628 template <class T>
00629 class WriteMultiData
00630 {
00631
00632 public:
00633
00634
00635 WriteMultiData(HDF5Handle& a_handle,
00636 const BoxLayout& a_layout,
00637 const int a_numIntv,
00638 const string& a_name,
00639 const int a_policyFlags = CH_HDF5::IOPolicyDefault,
00640 const IntVect& a_outputGhost = IntVect::Zero,
00641 const bool a_newForm = false);
00642
00643
00644 ~WriteMultiData()
00645 {
00646 H5Dclose(m_dataSet);
00647 }
00648
00649 private:
00650
00651
00652 WriteMultiData(const WriteMultiData&);
00653 WriteMultiData& operator=(const WriteMultiData);
00654
00655 public:
00656
00657
00658 int writeData(const BoxLayoutData<T>& a_data,
00659 const Interval& a_intvMem,
00660 const Interval& a_intvFile);
00661
00662 protected:
00663
00664 const int m_policyFlags;
00665 const IntVect m_outputGhost;
00666 const bool m_newForm;
00667
00668 char m_dataname[128];
00669 hid_t m_dataSet;
00670 Interval m_allIntvFile;
00671 long m_maxBoxPerProc;
00672
00673
00674
00675 Vector<Vector<long long> > m_offsets;
00676
00677
00678 Vector<hid_t> m_types;
00679 };
00680
00681
00682
00683
00684
00685
00686
00687 #if ( H5_VERS_MAJOR == 1 && H5_VERS_MINOR > 6 )
00688 typedef hsize_t ch_offset_t;
00689 #else
00690 #if ( H5_VERS_MAJOR == 1 && H5_VERS_MINOR == 6 && H5_VERS_RELEASE >= 4 )
00691 typedef hsize_t ch_offset_t;
00692 #else
00693 typedef hssize_t ch_offset_t;
00694 #endif
00695 #endif
00696
00697 template<class T>
00698 hid_t H5Type(const T* dummy);
00699
00700 template< >
00701 hid_t H5Type(const int* dummy);
00702
00703 template< >
00704 hid_t H5Type(const long long* dummy);
00705
00706 template< >
00707 hid_t H5Type(const float* dummy);
00708
00709 template< >
00710 hid_t H5Type(const double* dummy);
00711
00712 template< >
00713 hid_t H5Type(const Box* dummy);
00714
00715 template< >
00716 hid_t H5Type(const RealVect* dummy);
00717
00718 template< >
00719 hid_t H5Type(const IntVect* dummy);
00720
00721 template<class T>
00722 hid_t H5Type(const T* dummy)
00723 {
00724
00725 MayDay::Error(" H5Type(const T* dummy)");
00726 return -4;
00727 }
00728
00729 void createData(hid_t& a_dataset,
00730 hid_t& a_dataspace,
00731 HDF5Handle& handle,
00732 const std::string& name,
00733 hid_t type,
00734 hsize_t size);
00735
00736 template <class T>
00737 void createDataset(hid_t& a_dataset,
00738 hid_t& a_dataspace,
00739 HDF5Handle& handle,
00740 const std::string& name,
00741 const T* dummy,
00742 hsize_t size)
00743 {
00744 createData(a_dataset, a_dataspace, handle, name, H5Type(dummy), size);
00745
00746 }
00747
00748 void writeDataset(hid_t a_dataset,
00749 hid_t a_dataspace,
00750 const void* start,
00751 ch_offset_t off,
00752 hsize_t count);
00753
00754 void readDataset(hid_t a_dataset,
00755 hid_t a_dataspace,
00756 void* start,
00757 ch_offset_t off,
00758 hsize_t count);
00759
00760
00761
00762 struct OffsetBuffer
00763 {
00764 Vector<int> index;
00765 Vector<Vector<int> > offsets;
00766 void operator=(const OffsetBuffer& rhs);
00767 };
00768
00769 ostream& operator<<(ostream& os, const OffsetBuffer& ob);
00770
00771 #include "NamespaceFooter.H"
00772
00773 #include "BaseNamespaceHeader.H"
00774
00775 #include "NamespaceVar.H"
00776
00777
00778 template < >
00779 int linearSize(const CH_XDIR::OffsetBuffer& a_input);
00780
00781
00782 template < >
00783 void linearIn(CH_XDIR::OffsetBuffer& a_outputT, const void* const a_inBuf);
00784
00785
00786 template < >
00787 void linearOut(void* const a_outBuf, const CH_XDIR::OffsetBuffer& a_inputT);
00788
00789 template < > int linearSize(const Vector<CH_XDIR::OffsetBuffer>& a_input);
00790 template < > void linearIn(Vector<CH_XDIR::OffsetBuffer>& a_outputT, const void* const inBuf);
00791 template < > void linearOut(void* const a_outBuf, const Vector<CH_XDIR::OffsetBuffer>& a_inputT);
00792
00793 #include "BaseNamespaceFooter.H"
00794 #include "NamespaceHeader.H"
00795
00796
00797
00798 template <>
00799 inline void dataTypes(Vector<hid_t>& a_types, const BaseFab<int>& dummy)
00800 {
00801 a_types.resize(1);
00802 a_types[0] = H5T_NATIVE_INT;
00803 }
00804
00805 template <>
00806 inline void dataTypes(Vector<hid_t>& a_types, const BaseFab<char>& dummy)
00807 {
00808 a_types.resize(1);
00809 a_types[0] = H5T_NATIVE_CHAR;
00810 }
00811
00812
00813
00814
00815 template <>
00816 inline void dataTypes(Vector<hid_t>& a_types, const FArrayBox& dummy)
00817 {
00818 a_types.resize(1);
00819 a_types[0] = H5T_NATIVE_REAL;
00820 }
00821
00822 template <>
00823 inline void dataTypes(Vector<hid_t>& a_types, const FluxBox& dummy)
00824 {
00825 a_types.resize(1);
00826 a_types[0] = H5T_NATIVE_REAL;
00827 }
00828
00829 template <>
00830 inline void dataSize(const BaseFab<int>& item, Vector<int>& a_sizes,
00831 const Box& box, const Interval& comps)
00832 {
00833 a_sizes[0] = box.numPts() * comps.size();
00834 }
00835
00836 template <>
00837 inline void dataSize(const BaseFab<char>& item, Vector<int>& a_sizes,
00838 const Box& box, const Interval& comps)
00839 {
00840 a_sizes[0] = box.numPts() * comps.size();
00841 }
00842
00843 template <>
00844 inline void dataSize(const FArrayBox& item, Vector<int>& a_sizes,
00845 const Box& box, const Interval& comps)
00846 {
00847 a_sizes[0] = box.numPts() * comps.size();
00848 }
00849
00850 template <>
00851 inline void dataSize(const FluxBox& item, Vector<int>& a_sizes,
00852 const Box& box, const Interval& comps)
00853 {
00854 int& t = a_sizes[0];
00855 t = 0;
00856 for (int dir =0; dir<CH_SPACEDIM; dir++)
00857 {
00858 Box edgeBox(surroundingNodes(box,dir));
00859 t += edgeBox.numPts()*comps.size();
00860 }
00861 }
00862
00863 template <>
00864 inline const char* name(const FArrayBox& a_dummySpecializationArg)
00865 {
00866
00867
00868 const char* name = "FArrayBox";
00869 return name;
00870 }
00871
00872 template <>
00873 inline const char* name(const BaseFab<int>& a_dummySpecializationArg)
00874 {
00875
00876
00877 const char* name = "BaseFab<int>";
00878 return name;
00879 }
00880
00881 template <>
00882 inline const char* name(const BaseFab<char>& a_dummySpecializationArg)
00883 {
00884
00885
00886 const char* name = "BaseFab<char>";
00887 return name;
00888 }
00889
00890
00891
00892
00893 template <class T>
00894 inline void dataTypes(Vector<hid_t>& a_types, const T& dummy)
00895 {
00896 a_types.resize(1);
00897 a_types[0] = H5T_NATIVE_CHAR;
00898 }
00899
00900 template <class T>
00901 inline void dataSize(const T& item, Vector<int>& a_sizes,
00902 const Box& box, const Interval& comps)
00903 {
00904 a_sizes[0] = item.size(box, comps);
00905 }
00906
00907 template <class T>
00908 inline void write(const T& item, Vector<void*>& a_allocatedBuffers,
00909 const Box& box, const Interval& comps)
00910 {
00911 item.linearOut(a_allocatedBuffers[0], box, comps);
00912 }
00913
00914 template <class T>
00915 inline void read(T& item, Vector<void*>& a_allocatedBuffers,
00916 const Box& box, const Interval& comps)
00917 {
00918 item.linearIn(a_allocatedBuffers[0], box, comps);
00919 }
00920
00921 template <class T>
00922 inline const char* name(const T& a_dummySpecializationArg)
00923 {
00924 static const char* name = "unknown";
00925 return name;
00926 }
00927
00928 template <class T>
00929 void getOffsets(Vector<Vector<long long> >& offsets, const BoxLayoutData<T>& a_data,
00930 int types, const Interval& comps, const IntVect& outputGhost)
00931 {
00932 CH_TIME("getOffsets");
00933 const BoxLayout& layout = a_data.boxLayout();
00934 {
00935 CH_TIMELEAF("offsets.resize");
00936 offsets.resize(types, Vector<long long>(layout.size()+1));
00937
00938 }
00939 for (int t=0; t<types; t++) offsets[t][0] = 0;
00940 Vector<int> thisSize(types);
00941 if (T::preAllocatable()==0)
00942 {
00943 T dummy;
00944 unsigned int index = 1;
00945 for (LayoutIterator it(layout.layoutIterator()); it.ok(); ++it)
00946 {
00947 Box region = layout[it()];
00948 region.grow(outputGhost);
00949 {
00950 CH_TIMELEAF("dataSize");
00951 dataSize(dummy, thisSize, region, comps);
00952 }
00953 for (unsigned int i=0; i<thisSize.size(); ++i)
00954 {
00955 CH_TIMELEAF("offsets calc");
00956
00957 offsets[i][index] = offsets[i][index-1] + thisSize[i];
00958 }
00959 ++index;
00960 }
00961 }
00962 else
00963 {
00964 OffsetBuffer buff;
00965
00966 for (DataIterator dit(a_data.dataIterator()); dit.ok(); ++dit)
00967 {
00968 int index = a_data.boxLayout().index(dit());
00969
00970 buff.index.push_back(index);
00971 Box region = layout[dit()];
00972 region.grow(outputGhost);
00973 {
00974 CH_TIMELEAF("dataSize");
00975 dataSize(a_data[dit()], thisSize, region, comps);
00976 }
00977 buff.offsets.push_back(thisSize);
00978 }
00979 Vector<OffsetBuffer> gathering(numProc());
00980 {
00981 CH_TIMELEAF("gather");
00982 gather(gathering, buff, uniqueProc(SerialTask::compute));
00983 }
00984 {
00985 CH_TIMELEAF("broadcast");
00986 broadcast(gathering, uniqueProc(SerialTask::compute));
00987 }
00988
00989 for (int i=0; i<numProc(); ++i)
00990 {
00991 OffsetBuffer& offbuf = gathering[i];
00992 for (int num=0; num<offbuf.index.size(); num++)
00993 {
00994 int index = offbuf.index[num];
00995 for (unsigned int j=0; j<types; ++j)
00996 {
00997 CH_TIMELEAF("offsets calc");
00998
00999 offsets[j][index+1] = offbuf.offsets[num][j];
01000 }
01001 }
01002 }
01003 for (int i=0; i<layout.size(); i++)
01004 {
01005 for (unsigned int j=0; j<types; ++j)
01006 {
01007 CH_TIMELEAF("offsets calc");
01008
01009 offsets[j][i+1] += offsets[j][i];
01010 }
01011 }
01012 }
01013
01014
01015 }
01016
01017
01018
01019
01020
01021
01022 template <class T>
01023 void getOffsets(Vector<Vector<long long> >& a_offsets,
01024 const BoxLayout a_layout,
01025 int a_numTypes,
01026 const Interval& a_comps,
01027 const IntVect& a_outputGhost)
01028 {
01029 CH_TIME("getOffsets (prealloc)");
01030 CH_assert(T::preAllocatable() == 0);
01031
01032 a_offsets.resize(a_numTypes, Vector<long long>(a_layout.size()+1));
01033 for (int t = 0; t != a_numTypes; ++t)
01034 {
01035 a_offsets[t][0] = 0;
01036 }
01037 Vector<int> thisSize(a_numTypes);
01038 T dummy;
01039 unsigned int index = 1;
01040 for (LayoutIterator it(a_layout.layoutIterator()); it.ok(); ++it)
01041 {
01042 Box region = a_layout[it()];
01043 region.grow(a_outputGhost);
01044 dataSize(dummy, thisSize, region, a_comps);
01045 for (int i = 0; i != thisSize.size(); ++i)
01046 {
01047
01048 a_offsets[i][index] = a_offsets[i][index-1] + thisSize[i];
01049 }
01050 ++index;
01051 }
01052 }
01053
01054
01055
01056
01057
01058 template <class T>
01059 int write(HDF5Handle& a_handle, const BoxLayoutData<T>& a_data,
01060 const std::string& a_name, IntVect outputGhost,
01061 const Interval& in_comps, bool newForm)
01062 {
01063 CH_TIME("write_Level");
01064 int ret = 0;
01065
01066 Interval comps(in_comps);
01067 if ( comps.size() == 0) comps = a_data.interval();
01068 T dummy;
01069 Vector<hid_t> types;
01070 dataTypes(types, dummy);
01071
01072 Vector<Vector<long long> > offsets;
01073 Vector<long long> bufferCapacity(types.size(), 1);
01074 Vector<void*> buffers(types.size(), NULL);
01075
01076 getOffsets(offsets, a_data, types.size(), comps, outputGhost);
01077
01078
01079 hsize_t flatdims[1];
01080 char dataname[100];
01081 Vector<hid_t> dataspace(types.size());
01082 Vector<hid_t> dataset(types.size());
01083
01084 herr_t err;
01085 hsize_t count[1];
01086 ch_offset_t offset[1];
01087 CH_assert(!(newForm && types.size() != 1));
01088
01089 for (unsigned int i=0; i<types.size(); ++i)
01090 {
01091 flatdims[0] = offsets[i][offsets[i].size()-1];
01092 if (newForm)
01093 {
01094 if (a_name == "M")
01095 {
01096 strcpy(dataname, "Mask");
01097 }
01098 else
01099 {
01100 sprintf(dataname, "%sRegular", a_name.c_str());
01101 }
01102 }
01103 else
01104 {
01105 sprintf(dataname, "%s:datatype=%i",a_name.c_str(), i);
01106 }
01107 {
01108 CH_TIME("H5Screate");
01109 dataspace[i] = H5Screate_simple(1, flatdims, NULL);
01110 }
01111 CH_assert(dataspace[i] >=0);
01112 {
01113 CH_TIME("H5Dcreate");
01114 #ifdef H516
01115 dataset[i] = H5Dcreate(a_handle.groupID(), dataname,
01116 types[i],
01117 dataspace[i], H5P_DEFAULT);
01118 #else
01119 dataset[i] = H5Dcreate2(a_handle.groupID(), dataname,
01120 types[i],
01121 dataspace[i], H5P_DEFAULT,
01122 H5P_DEFAULT, H5P_DEFAULT);
01123 #endif
01124 }
01125 CH_assert(dataset[i] >= 0);
01126 }
01127
01128 hid_t offsetspace, offsetData;
01129 for (unsigned int i=0; i<types.size(); ++i)
01130 {
01131 flatdims[0] = offsets[i].size();
01132 if (newForm)
01133 {
01134 if (a_name == "M")
01135 {
01136 strcpy(dataname, "MaskOffsets");
01137 }
01138 else
01139 {
01140 sprintf(dataname, "%sOffsets",a_name.c_str());
01141 }
01142 }
01143 else
01144 {
01145 sprintf(dataname, "%s:offsets=%i",a_name.c_str(), i);
01146 }
01147 {
01148 CH_TIME("H5S_H5D_offsets_create");
01149 offsetspace = H5Screate_simple(1, flatdims, NULL);
01150 CH_assert(offsetspace >= 0);
01151 #ifdef H516
01152 offsetData = H5Dcreate(a_handle.groupID(), dataname,
01153 H5T_NATIVE_LLONG, offsetspace,
01154 H5P_DEFAULT);
01155 #else
01156 offsetData = H5Dcreate2(a_handle.groupID(), dataname,
01157 H5T_NATIVE_LLONG, offsetspace,
01158 H5P_DEFAULT,
01159 H5P_DEFAULT,H5P_DEFAULT);
01160 #endif
01161 CH_assert(offsetData >= 0);
01162 }
01163 if (procID() == 0)
01164 {
01165 CH_TIME("WriteOffsets");
01166 hid_t memdataspace = H5Screate_simple(1, flatdims, NULL);
01167 CH_assert(memdataspace >= 0);
01168 err = H5Dwrite(offsetData, H5T_NATIVE_LLONG, memdataspace, offsetspace,
01169 H5P_DEFAULT, &(offsets[i][0]));
01170 CH_assert(err >= 0);
01171 H5Sclose(memdataspace);
01172 }
01173 {
01174 CH_TIME("H5S_H5D_offsets_close");
01175 H5Sclose(offsetspace);
01176 H5Dclose(offsetData);
01177 }
01178 }
01179
01180
01181 if (!newForm)
01182 {
01183 CH_TIME("WriteAttributes");
01184 HDF5HeaderData info;
01185 info.m_int["comps"] = comps.size();
01186 info.m_string["objectType"] = name(dummy);
01187 info.m_intvect["outputGhost"] = outputGhost;
01188 info.m_intvect["ghost"] = outputGhost;
01189 std::string group = a_handle.getGroup();
01190 a_handle.setGroup(group+"/"+a_name+"_attributes");
01191 info.writeToFile(a_handle);
01192 a_handle.setGroup(group);
01193 }
01194
01195
01196
01197
01198 Vector<size_t> type_size(types.size());
01199 for (unsigned int i=0; i<types.size(); ++i)
01200 {
01201 type_size[i] = H5Tget_size(types[i]);
01202 }
01203
01204 Vector<int> thisSize(types.size());
01205
01206
01207
01208
01209
01210
01211
01212
01213
01214
01215 #ifdef CH_MPI
01216
01217 #endif
01218 #ifdef AGGREGATE_BOXES_
01219
01220
01221
01222 Vector<hsize_t> aggCount(types.size(), 0);
01223
01224
01225 Vector<long long> aggBufferSize(types.size(), 0);
01226 Vector<void*> aggBuffers(types.size(), NULL);
01227 #endif
01228
01229
01230
01231
01232
01233
01234
01235
01236
01237 {
01238 CH_TIME("ComputeBufferCapacity");
01239 for (DataIterator it = a_data.dataIterator(); it.ok(); ++it)
01240 {
01241 unsigned int index = a_data.boxLayout().index(it());
01242 for (unsigned int i=0; i<types.size(); ++i)
01243 {
01244 long long size = (offsets[i][index+1] - offsets[i][index])
01245 * type_size[i];
01246
01247 CH_assert(size >= 0);
01248 if (size > bufferCapacity[i])
01249 {
01250 bufferCapacity[i] = size;
01251 }
01252 #ifdef AGGREGATE_BOXES_
01253 aggBufferSize[i] += size;
01254 #endif
01255 }
01256 }
01257 }
01258
01259 for (unsigned int i=0; i<types.size(); ++i)
01260 {
01261 CH_TIME("mallocMT_buffers");
01262 buffers[i] = mallocMT(bufferCapacity[i]);
01263 if (buffers[i] == NULL)
01264 {
01265 pout() << " i=" << i
01266 << " types.size() = " << types.size()
01267 << " bufferCapacity[i] = " << (int)bufferCapacity[i]
01268 << endl;
01269 MayDay::Error("memory error in buffer allocation write");
01270 }
01271 #ifdef AGGREGATE_BOXES_
01272 aggBuffers[i] = mallocMT(aggBufferSize[i]);
01273 if (aggBuffers[i] == NULL)
01274 {
01275 MayDay::Error("memory error in aggregate buffer allocation in write");
01276 }
01277 #endif
01278 }
01279
01280 #ifdef CH_MPI
01281
01282 #endif
01283 #ifdef TRY_MPI_COLLECTIVES_
01284
01285
01286
01287
01288
01289
01290
01291 DataIterator dit = a_data.dataIterator();
01292 int maxNumBoxes = dit.size();
01293 int sendBuf = maxNumBoxes;
01294 int result = MPI_Allreduce(&sendBuf, &maxNumBoxes, 1, MPI_INT,MPI_MAX, Chombo_MPI::comm);
01295 if (result != MPI_SUCCESS)
01296 {
01297 MayDay::Error("Couldn't collect maximum number of boxes!");
01298 }
01299
01300
01301 hid_t DXPL = H5Pcreate(H5P_DATASET_XFER);
01302 if(!(a_handle.openMode()==HDF5Handle::CREATE_SERIAL))
01303 H5Pset_dxpl_mpio(DXPL, H5FD_MPIO_COLLECTIVE);
01304 #endif
01305
01306
01307
01308
01309 {
01310 CH_TIME("linearize_H5Dwrite");
01311 #ifdef AGGREGATE_BOXES_
01312
01313
01314 bool isFirstHyperslab = true;
01315
01316
01317
01318 Vector<void*> bufferLoc(types.size(), NULL);
01319
01320 for(unsigned int i=0; i<types.size(); ++i)
01321 {
01322 bufferLoc[i] = aggBuffers[i];
01323 }
01324 #endif
01325 for (DataIterator it = a_data.dataIterator(); it.ok(); ++it)
01326 {
01327 const T& data = a_data[it()];
01328 unsigned int index = a_data.boxLayout().index(it());
01329 Box box = a_data.box(it());
01330 box.grow(outputGhost);
01331
01332 {
01333 CH_TIMELEAF("linearize");
01334 #ifdef AGGREGATE_BOXES_
01335
01336
01337 for(unsigned int i=0; i<types.size(); i++)
01338 {
01339 data.linearOut(bufferLoc[i], box, comps);
01340 char* tempLoc = ((char*)bufferLoc[i]) + data.size(box,comps);
01341 bufferLoc[i] = (void*)tempLoc;
01342
01343 }
01344 #else
01345 write(data, buffers, box, comps);
01346 #endif
01347 }
01348
01349
01350
01351
01352
01353 for (unsigned int i=0; i<types.size(); ++i)
01354 {
01355 offset[0] = offsets[i][index];
01356 count[0] = offsets[i][index+1] - offset[0];
01357 #ifdef AGGREGATE_BOXES_
01358
01359
01360
01361
01362
01363
01364 aggCount[i] += count[0];
01365 if (isFirstHyperslab)
01366 {
01367
01368
01369
01370
01371 if(count[0] > 0)
01372 {
01373 err = H5Sselect_hyperslab(dataspace[i], H5S_SELECT_SET,
01374 offset, NULL,
01375 count, NULL);
01376 CH_assert(err >= 0);
01377 isFirstHyperslab = false;
01378 }
01379 else
01380 {
01381 H5Sselect_none(dataspace[i]);
01382 }
01383 }
01384 else
01385 {
01386 if(count[0] > 0)
01387 {
01388
01389
01390
01391 err = H5Sselect_hyperslab(dataspace[i], H5S_SELECT_OR,
01392 offset, NULL,
01393 count, NULL);
01394 CH_assert(err >= 0);
01395 }
01396
01397
01398 }
01399 #else
01400
01401 hid_t memdataspace=0;
01402 if (count[0] > 0)
01403 {
01404 err = H5Sselect_hyperslab(dataspace[i], H5S_SELECT_SET,
01405 offset, NULL,
01406 count, NULL);
01407 CH_assert(err >= 0);
01408 memdataspace = H5Screate_simple(1, count, NULL);
01409 CH_assert(memdataspace >= 0);
01410 }
01411 else
01412 {
01413 H5Sselect_none(dataspace[i]);
01414 H5Sselect_none(memdataspace);
01415 }
01416
01417
01418
01419
01420 {
01421 CH_TIMELEAF("H5Dwrite");
01422 #ifdef TRY_MPI_COLLECTIVES_
01423 err = H5Dwrite(dataset[i], types[i], memdataspace, dataspace[i],
01424 DXPL, buffers[i]);
01425 #else
01426 err = H5Dwrite(dataset[i], types[i], memdataspace, dataspace[i],
01427 H5P_DEFAULT, buffers[i]);
01428 #endif
01429 }
01430 CH_assert(err >= 0);
01431 H5Sclose(memdataspace);
01432 if (err < 0)
01433 {
01434 ret = err;
01435 pout() << "Before goto cleanup" << endl;
01436 goto cleanup;
01437 }
01438 #endif // end of ifdef AGGREGATE_BOXES_
01439 }
01440 }
01441
01442
01443
01444
01445
01446 #ifdef AGGREGATE_BOXES_
01447 for(unsigned int i=0; i<types.size(); ++i)
01448 {
01449 if(aggCount[i] > 0)
01450 {
01451 hid_t memdataspace = 0;
01452 memdataspace = H5Screate_simple(1, &(aggCount[i]), NULL);
01453 CH_assert(memdataspace >= 0);
01454 {
01455 CH_TIMELEAF("H5Dwrite");
01456 #ifdef TRY_MPI_COLLECTIVES_
01457 err = H5Dwrite(dataset[i], types[i], memdataspace, dataspace[i],
01458 DXPL, aggBuffers[i]);
01459 #else
01460 err = H5Dwrite(dataset[i], types[i], memdataspace, dataspace[i],
01461 H5P_DEFAULT, aggBuffers[i]);
01462 #endif
01463 }
01464 H5Sclose(memdataspace);
01465 if (err < 0)
01466 {
01467 ret = err;
01468 pout() << "Error! goto cleanup" << endl;
01469 goto cleanup;
01470 }
01471 }
01472
01473
01474 #ifdef TRY_MPI_COLLECTIVES_
01475 else
01476 {
01477 hid_t memdataspace = 0;
01478 memdataspace = H5Screate_simple(1, &(aggCount[i]), NULL);
01479 H5Sselect_none(memdataspace);
01480 H5Sselect_none(dataspace[i]);
01481 err = H5Dwrite(dataset[i], types[i], memdataspace, dataspace[i],
01482 DXPL, aggBuffers[i]);
01483 if (err < 0)
01484 {
01485 ret = err;
01486 pout() << "Before goto cleanup" << endl;
01487 goto cleanup;
01488 }
01489 }
01490 #endif
01491 }
01492 #else // not using aggregated collective buffering
01493 #ifdef TRY_MPI_COLLECTIVES_
01494
01495
01496
01497
01498 hid_t memdataspace = 0;
01499
01500
01501 memdataspace = H5Screate_simple(1, count, NULL);
01502 H5Sselect_none(memdataspace);
01503 int nBoxes = a_data.dataIterator().size();
01504 for(int iwrite = nBoxes; iwrite < maxNumBoxes; iwrite++)
01505 {
01506 for (unsigned int i=0; i<types.size(); ++i)
01507 {
01508 H5Sselect_none(dataspace[i]);
01509
01510
01511 err = H5Dwrite(dataset[i], types[i], memdataspace, dataspace[i],
01512 DXPL, buffers[i]);
01513 if (err < 0)
01514 {
01515 ret = err;
01516 pout() << "Before goto cleanup" << endl;
01517 goto cleanup;
01518 }
01519 }
01520 }
01521 H5Sclose(memdataspace);
01522
01523 #endif // end of #ifdef TRY_MPI_COLLECTIVES_
01524 #endif // end of #ifdef AGGREGATE_BOXES_
01525
01526 }
01527
01528 #ifdef TRY_MPI_COLLECTIVES_
01529 H5Pclose(DXPL);
01530 #endif
01531
01532
01533
01534 cleanup:
01535 for (unsigned int i=0; i<types.size(); ++i)
01536 {
01537 {
01538 CH_TIME("freeMT");
01539 freeMT(buffers[i]);
01540 #ifdef AGGREGATE_BOXES_
01541 freeMT(aggBuffers[i]);
01542 #endif
01543 }
01544 {
01545 CH_TIME("H5Sclose");
01546 H5Sclose(dataspace[i]);
01547 }
01548 {
01549 CH_TIME("H5Dclose");
01550 H5Dclose(dataset[i]);
01551 }
01552 }
01553 return ret;
01554
01555 }
01556
01557 template <class T>
01558 int write(HDF5Handle& a_handle, const LevelData<T>& a_data,
01559 const std::string& a_name, const IntVect& outputGhost, const Interval& in_comps)
01560 {
01561 CH_TIMERS("Write Level");
01562 CH_TIMER("calc minimum in outputGhost",t1);
01563 CH_TIMER("writeToFile",t2);
01564 CH_TIMER("setGroup",t3);
01565 HDF5HeaderData info;
01566 info.m_intvect["ghost"] = a_data.ghostVect();
01567 IntVect og(outputGhost);
01568 CH_START(t1);
01569 og.min(a_data.ghostVect());
01570 CH_STOP(t1);
01571 info.m_intvect["outputGhost"] = og;
01572 std::string group = a_handle.getGroup();
01573 a_handle.setGroup(group+"/"+a_name+"_attributes");
01574 CH_START(t2);
01575 info.writeToFile(a_handle);
01576 CH_STOP(t2);
01577 CH_START(t3);
01578 a_handle.setGroup(group);
01579 CH_STOP(t3);
01580 return write(a_handle, (const BoxLayoutData<T>&)a_data, a_name, og, in_comps);
01581 }
01582
01583 template <class T>
01584 int read(HDF5Handle& a_handle, LevelData<T>& a_data, const std::string& a_name,
01585 const DisjointBoxLayout& a_layout, const Interval& a_comps, bool a_redefineData)
01586 {
01587 if (a_redefineData)
01588 {
01589 HDF5HeaderData info;
01590 std::string group = a_handle.getGroup();
01591 if (a_handle.setGroup(group+"/"+a_name+"_attributes"))
01592 {
01593 std::string message = "error opening "
01594 +a_handle.getGroup()+"/"+a_name+"_attributes" ;
01595 MayDay::Warning(message.c_str());
01596 return 1;
01597 }
01598 info.readFromFile(a_handle);
01599 a_handle.setGroup(group);
01600 int ncomp = info.m_int["comps"];
01601 IntVect ghost = info.m_intvect["ghost"];
01602 if (a_comps.end() > 0 && ncomp < a_comps.end())
01603 {
01604 MayDay::Error("attempt to read component interval that is not available");
01605 }
01606 if (a_comps.size() == 0)
01607 a_data.define(a_layout, ncomp, ghost);
01608 else
01609 a_data.define(a_layout, a_comps.size(), ghost);
01610 }
01611 return read(a_handle, (BoxLayoutData<T>&)a_data, a_name, a_layout, a_comps, false);
01612
01613 }
01614
01615 template <class T>
01616 int read(HDF5Handle& a_handle, BoxLayoutData<T>& a_data, const std::string& a_name,
01617 const BoxLayout& a_layout, const Interval& a_comps, bool a_redefineData)
01618 {
01619 CH_TIME("read");
01620 int ret = 0;
01621
01622 herr_t err;
01623
01624 char dataname[100];
01625 hsize_t count[1];
01626 ch_offset_t offset[1];
01627 Vector<Vector<long long> > offsets;
01628
01629 T dummy;
01630 Vector<hid_t> types;
01631 dataTypes(types, dummy);
01632 Vector<hid_t> dataspace(types.size());
01633 Vector<hid_t> dataset(types.size());
01634 offsets.resize(types.size(), Vector<long long>(a_layout.size() +1));
01635
01636
01637
01638 Vector<Vector<char> > buffers(types.size(), 500);
01639
01640 for (unsigned int i=0; i<types.size(); ++i)
01641 {
01642 sprintf(dataname, "%s:datatype=%i",a_name.c_str(), i);
01643 #ifdef H516
01644 dataset[i] = H5Dopen(a_handle.groupID(), dataname);
01645 #else
01646 dataset[i] = H5Dopen2(a_handle.groupID(), dataname, H5P_DEFAULT);
01647 #endif
01648 if (dataset[i] < 0)
01649 {
01650 MayDay::Warning("dataset open failure"); return dataset[i];
01651 }
01652 dataspace[i] = H5Dget_space(dataset[i]);
01653 if (dataspace[i] < 0)
01654 {
01655 MayDay::Warning("dataspace open failure"); return dataspace[i];
01656 }
01657 }
01658
01659 hid_t offsetspace, offsetData;
01660 hsize_t flatdims[1];
01661 for (unsigned int i=0; i<types.size(); ++i)
01662 {
01663 flatdims[0] = offsets[i].size();
01664 sprintf(dataname, "%s:offsets=%i",a_name.c_str(), i);
01665 offsetspace = H5Screate_simple(1, flatdims, NULL);
01666 CH_assert(offsetspace >= 0);
01667 #ifdef H516
01668 offsetData = H5Dopen(a_handle.groupID(), dataname);
01669 #else
01670 offsetData = H5Dopen2(a_handle.groupID(), dataname,H5P_DEFAULT);
01671 #endif
01672 CH_assert(offsetData >= 0);
01673 hid_t memdataspace = H5Screate_simple(1, flatdims, NULL);
01674 CH_assert(memdataspace >= 0);
01675 err = H5Dread(offsetData, H5T_NATIVE_LLONG, memdataspace, offsetspace,
01676 H5P_DEFAULT, &(offsets[i][0]));
01677 CH_assert(err >=0);
01678 H5Sclose(memdataspace);
01679 H5Sclose(offsetspace);
01680 H5Dclose(offsetData);
01681 }
01682
01683 HDF5HeaderData info;
01684 std::string group = a_handle.getGroup();
01685 if (a_handle.setGroup(a_handle.getGroup()+"/"+a_name+"_attributes"))
01686 {
01687 std::string message = "error opening "+a_handle.getGroup()+"/"+a_name ;
01688 MayDay::Warning(message.c_str());
01689 return 1;
01690 }
01691
01692 info.readFromFile(a_handle);
01693 a_handle.setGroup(group);
01694 int ncomps = info.m_int["comps"];
01695 IntVect outputGhost(IntVect::Zero);
01696 if (info.m_intvect.find("outputGhost") != info.m_intvect.end())
01697 {
01698 outputGhost = info.m_intvect["outputGhost"];
01699 }
01700 if (ncomps <= 0)
01701 {
01702 MayDay::Warning("ncomps <= 0 in read");
01703 return ncomps;
01704 }
01705
01706 if (a_redefineData)
01707 {
01708 if (a_comps.size() != 0)
01709 {
01710 a_data.define(a_layout, a_comps.size());
01711 }
01712 else
01713 {
01714 a_data.define(a_layout, ncomps);
01715 }
01716 }
01717
01718 Interval comps(0, ncomps-1);
01719
01720
01721
01722
01723
01724 Vector<size_t> type_size(types.size());
01725 for (unsigned int i=0; i<types.size(); ++i)
01726 {
01727 type_size[i] = H5Tget_size(types[i]);
01728
01729 }
01730
01731 Vector<int> thisSize(types.size());
01732 DataIterator it = a_data.dataIterator();
01733 for ( ; it.ok(); ++it)
01734 {
01735 CH_TIMELEAF("H5Dread");
01736 T& data = a_data[it()];
01737 unsigned int index = a_data.boxLayout().index(it());
01738 Box box = a_data.box(it());
01739
01740 for (unsigned int i=0; i<types.size(); ++i)
01741 {
01742
01743 offset[0] = offsets[i][index];
01744 count[0] = offsets[i][index+1] - offset[0];
01745 if (count[0] > 0)
01746 {
01747 size_t size = count[0] * type_size[i];
01748 while (size > buffers[i].size())
01749 {
01750 buffers[i].resize(2*buffers[i].size());
01751 }
01752
01753 err = H5Sselect_hyperslab(dataspace[i], H5S_SELECT_SET,
01754 offset, NULL,
01755 count, NULL);
01756 CH_assert(err >= 0);
01757 hid_t memdataspace = H5Screate_simple(1, count, NULL);
01758 CH_assert(memdataspace >= 0);
01759 err = H5Dread(dataset[i], types[i], memdataspace, dataspace[i],
01760 H5P_DEFAULT, &(buffers[i][0]));
01761 CH_assert(err >= 0);
01762 H5Sclose(memdataspace);
01763 if (err < 0)
01764 {
01765 ret = err;
01766 goto cleanup;
01767 }
01768 }
01769 }
01770 box.grow(outputGhost);
01771 read(data, buffers, box, comps);
01772 }
01773
01774
01775
01776
01777
01778
01779
01780
01781
01782
01783
01784
01785
01786
01787
01788 cleanup:
01789 for (unsigned int i=0; i<types.size(); ++i)
01790 {
01791
01792 H5Sclose(dataspace[i]);
01793 H5Dclose(dataset[i]);
01794 }
01795 return ret;
01796 }
01797
01798 template <class T>
01799 int writeLevel(HDF5Handle& a_handle,
01800 const int& a_level,
01801 const T& a_data,
01802 const Real& a_dx,
01803 const Real& a_dt,
01804 const Real& a_time,
01805 const Box& a_domain,
01806 const int& a_refRatio,
01807 const IntVect& outputGhost,
01808 const Interval& comps)
01809 {
01810 int error;
01811 char levelName[10];
01812 std::string currentGroup = a_handle.getGroup();
01813 sprintf(levelName, "/level_%i",a_level);
01814 error = a_handle.setGroup(currentGroup + levelName);
01815 if (error != 0) return 1;
01816
01817 HDF5HeaderData meta;
01818 meta.m_real["dx"] = a_dx;
01819 meta.m_real["dt"] = a_dt;
01820 meta.m_real["time"] = a_time;
01821 meta.m_box["prob_domain"] = a_domain;
01822 meta.m_int["ref_ratio"] = a_refRatio;
01823
01824 error = meta.writeToFile(a_handle);
01825 if (error != 0) return 2;
01826
01827 error = write(a_handle, a_data.boxLayout());
01828 if (error != 0) return 3;
01829
01830 error = write(a_handle, a_data, "data", outputGhost, comps);
01831 if (error != 0) return 4;
01832
01833 a_handle.setGroup(currentGroup);
01834
01835 return 0;
01836 }
01837
01838 template <class T>
01839 int writeLevel(HDF5Handle& a_handle,
01840 const int& a_level,
01841 const T& a_data,
01842 const RealVect& a_dx,
01843 const Real& a_dt,
01844 const Real& a_time,
01845 const Box& a_domain,
01846 const IntVect& a_refRatios,
01847 const IntVect& outputGhost,
01848 const Interval& comps)
01849 {
01850 int error;
01851 char levelName[10];
01852 std::string currentGroup = a_handle.getGroup();
01853 sprintf(levelName, "/level_%i",a_level);
01854 error = a_handle.setGroup(currentGroup + levelName);
01855 if (error != 0) return 1;
01856
01857 HDF5HeaderData meta;
01858 meta.m_realvect["vec_dx"] = a_dx;
01859 meta.m_real["dt"] = a_dt;
01860 meta.m_real["time"] = a_time;
01861 meta.m_box["prob_domain"] = a_domain;
01862 meta.m_intvect["vec_ref_ratio"] = a_refRatios;
01863
01864 error = meta.writeToFile(a_handle);
01865 if (error != 0) return 2;
01866
01867 error = write(a_handle, a_data.boxLayout());
01868 if (error != 0) return 3;
01869
01870 error = write(a_handle, a_data, "data", outputGhost, comps);
01871 if (error != 0) return 4;
01872
01873 a_handle.setGroup(currentGroup);
01874
01875 return 0;
01876 }
01877
01878 template <class T>
01879 int readLevel(HDF5Handle& a_handle,
01880 const int& a_level,
01881 LevelData<T>& a_data,
01882 Real& a_dx,
01883 Real& a_dt,
01884 Real& a_time,
01885 Box& a_domain,
01886 int& a_refRatio,
01887 const Interval& a_comps,
01888 bool setGhost)
01889 {
01890 HDF5HeaderData header;
01891 header.readFromFile(a_handle);
01892
01893
01894
01895 int error;
01896 char levelName[10];
01897 std::string currentGroup = a_handle.getGroup();
01898 sprintf(levelName, "/level_%i",a_level);
01899 error = a_handle.setGroup(currentGroup + levelName);
01900 if (error != 0) return 1;
01901
01902 HDF5HeaderData meta;
01903 error = meta.readFromFile(a_handle);
01904 if (error != 0) return 2;
01905 a_dx = meta.m_real["dx"];
01906 a_dt = meta.m_real["dt"];
01907 a_time = meta.m_real["time"];
01908 a_domain = meta.m_box["prob_domain"];
01909 a_refRatio = meta.m_int["ref_ratio"];
01910 Vector<Box> boxes;
01911 error = read(a_handle, boxes);
01912 Vector<int> procIDs;
01913 LoadBalance(procIDs, boxes);
01914
01915 DisjointBoxLayout layout(boxes, procIDs, a_domain);
01916
01917 layout.close();
01918 if (error != 0) return 3;
01919
01920 error = read(a_handle, a_data, "data", layout, a_comps, true);
01921
01922 if (error != 0) return 4;
01923
01924 a_handle.setGroup(currentGroup);
01925
01926 return 0;
01927 }
01928
01929 template <class T>
01930 int readLevel(HDF5Handle& a_handle,
01931 const int& a_level,
01932 LevelData<T>& a_data,
01933 RealVect& a_dx,
01934 Real& a_dt,
01935 Real& a_time,
01936 Box& a_domain,
01937 IntVect& a_refRatio,
01938 const Interval& a_comps,
01939 bool setGhost)
01940 {
01941 HDF5HeaderData header;
01942 header.readFromFile(a_handle);
01943
01944
01945
01946 int error;
01947 char levelName[10];
01948 std::string currentGroup = a_handle.getGroup();
01949 sprintf(levelName, "/level_%i",a_level);
01950 error = a_handle.setGroup(currentGroup + levelName);
01951 if (error != 0) return 1;
01952
01953 HDF5HeaderData meta;
01954 error = meta.readFromFile(a_handle);
01955 if (error != 0) return 2;
01956
01957
01958 if (meta.m_realvect.find("vec_dx") != meta.m_realvect.end())
01959 {
01960 a_dx = meta.m_realvect["vec_dx"];
01961 }
01962 else
01963 {
01964 Real dxScalar = meta.m_real["dx"];
01965 a_dx = dxScalar * RealVect::Unit;
01966 }
01967 a_dt = meta.m_real["dt"];
01968 a_time = meta.m_real["time"];
01969 a_domain = meta.m_box["prob_domain"];
01970
01971
01972 if (meta.m_intvect.find("vec_ref_ratio") != meta.m_intvect.end())
01973 {
01974 a_refRatio = meta.m_intvect["vec_ref_ratio"];
01975 }
01976 else
01977 {
01978 int refRatioScalar = meta.m_int["ref_ratio"];
01979 a_refRatio = refRatioScalar * IntVect::Unit;
01980 }
01981 Vector<Box> boxes;
01982 error = read(a_handle, boxes);
01983 Vector<int> procIDs;
01984 LoadBalance(procIDs, boxes);
01985
01986 DisjointBoxLayout layout(boxes, procIDs, a_domain);
01987
01988 layout.close();
01989 if (error != 0) return 3;
01990
01991 error = read(a_handle, a_data, "data", layout, a_comps, true);
01992
01993 if (error != 0) return 4;
01994
01995 a_handle.setGroup(currentGroup);
01996
01997 return 0;
01998 }
01999
02000
02001
02002
02003
02004
02005
02006
02007
02008
02009
02010
02011
02012
02013
02014
02015
02016
02017
02018
02019
02020
02021
02022
02023
02024
02025
02026
02027
02028
02029
02030
02031
02032
02033
02034
02035
02036
02037
02038
02039
02040
02041
02042
02043 template <class T>
02044 WriteMultiData<T>::WriteMultiData(HDF5Handle& a_handle,
02045 const BoxLayout& a_layout,
02046 const int a_numIntv,
02047 const string& a_name,
02048 const int a_policyFlags,
02049 const IntVect& a_outputGhost,
02050 const bool a_newForm)
02051 :
02052 m_policyFlags(a_policyFlags),
02053 m_outputGhost(a_outputGhost),
02054 m_newForm(a_newForm)
02055 {
02056 CH_TIME("WriteMultiData::constructor");
02057 CH_assert(T::preAllocatable() == 0);
02058 #ifdef H516
02059 CH_assert(false);
02060 #endif
02061
02062
02063
02064 write(a_handle, a_layout);
02065
02066
02067
02068 T dummy;
02069 dataTypes<T>(m_types, dummy);
02070
02071 CH_assert(m_types.size() == 1);
02072
02073 m_allIntvFile.define(0, a_numIntv-1);
02074 getOffsets<T>(m_offsets,
02075 a_layout,
02076 m_types.size(),
02077 m_allIntvFile,
02078 m_outputGhost);
02079
02080 hsize_t flatdims[1];
02081
02082 {
02083 hid_t dataspace;
02084 flatdims[0] = m_offsets[0].back();
02085 if (m_newForm)
02086 {
02087 if (a_name == "M")
02088 {
02089 strcpy(m_dataname, "Mask");
02090 }
02091 else
02092 {
02093 sprintf(m_dataname, "%sRegular", a_name.c_str());
02094 }
02095 }
02096 else
02097 {
02098 sprintf(m_dataname, "%s:datatype=%i", a_name.c_str(), 0);
02099 }
02100 {
02101 CH_TIME("H5Screate");
02102 dataspace = H5Screate_simple(1, flatdims, NULL);
02103 }
02104 CH_assert(dataspace >=0);
02105 {
02106 CH_TIME("H5Dcreate");
02107 m_dataSet = H5Dcreate2(a_handle.groupID(), m_dataname,
02108 m_types[0],
02109 dataspace, H5P_DEFAULT,
02110 H5P_DEFAULT, H5P_DEFAULT);
02111 }
02112 CH_assert(m_dataSet >= 0);
02113 {
02114 CH_TIME("H5S_H5D_data_close");
02115 H5Sclose(dataspace);
02116 }
02117 }
02118
02119
02120 {
02121 hid_t offsetspace, offsetdataset;
02122 char offsetname[128];
02123 flatdims[0] = m_offsets[0].size();
02124 if (m_newForm)
02125 {
02126 if (a_name == "M")
02127 {
02128 strcpy(offsetname, "MaskOffsets");
02129 }
02130 else
02131 {
02132 sprintf(offsetname, "%sOffsets", a_name.c_str());
02133 }
02134 }
02135 else
02136 {
02137 sprintf(offsetname, "%s:offsets=%i", a_name.c_str(), 0);
02138 }
02139 {
02140 CH_TIME("H5S_H5D_offsets_create");
02141 offsetspace = H5Screate_simple(1, flatdims, NULL);
02142 CH_assert(offsetspace >= 0);
02143 offsetdataset = H5Dcreate2(a_handle.groupID(), offsetname,
02144 H5T_NATIVE_LLONG, offsetspace,
02145 H5P_DEFAULT,
02146 H5P_DEFAULT,H5P_DEFAULT);
02147 CH_assert(offsetdataset >= 0);
02148 }
02149 if (procID() == 0)
02150
02151 {
02152 herr_t err;
02153 CH_TIME("WriteOffsets");
02154 hid_t memdataspace = H5Screate_simple(1, flatdims, NULL);
02155 CH_assert(memdataspace >= 0);
02156 err = H5Dwrite(offsetdataset, H5T_NATIVE_LLONG, memdataspace,
02157 offsetspace, H5P_DEFAULT, &(m_offsets[0][0]));
02158 CH_assert(err >= 0);
02159 H5Sclose(memdataspace);
02160 }
02161 {
02162 CH_TIME("H5S_H5D_offsets_close");
02163 H5Sclose(offsetspace);
02164 H5Dclose(offsetdataset);
02165 }
02166 }
02167
02168
02169 {
02170 CH_TIME("WriteAttributes");
02171 HDF5HeaderData info;
02172 info.m_intvect["ghost"] = m_outputGhost;
02173 info.m_intvect["outputGhost"] = m_outputGhost;
02174 if (!m_newForm)
02175 {
02176 info.m_int["comps"] = m_allIntvFile.size();
02177 info.m_string["objectType"] = name(dummy);
02178 }
02179 std::string group = a_handle.getGroup();
02180 a_handle.setGroup(group+"/"+a_name+"_attributes");
02181 info.writeToFile(a_handle);
02182 a_handle.setGroup(group);
02183 }
02184
02185
02186 {
02187 m_maxBoxPerProc = a_layout.dataIterator().size();
02188 #ifdef CH_MPI
02189 long myCountBoxes = m_maxBoxPerProc;
02190 MPI_Allreduce(&myCountBoxes,
02191 &m_maxBoxPerProc,
02192 1, MPI_LONG, MPI_MAX, Chombo_MPI::comm);
02193 #endif
02194 }
02195 }
02196
02197
02198
02199
02200
02201
02202
02203
02204
02205
02206
02207
02208
02209 template <class T>
02210 int
02211 WriteMultiData<T>::writeData(const BoxLayoutData<T>& a_data,
02212 const Interval& a_intvMem,
02213 const Interval& a_intvFile)
02214 {
02215 CH_TIME("WriteMultiData::writeData");
02216 CH_assert(a_intvMem.size() > 0);
02217 CH_assert(a_intvMem.size() == a_intvFile.size());
02218 CH_assert(a_intvMem.begin() >= a_data.interval().begin() &&
02219 a_intvMem.end() <= a_data.interval().end());
02220 CH_assert(a_intvFile.begin() >= m_allIntvFile.begin() &&
02221 a_intvFile.end() <= m_allIntvFile.end());
02222
02223 int ret = 0;
02224
02225 hsize_t fileCount[1], memCount[SpaceDim+1];
02226 ch_offset_t fileOffset[1], memOffset[SpaceDim+1];
02227 herr_t err;
02228
02229
02230
02231 hid_t fileDataSpace = H5Dget_space(m_dataSet);
02232
02233
02234
02235 hid_t DXPL = H5Pcreate(H5P_DATASET_XFER);
02236 #ifdef CH_MPI
02237 if (m_policyFlags & CH_HDF5::IOPolicyCollectiveWrite)
02238 {
02239 H5Pset_dxpl_mpio(DXPL, H5FD_MPIO_COLLECTIVE);
02240 }
02241 #endif
02242
02243 if (m_policyFlags & CH_HDF5::IOPolicyMultiDimHyperslab)
02244
02245
02246
02247
02248 {
02249 CH_TIME("hyperslab_H5Dwrite");
02250 DataIterator dit = a_data.dataIterator();
02251 const void* buffer;
02252 for (int iBox = 0; iBox != m_maxBoxPerProc; ++iBox)
02253 {
02254 hid_t memDataSpace;
02255 if (dit.ok())
02256 {
02257 const T& data = a_data[dit];
02258 unsigned globalIdx = a_data.boxLayout().index(dit());
02259 Box dbox = data.box();
02260 Box rbox = a_data.box(dit());
02261 rbox.grow(m_outputGhost);
02262
02263 memCount[0] = a_data.nComp();
02264 D_TERM6(memCount[SpaceDim-0] = dbox.size(0);,
02265 memCount[SpaceDim-1] = dbox.size(1);,
02266 memCount[SpaceDim-2] = dbox.size(2);,
02267 memCount[SpaceDim-3] = dbox.size(3);,
02268 memCount[SpaceDim-4] = dbox.size(4);,
02269 memCount[SpaceDim-5] = dbox.size(5);)
02270 memDataSpace = H5Screate_simple(SpaceDim+1, memCount, NULL);
02271
02272 memCount[0] = a_intvMem.size();
02273 D_TERM6(memCount[SpaceDim-0] = rbox.size(0);,
02274 memCount[SpaceDim-1] = rbox.size(1);,
02275 memCount[SpaceDim-2] = rbox.size(2);,
02276 memCount[SpaceDim-3] = rbox.size(3);,
02277 memCount[SpaceDim-4] = rbox.size(4);,
02278 memCount[SpaceDim-5] = rbox.size(5);)
02279 memOffset[0] = a_intvMem.begin();
02280 D_TERM6(
02281 memOffset[SpaceDim-0] = rbox.smallEnd(0) - dbox.smallEnd(0);,
02282 memOffset[SpaceDim-1] = rbox.smallEnd(1) - dbox.smallEnd(1);,
02283 memOffset[SpaceDim-2] = rbox.smallEnd(2) - dbox.smallEnd(2);,
02284 memOffset[SpaceDim-3] = rbox.smallEnd(3) - dbox.smallEnd(3);,
02285 memOffset[SpaceDim-4] = rbox.smallEnd(4) - dbox.smallEnd(4);,
02286 memOffset[SpaceDim-5] = rbox.smallEnd(5) - dbox.smallEnd(5);)
02287 err = H5Sselect_hyperslab(memDataSpace, H5S_SELECT_SET,
02288 memOffset, NULL, memCount, NULL);
02289 CH_assert(err >= 0);
02290
02291 fileOffset[0] = m_offsets[0][globalIdx];
02292 fileCount[0] = m_offsets[0][globalIdx + 1] - fileOffset[0];
02293 if (fileCount[0] > 0)
02294 {
02295
02296 const hsize_t dpnts = rbox.numPts();
02297 fileOffset[0] += dpnts*a_intvFile.begin();
02298 fileCount[0] = dpnts*a_intvFile.size();
02299 CH_assert(fileOffset[0] + fileCount[0] <=
02300 m_offsets[0][globalIdx+1]);
02301 err = H5Sselect_hyperslab(fileDataSpace, H5S_SELECT_SET,
02302 fileOffset, NULL, fileCount, NULL);
02303 CH_assert(err >= 0);
02304 }
02305 else
02306 {
02307 H5Sselect_none(memDataSpace);
02308 H5Sselect_none(fileDataSpace);
02309 }
02310 buffer = data.dataPtr();
02311 ++dit;
02312 }
02313 else
02314 {
02315 std::memset(memCount, 0, (SpaceDim+1)*sizeof(hsize_t));
02316 memDataSpace = H5Screate_simple(SpaceDim+1, memCount, NULL);
02317 H5Sselect_none(memDataSpace);
02318 H5Sselect_none(fileDataSpace);
02319 buffer = 0;
02320 }
02321
02322 err = H5Dwrite(m_dataSet, m_types[0], memDataSpace, fileDataSpace,
02323 DXPL, buffer);
02324 CH_assert(err >= 0);
02325 H5Sclose(memDataSpace);
02326 if (err < 0)
02327 {
02328 ret = err;
02329 goto cleanup;
02330 }
02331 }
02332 }
02333 else
02334
02335
02336
02337
02338 {
02339 CH_TIME("linearize_H5Dwrite");
02340 void* linearBuffer;
02341
02342
02343 {
02344 CH_TIMELEAF("allocateLinearBuffer");
02345 long long bufferSize = 0;
02346 for (DataIterator dit = a_data.dataIterator(); dit.ok(); ++dit)
02347 {
02348 unsigned globalIdx = a_data.boxLayout().index(dit());
02349 {
02350 const long long allIntvSize =
02351 (m_offsets[0][globalIdx + 1] - m_offsets[0][globalIdx]);
02352 CH_assert(allIntvSize >= 0);
02353 if (allIntvSize > bufferSize)
02354 {
02355 bufferSize = allIntvSize;
02356 }
02357 }
02358 }
02359
02360 bufferSize = ((bufferSize/m_allIntvFile.size())*a_intvMem.size())*
02361 H5Tget_size(m_types[0]);
02362 linearBuffer = mallocMT(bufferSize);
02363 if (linearBuffer == NULL)
02364 {
02365 pout() << " bufferCapacity = " << (int)bufferSize << endl;
02366 MayDay::Error("Memory error in buffer allocation write");
02367 }
02368 }
02369
02370
02371 DataIterator dit = a_data.dataIterator();
02372 for (int iBox = 0; iBox != m_maxBoxPerProc; ++iBox)
02373 {
02374 hid_t memDataSpace;
02375 if (dit.ok())
02376 {
02377 const T& data = a_data[dit];
02378 unsigned globalIdx = a_data.boxLayout().index(dit());
02379 Box rbox = a_data.box(dit());
02380 rbox.grow(m_outputGhost);
02381 {
02382 CH_TIMELEAF("linearize");
02383 data.linearOut(linearBuffer, rbox, a_intvMem);
02384 }
02385 const hsize_t dpnts = rbox.numPts();
02386
02387 memCount[0] = dpnts*a_intvMem.size();
02388 memDataSpace = H5Screate_simple(1, memCount, NULL);
02389
02390 fileOffset[0] = m_offsets[0][globalIdx];
02391 fileCount[0] = m_offsets[0][globalIdx + 1] - fileOffset[0];
02392 if (fileCount[0] > 0)
02393 {
02394
02395 fileOffset[0] += dpnts*a_intvFile.begin();
02396 fileCount[0] = dpnts*a_intvFile.size();
02397 CH_assert(fileOffset[0] + fileCount[0] <=
02398 m_offsets[0][globalIdx+1]);
02399 err = H5Sselect_hyperslab(fileDataSpace, H5S_SELECT_SET,
02400 fileOffset, NULL, fileCount, NULL);
02401 CH_assert(err >= 0);
02402 }
02403 else
02404 {
02405 H5Sselect_none(memDataSpace);
02406 H5Sselect_none(fileDataSpace);
02407 }
02408 ++dit;
02409 }
02410 else
02411 {
02412 memCount[0] = 0;
02413 memDataSpace = H5Screate_simple(SpaceDim+1, memCount, NULL);
02414 H5Sselect_none(memDataSpace);
02415 H5Sselect_none(fileDataSpace);
02416 }
02417
02418 err = H5Dwrite(m_dataSet, m_types[0], memDataSpace, fileDataSpace,
02419 DXPL, linearBuffer);
02420 CH_assert(err >= 0);
02421 H5Sclose(memDataSpace);
02422 if (err < 0)
02423 {
02424 ret = err;
02425 freeMT(linearBuffer);
02426 goto cleanup;
02427 }
02428 }
02429 freeMT(linearBuffer);
02430 }
02431
02432 cleanup: ;
02433 H5Pclose(DXPL);
02434 H5Sclose(fileDataSpace);
02435
02436 return ret;
02437 }
02438
02439 #include "NamespaceFooter.H"
02440
02441 #else // CH_USE_HDF5
02442
02443
02444 #define HOFFSET(S,M) (offsetof(S,M))
02445
02446 #endif // CH_USE_HDF5 not defined
02447 #endif // CH_HDF5_H