AlmaBTE  1.3
A solver of the space- and time-dependent Boltzmann transport equation for phonons
bulk_hdf5.hpp
Go to the documentation of this file.
1 // Copyright 2015-2018 The ALMA Project Developers
2 //
3 // Licensed under the Apache License, Version 2.0 (the "License");
4 // you may not use this file except in compliance with the License.
5 // You may obtain a copy of the License at
6 //
7 // http://www.apache.org/licenses/LICENSE-2.0
8 //
9 // Unless required by applicable law or agreed to in writing, software
10 // distributed under the License is distributed on an "AS IS" BASIS,
11 // WITHOUT WARRANTIES OR CONDITIONS OF ANY KIND, either express or
12 // implied. See the License for the specific language governing
13 // permissions and limitations under the License.
14 
15 #pragma once
16 
20 
21 #include <string>
22 #include <vector>
23 #include <tuple>
24 #include <boost/mpi.hpp>
25 #include <structures.hpp>
26 #include <symmetry.hpp>
27 #include <qpoint_grid.hpp>
28 #include <processes.hpp>
29 
30 namespace alma {
43 void save_bulk_hdf5(const char* filename,
44  const std::string& description,
45  const Crystal_structure& cell,
46  const Symmetry_operations& symmetries,
47  const Gamma_grid& grid,
48  const std::vector<Threeph_process>& processes,
49  const boost::mpi::communicator& comm);
50 
64 std::tuple<std::string,
65  std::unique_ptr<Crystal_structure>,
66  std::unique_ptr<Symmetry_operations>,
67  std::unique_ptr<Gamma_grid>,
68  std::unique_ptr<std::vector<Threeph_process>>>
69 load_bulk_hdf5(const char* filename, const boost::mpi::communicator& comm);
70 
81 std::vector<std::string> list_scattering_subgroups(
82  const char* filename,
83  const boost::mpi::communicator& comm);
84 
88 public:
90  const char* filename,
91  const std::string& groupname,
92  const boost::mpi::communicator& comm);
93 
95  const std::string name;
98  const bool preserves_symmetry;
100  const std::string description;
102  const Eigen::ArrayXXd w0;
105  const std::vector<std::string> attributes;
108  const std::vector<std::string> datasets;
110  const std::vector<std::string> groups;
113  Scattering_subgroup(const std::string& _name,
114  const bool _preserves_symmetry,
115  const std::string& _description,
116  const Eigen::Ref<const Eigen::ArrayXXd>& _w0)
117  : name(_name), preserves_symmetry(_preserves_symmetry),
118  description(_description), w0(_w0), attributes({}), datasets({}),
119  groups({}) {
120  }
121 
122 private:
124  Scattering_subgroup(const std::string& _name,
125  const bool _preserves_symmetry,
126  const std::string& _description,
127  const Eigen::Ref<const Eigen::ArrayXXd>& _w0,
128  const std::vector<std::string>& _attributes,
129  const std::vector<std::string>& _datasets,
130  const std::vector<std::string>& _groups)
131  : name(_name), preserves_symmetry(_preserves_symmetry),
132  description(_description), w0(_w0), attributes(_attributes),
133  datasets(_datasets), groups(_groups) {
134  }
135 };
146  const char* filename,
147  const std::string& groupname,
148  const boost::mpi::communicator& comm);
149 
157 void write_scattering_subgroup(const char* filename,
158  const Scattering_subgroup& group,
159  const boost::mpi::communicator& comm);
160 } // namespace alma
const std::vector< std::string > datasets
List with the names of all non-mandatory datasets in the subgroup.
Definition: bulk_hdf5.hpp:108
Definition: analytic1d.hpp:26
Scattering_subgroup(const std::string &_name, const bool _preserves_symmetry, const std::string &_description, const Eigen::Ref< const Eigen::ArrayXXd > &_w0)
Public constructor that does not allow extra attributes, datasets or groups.
Definition: bulk_hdf5.hpp:113
void save_bulk_hdf5(const char *filename, const std::string &description, const Crystal_structure &cell, const Symmetry_operations &symmetries, const Gamma_grid &grid, const std::vector< Threeph_process > &processes, const boost::mpi::communicator &comm)
Save all relevant information about a bulk material to an HDF5 file.
Definition: bulk_hdf5.cpp:125
const Eigen::ArrayXXd w0
Set of scattering rates [1 / ps].
Definition: bulk_hdf5.hpp:102
Classes and functions used to manipulate grids in reciprocal space.
std::vector< std::string > list_scattering_subgroups(const char *filename, const boost::mpi::communicator &comm)
Return the names of all subgroups in the "/scattering" group of an HDF5 file.
Definition: bulk_hdf5.cpp:785
const bool preserves_symmetry
True unless the scattering source alters the q-point equivalence classes.
Definition: bulk_hdf5.hpp:98
const std::string description
Free-format description.
Definition: bulk_hdf5.hpp:100
Detection and representation of allowed three-phonon processes.
Definitions of the basic data-handling classes in ALMA.
C++ interface to Atsushi Togo&#39;s spglib.
const std::string name
Name of the subgroup.
Definition: bulk_hdf5.hpp:95
void write_scattering_subgroup(const char *filename, const Scattering_subgroup &group, const boost::mpi::communicator &comm)
Add a subgroup to the "/scattering" group of the HDF5 file.
Definition: bulk_hdf5.cpp:938
const std::vector< std::string > groups
List with the names of all subsubgroups in the subgroup.
Definition: bulk_hdf5.hpp:110
POD class describing a subgroup of the "/scattering" group of an HDF5 file.
Definition: bulk_hdf5.hpp:87
const std::vector< std::string > attributes
List with the names of all non-mandatory attributes in the subgroup.
Definition: bulk_hdf5.hpp:105
friend Scattering_subgroup load_scattering_subgroup(const char *filename, const std::string &groupname, const boost::mpi::communicator &comm)
Read the scattering rates stored in a subgroup from the "/scattering" group.
Definition: bulk_hdf5.cpp:828
std::tuple< std::string, std::unique_ptr< Crystal_structure >, std::unique_ptr< Symmetry_operations >, std::unique_ptr< Gamma_grid >, std::unique_ptr< std::vector< Threeph_process > > > load_bulk_hdf5(const char *filename, const boost::mpi::communicator &comm)
Reconstruct all data structures from a file saved with save_bulk_hdf5().
Definition: bulk_hdf5.cpp:442