Hi Dario, To add to what’s already been suggested, you can also look at the TransferableQuadraturePointData <https://dealii.org/current/doxygen/deal.II/classTransferableQuadraturePointData.html> and ContinuousQuadratureDataTransfer <https://dealii.org/current/doxygen/deal.II/classparallel_1_1distributed_1_1ContinuousQuadratureDataTransfer.html> classes (and the corresponding tests <https://github.com/dealii/dealii/search?q=ContinuousQuadratureDataTransfer>). These as mentioned in the CellDataStorage <https://dealii.org/current/doxygen/deal.II/classCellDataStorage.html> class documentation, and were implemented to support the transfer of QP data upon refinement / coarsening with p::d::SolutionTransfer. The TransferableQuadraturePointData::pack_values()/unpack_values() functions can also be used to help with the serialisation process.
As a side note, why there is no direct support for serialisation in the CellDataStorage is that the data stored there is very general and user-controlled, so its not straightforward to recreate polymorphic types upon restarting. For instance, there is no guarantee that you keep the data types stored in the CellDataStorage consistent at save and load time. So the idea is that the user is delegated responsibility for creating and initialising these classes, and then the unpack_values() function can be used to unroll some serialised data in to the fully initialised user-defined class instance. How you could to this is to create a save()/load() function (see how we use BOOST_SERIALIZATION_SPLIT_MEMBER() in the various classes that have this) and call the pack_values() and unpack_values() in save() and load() respectively. In save(), for instance, you could call pack_data() and unroll all of your data into the input std::vector that would have scope in the save() function. You would then serialise the std::vector. Loading would be the opposite — deserialise the vector, and pass it into the unpack_values() function that then extracts the data and initialises your class member variables with the loaded data. Something like that… The benefit of doing it this way is that you then get the interpolation-upon-refinement capabilities for the ContinuousQuadratureDataTransfer class for free! I see that we don’t actually have a test that shows how to do this, but if you’d be interested in adding one, that would be much appreciated. We could then augment the documentation with a working example to show how one can do it. Best, Jean-Paul > On 2. Jul 2021, at 18:19, 'Dario Abbondanza' via deal.II User Group > <dealii@googlegroups.com> wrote: > > Thank you for your reply. > > > I don't have a particularly good solution in my head, but if you look at how > > programs such as ASPECT do checkpoint/restart, then they use a mechanism > > such > > as p::d::SolutionTransfer to attach data to the triangulation, save/restore, > > and unattach. > > This is exactly what I am doing for the solution vectors. > > > > You might be able to achieve the same effect by using the > > CellDataTransfer class, though I have no idea if anyone has ever tried that. > > I will have a look at the CellDataTransfer class and see what I am able to > achieve. > Eventually I will reply to this discussion. > > Best, > Dario > > Il giorno venerdì 2 luglio 2021 alle 05:51:58 UTC+2 Wolfgang Bangerth ha > scritto: > On 7/1/21 4:48 AM, 'Dario Abbondanza' via deal.II User Group wrote: > > > > I'm studying a problem in which I need to store quantities (symmetric > > tensors > > and doubles) at quadrature points, and I do that by using a PointHistory > > class > > as done in step 44 > > (https://www.dealii.org/current/doxygen/deal.II/step_44.html#Quadraturepointhistory > > > > <https://www.dealii.org/current/doxygen/deal.II/step_44.html#Quadraturepointhistory>). > > > > The only difference is that I'm defining the CellDataStorage as > > > > CellDataStorage<typename > > parallel::distributed::Triangulation<dim>::cell_iterator, > > PointHistory<dim> > > > quadrature_point_history; > > > > because I use a parallel distributed triangulation. > > I don't have a particularly good solution in my head, but if you look at how > programs such as ASPECT do checkpoint/restart, then they use a mechanism such > as p::d::SolutionTransfer to attach data to the triangulation, save/restore, > and unattach. You might be able to achieve the same effect by using the > CellDataTransfer class, though I have no idea if anyone has ever tried that. > > Best > W. > > -- > ------------------------------------------------------------------------ > Wolfgang Bangerth email: bang...@colostate.edu > <applewebdata://6E6DEF90-6D58-42D0-B05A-47918F563960> > www: http://www.math.colostate.edu/~bangerth/ > <http://www.math.colostate.edu/~bangerth/> > > > ________________________________________________________ > Le informazioni contenute in questo messaggio di posta elettronica sono > strettamente riservate e indirizzate esclusivamente al destinatario. Si prega > di non leggere, fare copia, inoltrare a terzi o conservare tale messaggio se > non si è il legittimo destinatario dello stesso. Qualora tale messaggio sia > stato ricevuto per errore, si prega di restituirlo al mittente e di > cancellarlo permanentemente dal proprio computer. > The information contained in this e mail message is strictly confidential and > intended for the use of the addressee only. If you are not the intended > recipient, please do not read, copy, forward or store it on your computer. If > you have received the message in error, please forward it back to the sender > and delete it permanently from your computer system. > > > Fai crescere i nostri giovani ricercatori > dona il 5 per mille alla Sapienza > codice fiscale 80209930587 > > -- > The deal.II project is located at http://www.dealii.org/ > <http://www.dealii.org/> > For mailing list/forum options, see > https://groups.google.com/d/forum/dealii?hl=en > <https://groups.google.com/d/forum/dealii?hl=en> > --- > You received this message because you are subscribed to the Google Groups > "deal.II User Group" group. > To unsubscribe from this group and stop receiving emails from it, send an > email to dealii+unsubscr...@googlegroups.com > <mailto:dealii+unsubscr...@googlegroups.com>. > To view this discussion on the web visit > https://groups.google.com/d/msgid/dealii/d808015d-07ed-46d6-874b-4877958c222an%40googlegroups.com > > <https://groups.google.com/d/msgid/dealii/d808015d-07ed-46d6-874b-4877958c222an%40googlegroups.com?utm_medium=email&utm_source=footer>. -- The deal.II project is located at http://www.dealii.org/ For mailing list/forum options, see https://groups.google.com/d/forum/dealii?hl=en --- You received this message because you are subscribed to the Google Groups "deal.II User Group" group. To unsubscribe from this group and stop receiving emails from it, send an email to dealii+unsubscr...@googlegroups.com. To view this discussion on the web visit https://groups.google.com/d/msgid/dealii/7660A5BA-951B-4124-A75B-19D411A0A540%40gmail.com.