void
restart();
+ /**
+ * Determine an estimate for the memory consumption (in bytes) of this
+ * object.
+ */
+ std::size_t
+ memory_consumption() const;
+
+ /**
+ * Write or read the data of this object to or from a stream for the purpose
+ * of serialization using the [BOOST serialization
+ * library](https://www.boost.org/doc/libs/1_74_0/libs/serialization/doc/index.html).
+ */
+ template <class Archive>
+ void
+ serialize(Archive &ar, const unsigned int version);
+
private:
/**
* The beginning of the time interval.
}
+
+template <class Archive>
+inline void
+DiscreteTime::serialize(Archive &ar, const unsigned int)
+{
+ ar &start_time;
+ ar &end_time;
+ ar ¤t_time;
+ ar &next_time;
+ ar &previous_time;
+ ar &start_step_size;
+ ar &step_number;
+}
+
DEAL_II_NAMESPACE_CLOSE
#endif
#include <deal.II/base/discrete_time.h>
#include <deal.II/base/exceptions.h>
+#include <deal.II/base/memory_consumption.h>
DEAL_II_NAMESPACE_OPEN
step_number = 0;
}
+
+
+std::size_t
+DiscreteTime::memory_consumption() const
+{
+ return (MemoryConsumption::memory_consumption(start_time) +
+ MemoryConsumption::memory_consumption(end_time) +
+ MemoryConsumption::memory_consumption(current_time) +
+ MemoryConsumption::memory_consumption(next_time) +
+ MemoryConsumption::memory_consumption(previous_time) +
+ MemoryConsumption::memory_consumption(start_step_size) +
+ MemoryConsumption::memory_consumption(step_number));
+}
+
DEAL_II_NAMESPACE_CLOSE
--- /dev/null
+// ------------------------------------------------------------------------
+//
+// SPDX-License-Identifier: LGPL-2.1-or-later
+// Copyright (C) 2020 - 2022 by the deal.II authors
+//
+// This file is part of the deal.II library.
+//
+// Part of the source code is dual licensed under Apache-2.0 WITH
+// LLVM-exception OR LGPL-2.1-or-later. Detailed license information
+// governing the source code and code contributions can be found in
+// LICENSE.md and CONTRIBUTING.md at the top level directory of deal.II.
+//
+// ------------------------------------------------------------------------
+
+
+// check and illustrate the serialization process for DiscreteTime class
+
+#include <deal.II/base/discrete_time.h>
+
+#include "serialization.h"
+
+
+
+void
+test()
+{
+ // save data to archive
+ std::ostringstream oss;
+ {
+ DiscreteTime time(0.0, 1.0, 0.6);
+ time.advance_time();
+
+ deallog << "Before serialization start time = " << time.get_start_time()
+ << ", end time = " << time.get_end_time()
+ << ", current time = " << time.get_current_time()
+ << ", next time = " << time.get_next_time()
+ << ", previous time = " << time.get_previous_time()
+ << ", next step size = " << time.get_next_step_size()
+ << ", previous step size = " << time.get_previous_step_size()
+ << ", step number = " << time.get_step_number() << std::endl;
+
+ boost::archive::text_oarchive oa(oss, boost::archive::no_header);
+ oa << time;
+
+ // archive and stream closed when
+ // destructors are called
+ }
+
+ // verify correctness of the serialization.
+ {
+ std::istringstream iss(oss.str());
+ boost::archive::text_iarchive ia(iss, boost::archive::no_header);
+
+ DiscreteTime time(1.0, 2.0, 0.5); // Initialize with garbage values.
+ ia >> time;
+
+ deallog << "After serialization start time = " << time.get_start_time()
+ << ", end time = " << time.get_end_time()
+ << ", current time = " << time.get_current_time()
+ << ", next time = " << time.get_next_time()
+ << ", previous time = " << time.get_previous_time()
+ << ", next step size = " << time.get_next_step_size()
+ << ", previous step size = " << time.get_previous_step_size()
+ << ", step number = " << time.get_step_number() << std::endl;
+ }
+
+ deallog << "OK" << std::endl << std::endl;
+}
+
+
+int
+main(int argc, char *argv[])
+{
+ initlog();
+ deallog << std::setprecision(4);
+
+ test();
+}