From: Wolfgang Bangerth Date: Sat, 3 Dec 2016 23:03:08 +0000 (-0700) Subject: Move DataOutInterface::write_pvd/visit_record() to namespace DataOutBase. X-Git-Tag: v8.5.0-rc1~343^2~2 X-Git-Url: https://gitweb.dealii.org/cgi-bin/gitweb.cgi?a=commitdiff_plain;h=09ab757dc49b2afe8b4bf70d23857a89bc3e1195;p=dealii.git Move DataOutInterface::write_pvd/visit_record() to namespace DataOutBase. These functions did not depend on the state of the DataOutInterface object and could therefore be made 'static'. On the other hand, we have traditionally kept such functions in namespace DataOutBase. Move them there, and deprecate the old versions of these functions. --- diff --git a/doc/news/changes.h b/doc/news/changes.h index 4aeeeed31f..401747de27 100644 --- a/doc/news/changes.h +++ b/doc/news/changes.h @@ -38,6 +38,18 @@ inconvenience this causes.

    +
  1. Deprecated: The DataOutInterface::write_pvd_record() and + DataOutInterface::write_visit_record() functions were actually + independent of any kind of data being written. As a consequence, + they did not depend on the state of the DataOutInterface object + to which they belonged (or any object of a derived class). Such + functions typically reside in the DataOutBase namespace instead + where they have now been moved. The functions in DataOutInterface + are now deprecated. +
    + (Wolfgang Bangerth, 2016/12/03) +
  2. +
  3. Changed: VectorTools::create_right_hand_side and VectorTools::create_boundary_right_hand_side now take an additional template parameter VectorType. diff --git a/include/deal.II/base/data_out_base.h b/include/deal.II/base/data_out_base.h index 2ca81c49f2..2796d5df92 100644 --- a/include/deal.II/base/data_out_base.h +++ b/include/deal.II/base/data_out_base.h @@ -1726,6 +1726,136 @@ namespace DataOutBase const VtkFlags &flags, std::ostream &out); + /** + * In ParaView it is possible to visualize time-dependent data tagged with + * the current integration time of a time dependent simulation. To use this + * feature you need a .pvd file that describes which VTU or + * PVTU file belongs to which timestep. This function writes a file that + * provides this mapping, i.e., it takes a list of pairs each of which + * indicates a particular time instant and the corresponding file that + * contains the graphical data for this time instant. + * + * A typical use case, in program that computes a time dependent solution, + * would be the following (time and time_step are + * member variables of the class with types double and + * unsigned int, respectively; the variable + * times_and_names is of type + * std::vector@ @>): + * + * @code + * template + * void MyEquation::output_results () const + * { + * DataOut data_out; + * + * data_out.attach_dof_handler (dof_handler); + * data_out.add_data_vector (solution, "U"); + * data_out.build_patches (); + * + * const std::string filename = "solution-" + + * Utilities::int_to_string (timestep_number, 3) + + * ".vtu"; + * std::ofstream output (filename.c_str()); + * data_out.write_vtu (output); + * + * times_and_names.push_back (std::pair (time, filename)); + * std::ofstream pvd_output ("solution.pvd"); + * DataOutBase::write_pvd_record (pvd_output, times_and_names); + * } + * @endcode + * + * @note See DataOutInterface::write_vtu or DataOutInterface::write_pvtu_record + * for writing solutions at each timestep. + * + * @note The second element of each pair, i.e., the file in which the + * graphical data for each time is stored, may itself be again a file that + * references other files. For example, it could be the name for a + * .pvtu file that references multiple parts of a parallel + * computation. + * + * @author Marco Engelhard, 2012 + */ + void write_pvd_record (std::ostream &out, + const std::vector > ×_and_names); + + /** + * This function is the exact equivalent of the write_pvtu_record() function + * but for older versions of the VisIt visualization program and for one + * visualization graph (or one time step only). See there for the purpose of + * this function. + * + * This function is documented in the "Creating a master file for parallel" + * section (section 5.7) of the "Getting data into VisIt" report that can be + * found here: + * https://wci.llnl.gov/codes/visit/2.0.0/GettingDataIntoVisIt2.0.0.pdf + */ + void write_visit_record (std::ostream &out, + const std::vector &piece_names); + + /** + * This function is equivalent to the write_visit_record() above but for + * multiple time steps. Here is an example of how the function would be + * used: + * @code + * const unsigned int number_of_time_steps = 3; + * std::vector > piece_names(number_of_time_steps); + * + * piece_names[0].push_back("subdomain_01.time_step_0.vtk"); + * piece_names[0].push_back("subdomain_02.time_step_0.vtk"); + * + * piece_names[1].push_back("subdomain_01.time_step_1.vtk"); + * piece_names[1].push_back("subdomain_02.time_step_1.vtk"); + * + * piece_names[2].push_back("subdomain_01.time_step_2.vtk"); + * piece_names[2].push_back("subdomain_02.time_step_2.vtk"); + * + * std::ofstream visit_output ("master_file.visit"); + * + * DataOutBase::write_visit_record(visit_output, piece_names); + * @endcode + * + * This function is documented in the "Creating a master file for parallel" + * section (section 5.7) of the "Getting data into VisIt" report that can be + * found here: + * https://wci.llnl.gov/codes/visit/2.0.0/GettingDataIntoVisIt2.0.0.pdf + */ + void write_visit_record (std::ostream &out, + const std::vector > &piece_names); + + /** + * This function is equivalent to the write_visit_record() above but for + * multiple time steps and with additional information about the time for + * each timestep. Here is an example of how the function would be + * used: + * @code + * const unsigned int number_of_time_steps = 3; + * std::vector > > times_and_piece_names(number_of_time_steps); + * + * times_and_piece_names[0].first = 0.0; + * times_and_piece_names[0].second.push_back("subdomain_01.time_step_0.vtk"); + * times_and_piece_names[0].second.push_back("subdomain_02.time_step_0.vtk"); + * + * times_and_piece_names[1].first = 0.5; + * times_and_piece_names[1].second.push_back("subdomain_01.time_step_1.vtk"); + * times_and_piece_names[1].second.push_back("subdomain_02.time_step_1.vtk"); + * + * times_and_piece_names[2].first = 1.0; + * times_and_piece_names[2].second.push_back("subdomain_01.time_step_2.vtk"); + * times_and_piece_names[2].second.push_back("subdomain_02.time_step_2.vtk"); + * + * std::ofstream visit_output ("master_file.visit"); + * + * DataOutBase::write_visit_record(visit_output, times_and_piece_names); + * @endcode + * + * This function is documented in the "Creating a master file for parallel" + * section (section 5.7) of the "Getting data into VisIt" report that can be + * found here: + * https://wci.llnl.gov/codes/visit/2.0.0/GettingDataIntoVisIt2.0.0.pdf + */ + void write_visit_record (std::ostream &out, + const std::vector > > ×_and_piece_names); + /** * Write the given list of patches to the output stream in SVG format. * @@ -2177,138 +2307,34 @@ public: const std::vector &piece_names) const; /** - * In ParaView it is possible to visualize time-dependent data tagged with - * the current integration time of a time dependent simulation. To use this - * feature you need a .pvd file that describes which VTU or - * PVTU file belongs to which timestep. This function writes a file that - * provides this mapping, i.e., it takes a list of pairs each of which - * indicates a particular time instant and the corresponding file that - * contains the graphical data for this time instant. - * - * A typical use case, in program that computes a time dependent solution, - * would be the following (time and time_step are - * member variables of the class with types double and - * unsigned int, respectively; the variable - * times_and_names is of type - * std::vector@ @>): - * - * @code - * template - * void MyEquation::output_results () const - * { - * DataOut data_out; - * - * data_out.attach_dof_handler (dof_handler); - * data_out.add_data_vector (solution, "U"); - * data_out.build_patches (); - * - * const std::string filename = "solution-" + - * Utilities::int_to_string (timestep_number, 3) + - * ".vtu"; - * std::ofstream output (filename.c_str()); - * data_out.write_vtu (output); - * - * times_and_names.push_back (std::pair (time, filename)); - * std::ofstream pvd_output ("solution.pvd"); - * data_out.write_pvd_record (pvd_output, times_and_names); - * } - * @endcode - * - * @note See DataOutBase::write_vtu or DataOutInterface::write_pvtu_record - * for writing solutions at each timestep. - * - * @note The second element of each pair, i.e., the file in which the - * graphical data for each time is stored, may itself be again a file that - * references other files. For example, it could be the name for a - * .pvtu file that references multiple parts of a parallel - * computation. - * - * @author Marco Engelhard, 2012 + * @deprecated Use DataOutBase::write_pvd_record() instead */ - void write_pvd_record (std::ostream &out, - const std::vector > ×_and_names) const; + static void write_pvd_record (std::ostream &out, + const std::vector > ×_and_names) DEAL_II_DEPRECATED; /** - * This function is the exact equivalent of the write_pvtu_record() function - * but for older versions of the VisIt visualization program and for one - * visualization graph (or one time step only). See there for the purpose of - * this function. - * - * This function is documented in the "Creating a master file for parallel" - * section (section 5.7) of the "Getting data into VisIt" report that can be - * found here: - * https://wci.llnl.gov/codes/visit/2.0.0/GettingDataIntoVisIt2.0.0.pdf + * @deprecated Use DataOutBase::write_visit_record() instead */ - void write_visit_record (std::ostream &out, - const std::vector &piece_names) const; + static + void + write_visit_record (std::ostream &out, + const std::vector &piece_names) DEAL_II_DEPRECATED; /** - * This function is equivalent to the write_visit_record() above but for - * multiple time steps. Here is an example of how the function would be - * used: - * @code - * DataOut data_out; - * - * const unsigned int number_of_time_steps = 3; - * std::vector > piece_names(number_of_time_steps); - * - * piece_names[0].push_back("subdomain_01.time_step_0.vtk"); - * piece_names[0].push_back("subdomain_02.time_step_0.vtk"); - * - * piece_names[1].push_back("subdomain_01.time_step_1.vtk"); - * piece_names[1].push_back("subdomain_02.time_step_1.vtk"); - * - * piece_names[2].push_back("subdomain_01.time_step_2.vtk"); - * piece_names[2].push_back("subdomain_02.time_step_2.vtk"); - * - * std::ofstream visit_output ("master_file.visit"); - * - * data_out.write_visit_record(visit_output, piece_names); - * @endcode - * - * This function is documented in the "Creating a master file for parallel" - * section (section 5.7) of the "Getting data into VisIt" report that can be - * found here: - * https://wci.llnl.gov/codes/visit/2.0.0/GettingDataIntoVisIt2.0.0.pdf + * @deprecated Use DataOutBase::write_visit_record() instead */ - void write_visit_record (std::ostream &out, - const std::vector > &piece_names) const; + static + void + write_visit_record (std::ostream &out, + const std::vector > &piece_names) DEAL_II_DEPRECATED; /** - * This function is equivalent to the write_visit_record() above but for - * multiple time steps and with additional information about the time for - * each timestep. Here is an example of how the function would be - * used: - * @code - * DataOut data_out; - * - * const unsigned int number_of_time_steps = 3; - * std::vector > > times_and_piece_names(number_of_time_steps); - * - * times_and_piece_names[0].first = 0.0; - * times_and_piece_names[0].second.push_back("subdomain_01.time_step_0.vtk"); - * times_and_piece_names[0].second.push_back("subdomain_02.time_step_0.vtk"); - * - * times_and_piece_names[1].first = 0.5; - * times_and_piece_names[1].second.push_back("subdomain_01.time_step_1.vtk"); - * times_and_piece_names[1].second.push_back("subdomain_02.time_step_1.vtk"); - * - * times_and_piece_names[2].first = 1.0; - * times_and_piece_names[2].second.push_back("subdomain_01.time_step_2.vtk"); - * times_and_piece_names[2].second.push_back("subdomain_02.time_step_2.vtk"); - * - * std::ofstream visit_output ("master_file.visit"); - * - * data_out.write_visit_record(visit_output, times_and_piece_names); - * @endcode - * - * This function is documented in the "Creating a master file for parallel" - * section (section 5.7) of the "Getting data into VisIt" report that can be - * found here: - * https://wci.llnl.gov/codes/visit/2.0.0/GettingDataIntoVisIt2.0.0.pdf + * @deprecated Use DataOutBase::write_visit_record() instead */ - void write_visit_record (std::ostream &out, - const std::vector > > ×_and_piece_names) const; + static + void + write_visit_record (std::ostream &out, + const std::vector > > ×_and_piece_names) DEAL_II_DEPRECATED; /** * Obtain data through get_patches() and write it to out in SVG diff --git a/source/base/data_out_base.cc b/source/base/data_out_base.cc index f137ec11c7..53ba05fc3c 100644 --- a/source/base/data_out_base.cc +++ b/source/base/data_out_base.cc @@ -5318,6 +5318,110 @@ namespace DataOutBase } + void + write_pvd_record (std::ostream &out, + const std::vector > ×_and_names) + { + AssertThrow (out, ExcIO()); + + out << "\n"; + + out << "\n"; + + out << "\n"; + out << " \n"; + + std::streamsize ss = out.precision(); + out.precision(12); + + for (unsigned int i=0; i\n"; + + out << " \n"; + out << "\n"; + + out.flush(); + out.precision(ss); + + AssertThrow (out, ExcIO()); + } + + + + void + write_visit_record (std::ostream &out, + const std::vector &piece_names) + { + out << "!NBLOCKS " << piece_names.size() << '\n'; + for (unsigned int i=0; i > &piece_names) + { + AssertThrow (out, ExcIO()); + + if (piece_names.size() == 0) + return; + + const double nblocks = piece_names[0].size(); + Assert(nblocks > 0, ExcMessage("piece_names should be a vector of nonempty vectors.") ) + + out << "!NBLOCKS " << nblocks << '\n'; + for (std::vector >::const_iterator domain = piece_names.begin(); domain != piece_names.end(); ++domain) + { + Assert(domain->size() == nblocks, ExcMessage("piece_names should be a vector of equal sized vectors.") ) + for (std::vector::const_iterator subdomain = domain->begin(); subdomain != domain->end(); ++subdomain) + out << *subdomain << '\n'; + } + + out << std::flush; + } + + + + void + write_visit_record (std::ostream &out, + const std::vector > > ×_and_piece_names) + { + AssertThrow (out, ExcIO()); + + if (times_and_piece_names.size() == 0) + return; + + const double nblocks = times_and_piece_names[0].second.size(); + Assert(nblocks > 0, ExcMessage("time_and_piece_names should contain nonempty vectors of filenames for every timestep.") ) + + for (std::vector > >::const_iterator domain = times_and_piece_names.begin(); + domain != times_and_piece_names.end(); ++domain) + out << "!TIME " << domain->first << '\n'; + + out << "!NBLOCKS " << nblocks << '\n'; + for (std::vector > >::const_iterator domain = times_and_piece_names.begin(); + domain != times_and_piece_names.end(); ++domain) + { + Assert(domain->second.size() == nblocks, ExcMessage("piece_names should be a vector of equal sized vectors.") ) + for (std::vector::const_iterator subdomain = domain->second.begin(); subdomain != domain->second.end(); ++subdomain) + out << *subdomain << '\n'; + } + + out << std::flush; + } + + + template void write_svg (const std::vector > &, const std::vector &, @@ -6223,43 +6327,6 @@ void DataOutInterface::write_vtu_in_parallel (const char *filename } -template -void -DataOutInterface:: -write_pvd_record (std::ostream &out, - const std::vector > ×_and_names) const -{ - AssertThrow (out, ExcIO()); - - out << "\n"; - - out << "\n"; - - out << "\n"; - out << " \n"; - - std::streamsize ss = out.precision(); - out.precision(12); - - for (unsigned int i=0; i\n"; - - out << " \n"; - out << "\n"; - - out.flush(); - out.precision(ss); - - AssertThrow (out, ExcIO()); -} - - template void DataOutInterface::write_pvtu_record (std::ostream &out, @@ -6360,40 +6427,20 @@ DataOutInterface::write_pvtu_record (std::ostream &out, template void -DataOutInterface::write_visit_record (std::ostream &out, - const std::vector &piece_names) const +DataOutInterface:: +write_pvd_record (std::ostream &out, + const std::vector > ×_and_names) { - out << "!NBLOCKS " << piece_names.size() << '\n'; - for (unsigned int i=0; i void DataOutInterface::write_visit_record (std::ostream &out, - const std::vector > &piece_names) const + const std::vector &piece_names) { - AssertThrow (out, ExcIO()); - - if (piece_names.size() == 0) - return; - - const double nblocks = piece_names[0].size(); - Assert(nblocks > 0, ExcMessage("piece_names should be a vector of nonempty vectors.") ) - - out << "!NBLOCKS " << nblocks << '\n'; - for (std::vector >::const_iterator domain = piece_names.begin(); domain != piece_names.end(); ++domain) - { - Assert(domain->size() == nblocks, ExcMessage("piece_names should be a vector of equal sized vectors.") ) - for (std::vector::const_iterator subdomain = domain->begin(); subdomain != domain->end(); ++subdomain) - out << *subdomain << '\n'; - } - - out << std::flush; + DataOutBase::write_visit_record (out, piece_names); } @@ -6401,30 +6448,19 @@ DataOutInterface::write_visit_record (std::ostream &out, template void DataOutInterface::write_visit_record (std::ostream &out, - const std::vector > > ×_and_piece_names) const + const std::vector > &piece_names) { - AssertThrow (out, ExcIO()); - - if (times_and_piece_names.size() == 0) - return; + DataOutBase::write_visit_record (out, piece_names); +} - const double nblocks = times_and_piece_names[0].second.size(); - Assert(nblocks > 0, ExcMessage("time_and_piece_names should contain nonempty vectors of filenames for every timestep.") ) - for (std::vector > >::const_iterator domain = times_and_piece_names.begin(); - domain != times_and_piece_names.end(); ++domain) - out << "!TIME " << domain->first << '\n'; - out << "!NBLOCKS " << nblocks << '\n'; - for (std::vector > >::const_iterator domain = times_and_piece_names.begin(); - domain != times_and_piece_names.end(); ++domain) - { - Assert(domain->second.size() == nblocks, ExcMessage("piece_names should be a vector of equal sized vectors.") ) - for (std::vector::const_iterator subdomain = domain->second.begin(); subdomain != domain->second.end(); ++subdomain) - out << *subdomain << '\n'; - } - - out << std::flush; +template +void +DataOutInterface::write_visit_record (std::ostream &out, + const std::vector > > ×_and_piece_names) +{ + DataOutBase::write_visit_record (out, times_and_piece_names); }