Mentions légales du service

Skip to content
Snippets Groups Projects
Commit 08d59b30 authored by Olivier COULAUD's avatar Olivier COULAUD
Browse files

tools compile with new writeHeader

parent 15038cb8
No related branches found
No related tags found
No related merge requests found
Pipeline #1107475 skipped
...@@ -682,6 +682,8 @@ namespace scalfmm::io ...@@ -682,6 +682,8 @@ namespace scalfmm::io
* Should be used if we write the particles with writeArrayOfReal method * Should be used if we write the particles with writeArrayOfReal method
* *
* @tparam PointType The type of a point * @tparam PointType The type of a point
* @tparam IntType1 The type of a nbParticles
* @tparam IntType1 The type of a theother integer variables
* @param centerOfBox The center of the Box (FPoint<FReal, dimension> class) * @param centerOfBox The center of the Box (FPoint<FReal, dimension> class)
* @param boxWidth The width of the box * @param boxWidth The width of the box
* @param nbParticles Number of particles in the box (or to save) * @param nbParticles Number of particles in the box (or to save)
...@@ -690,10 +692,10 @@ namespace scalfmm::io ...@@ -690,10 +692,10 @@ namespace scalfmm::io
* @param dimension. The dimension of the problem * @param dimension. The dimension of the problem
* @param nb_input_values. The number of input values in a particle * @param nb_input_values. The number of input values in a particle
*/ */
template<class PointType, typename IntType> template<class PointType, typename IntType1, typename IntType2>
auto writeHeader(const PointType& centerOfBox, const FReal& boxWidth, const IntType& nbParticles, auto writeHeader(const PointType& centerOfBox, const FReal& boxWidth, const IntType1& nbParticles,
const IntType& dataType, const IntType& nbDataPerRecord, const IntType& dimension, const IntType2& dataType, const IntType2& nbDataPerRecord, const IntType2& dimension,
const IntType& nb_input_values) -> void const IntType2& nb_input_values) -> void
{ {
std::array<unsigned int, 4> typeFReal = { std::array<unsigned int, 4> typeFReal = {
static_cast<unsigned int>(dataType), static_cast<unsigned int>(nbDataPerRecord), static_cast<unsigned int>(dataType), static_cast<unsigned int>(nbDataPerRecord),
...@@ -723,11 +725,10 @@ namespace scalfmm::io ...@@ -723,11 +725,10 @@ namespace scalfmm::io
} }
std::cout << std::endl << std::flush; std::cout << std::endl << std::flush;
} }
template<typename IntType1, typename IntType2>
template<typename IntType> auto writeHeader(const FReal* centerOfBox, const FReal& boxWidth, const IntType1& nbParticles,
auto writeHeader(const FReal* centerOfBox, const FReal& boxWidth, const IntType& nbParticles, const IntType2& dataType, const IntType2& nbDataPerRecord, const IntType2& dimension,
const IntType& dataType, const IntType& nbDataPerRecord, const IntType& dimension, const IntType2& nb_input_values) -> void
const IntType& nb_input_values) -> void
{ {
std::array<unsigned int, 4> typeFReal = { std::array<unsigned int, 4> typeFReal = {
static_cast<unsigned int>(dataType), static_cast<unsigned int>(nbDataPerRecord), static_cast<unsigned int>(dataType), static_cast<unsigned int>(nbDataPerRecord),
...@@ -895,8 +896,8 @@ namespace scalfmm::io ...@@ -895,8 +896,8 @@ namespace scalfmm::io
// //
// write the particles // write the particles
const auto& centre = tree.box_center(); const auto& centre = tree.box_center();
this->writeHeader(centre, tree.box_width(), number_particles, sizeof(value_type), nb_elt_per_par, this->writeHeader(centre, tree.box_width(), number_particles, static_cast<std::size_t>(sizeof(value_type)),
centre.dimension, nb_input_elements); nb_elt_per_par, centre.dimension, nb_input_elements);
this->writeArrayOfReal(particles.data()->data(), nb_elt_per_par, number_particles); this->writeArrayOfReal(particles.data()->data(), nb_elt_per_par, number_particles);
} }
...@@ -921,9 +922,11 @@ namespace scalfmm::io ...@@ -921,9 +922,11 @@ namespace scalfmm::io
{ {
// get the number of elements per particles in the container build with tuples. // get the number of elements per particles in the container build with tuples.
using particle_type = typename ContainerType::value_type; using particle_type = typename ContainerType::value_type;
constexpr int nb_elt_per_par = meta::tuple_size_v<particle_type>; constexpr std::size_t nb_elt_per_par = meta::tuple_size_v<particle_type>;
static constexpr std::size_t dimension = particle_type::dimension_size;
// Not good output_values are put in input_values // Not good output_values are put in input_values
constexpr int nb_input_per_par = particle_type::inputs_size; constexpr std::size_t nb_input_per_par = particle_type::inputs_size;
/// @todo check for different input and output types (double versus complexe) /// @todo check for different input and output types (double versus complexe)
using data_type = typename particle_type::outputs_value_type; using data_type = typename particle_type::outputs_value_type;
// //
...@@ -943,8 +946,8 @@ namespace scalfmm::io ...@@ -943,8 +946,8 @@ namespace scalfmm::io
}; };
// write the particles // write the particles
// Here we need to separate input from output variables - no tools yet // Here we need to separate input from output variables - no tools yet
this->writeHeader(center, box_width, number_particles, sizeof(data_type), nb_elt_per_par, center.dimension, this->writeHeader(center, box_width, number_particles, static_cast<std::size_t>(sizeof(data_type)),
nb_input_per_par); nb_elt_per_par, dimension, nb_input_per_par);
this->writeArrayOfReal(particles.data()->data(), nb_elt_per_par, number_particles); this->writeArrayOfReal(particles.data()->data(), nb_elt_per_par, number_particles);
} }
......
...@@ -538,15 +538,15 @@ auto main([[maybe_unused]] int argc, [[maybe_unused]] char* argv[]) -> int ...@@ -538,15 +538,15 @@ auto main([[maybe_unused]] int argc, [[maybe_unused]] char* argv[]) -> int
if(use_double) if(use_double)
{ {
scalfmm::io::FFmaGenericWriter<value_type> writer(output_file); scalfmm::io::FFmaGenericWriter<value_type> writer(output_file);
writer.writeHeader(boxCenter, boxWith, nb_particles, sizeof(value_type), n_data_per_particle, dimension, writer.writeHeader(boxCenter, boxWith, nb_particles, static_cast<int>(sizeof(value_type)), n_data_per_particle,
number_input_values); dimension, number_input_values);
writer.writeArrayOfReal(data.data(), n_data_per_particle, nb_particles); writer.writeArrayOfReal(data.data(), n_data_per_particle, nb_particles);
} }
else else
{ {
scalfmm::io::FFmaGenericWriter<float> writer(output_file); scalfmm::io::FFmaGenericWriter<float> writer(output_file);
writer.writeHeader(boxCenter, boxWith, nb_particles, sizeof(float), n_data_per_particle, dimension, writer.writeHeader(boxCenter, boxWith, nb_particles, static_cast<int>(sizeof(float)), n_data_per_particle,
number_input_values); dimension, number_input_values);
std::vector<float> dataFloat(data.size()); std::vector<float> dataFloat(data.size());
std::copy(data.begin(), data.end(), dataFloat.begin()); std::copy(data.begin(), data.end(), dataFloat.begin());
writer.writeArrayOfReal(dataFloat.data(), n_data_per_particle, nb_particles); writer.writeArrayOfReal(dataFloat.data(), n_data_per_particle, nb_particles);
......
...@@ -90,7 +90,7 @@ namespace loc_args ...@@ -90,7 +90,7 @@ namespace loc_args
}; };
} // namespace loc_args } // namespace loc_args
template<int Dimension, typename value_type> template<unsigned int Dimension, typename value_type>
auto run(const std::string& input_file, const std::string& output_file, const bool verbose) -> int auto run(const std::string& input_file, const std::string& output_file, const bool verbose) -> int
{ {
...@@ -99,7 +99,7 @@ auto run(const std::string& input_file, const std::string& output_file, const bo ...@@ -99,7 +99,7 @@ auto run(const std::string& input_file, const std::string& output_file, const bo
scalfmm::io::FFmaGenericLoader<value_type, Dimension> loader(input_file, verbose); scalfmm::io::FFmaGenericLoader<value_type, Dimension> loader(input_file, verbose);
// //
auto nb_particles = loader.getNumberOfParticles(); std::size_t nb_particles = loader.getNumberOfParticles();
const unsigned int n_data_per_particle = loader.getNbRecordPerline(); const unsigned int n_data_per_particle = loader.getNbRecordPerline();
if(Dimension != loader.get_dimension()) if(Dimension != loader.get_dimension())
...@@ -127,8 +127,9 @@ auto run(const std::string& input_file, const std::string& output_file, const bo ...@@ -127,8 +127,9 @@ auto run(const std::string& input_file, const std::string& output_file, const bo
// } // }
scalfmm::io::FFmaGenericWriter<value_type> writer(output_file); scalfmm::io::FFmaGenericWriter<value_type> writer(output_file);
writer.writeHeader(loader.getCenterOfBox(), loader.getBoxWidth(), nb_particles, sizeof(value_type), writer.writeHeader(loader.getCenterOfBox(), loader.getBoxWidth(), nb_particles,
n_data_per_particle, Dimension, loader.get_number_of_input_per_record()); static_cast<unsigned int>(sizeof(value_type)), n_data_per_particle, Dimension,
loader.get_number_of_input_per_record());
writer.writeArrayOfReal(particles1.data(), n_data_per_particle, nb_particles); writer.writeArrayOfReal(particles1.data(), n_data_per_particle, nb_particles);
std::cout << "End of writing" << std::endl; std::cout << "End of writing" << std::endl;
......
0% Loading or .
You are about to add 0 people to the discussion. Proceed with caution.
Finish editing this message first!
Please register or to comment