MFEM  v4.6.0
Finite element discretization library
fmsdatacollection.cpp
Go to the documentation of this file.
1 // Copyright (c) 2010-2023, Lawrence Livermore National Security, LLC. Produced
2 // at the Lawrence Livermore National Laboratory. All Rights reserved. See files
3 // LICENSE and NOTICE for details. LLNL-CODE-806117.
4 //
5 // This file is part of the MFEM library. For more information and source code
6 // availability visit https://mfem.org.
7 //
8 // MFEM is free software; you can redistribute it and/or modify it under the
9 // terms of the BSD-3 license. We welcome feedback and contributions, see file
10 // CONTRIBUTING.md for details.
11 
12 #include "../config/config.hpp"
13 
14 #ifdef MFEM_USE_FMS
15 
16 #include "fem.hpp"
17 #include "../general/text.hpp"
18 
19 #include <fmsio.h>
20 
21 #include <string>
22 #include <sstream>
23 
24 namespace mfem
25 {
26 
27 // class FMSDataCollection implementation
28 
29 FMSDataCollection::FMSDataCollection(const std::string& coll_name,
30  Mesh *mesh)
31  : DataCollection(coll_name, mesh),
32  fms_protocol("ascii")
33 {
34  appendRankToFileName = false; // always include rank in file names
35  cycle = 0; // always include cycle in directory names
36 }
37 
38 #ifdef MFEM_USE_MPI
39 
41  const std::string& coll_name,
42  Mesh *mesh)
43  : DataCollection(coll_name, mesh),
44  fms_protocol("ascii")
45 {
46  m_comm = comm;
47  MPI_Comm_rank(comm, &myid);
48  MPI_Comm_size(comm, &num_procs);
49  appendRankToFileName = true; // always include rank in file names
50  cycle = 0; // always include cycle in directory names
51 }
52 
53 #endif
54 
56 {
57  // empty
58 }
59 
61 {
62  // Convert this to FmsDataCollection.
63 
64  FmsDataCollection dc;
65  if (DataCollectionToFmsDataCollection(this, &dc) == 0)
66  {
67  std::string root(RootFileName());
68  int err = FmsIOWrite(root.c_str(), fms_protocol.c_str(), dc);
69  FmsDataCollectionDestroy(&dc);
70  if (err)
71  {
72  MFEM_ABORT("Error creating FMS file: " << root);
73  }
74  }
75  else
76  {
77  MFEM_ABORT("Error converting data collection");
78  }
79 }
80 
81 void FMSDataCollection::Load(int cycle)
82 {
83  DeleteAll();
84  this->cycle = cycle;
85 
86  FmsDataCollection dc;
87  std::string root(RootFileName());
88  int err = FmsIORead(root.c_str(), fms_protocol.c_str(), &dc);
89 
90  if (err == 0)
91  {
92  DataCollection *mdc = nullptr;
93  if (FmsDataCollectionToDataCollection(dc,&mdc) == 0)
94  {
95  // Tell the data collection we read that it does not own data.
96  // We will steal its data.
97  mdc->SetOwnData(false);
98 
99  SetCycle(mdc->GetCycle());
100  SetTime(mdc->GetTime());
101  SetTimeStep(mdc->GetTimeStep());
102  name = mdc->GetCollectionName();
103 
104  // Set mdc's mesh as our mesh.
105  SetMesh(mdc->GetMesh());
106 
107  // Set mdc's fields/qfields as ours.
108  std::vector<std::string> names;
109  for (const auto &pair : mdc->GetFieldMap())
110  {
111  names.push_back(pair.first);
112  RegisterField(pair.first, pair.second);
113  }
114  for (const auto &name : names)
115  {
116  mdc->DeregisterField(name);
117  }
118 
119  names.clear();
120  for (const auto &pair : mdc->GetQFieldMap())
121  {
122  names.push_back(pair.first);
123  RegisterQField(pair.first, pair.second);
124  }
125  for (const auto &name : names)
126  {
127  mdc->DeregisterField(name);
128  }
129 
130  // Indicate that we own the data.
131  SetOwnData(true);
132 
133  // Delete mdc. We stole its contents.
134  delete mdc;
135  }
136  FmsDataCollectionDestroy(&dc);
137  }
138  else
139  {
140  MFEM_ABORT("Error reading data collection: " << root);
141  }
142 }
143 
144 void FMSDataCollection::SetProtocol(const std::string &protocol)
145 {
146  fms_protocol = protocol;
147 }
148 
150 {
151  std::string res;
152  if (pad_digits_cycle)
153  {
154  res = prefix_path + name + "_" +
156  ".fms";
157  }
158  else
159  {
160  res = prefix_path + name + ".fms";
161  }
162  return res;
163 }
164 
165 } // namespace mfem
166 
167 #endif
void SetCycle(int c)
Set time cycle (for time-dependent simulations)
OutStream err(std::cerr)
Global stream used by the library for standard error output. Initially it uses the same std::streambu...
Definition: globals.hpp:71
bool appendRankToFileName
Append rank to any output file names.
virtual void Save()
Save the collection and a FMS blueprint root file.
virtual void RegisterQField(const std::string &q_field_name, QuadratureFunction *qf)
Add a QuadratureFunction to the collection.
void DeleteAll()
Delete data owned by the DataCollection including field information.
int pad_digits_cycle
Number of digits used for the cycle and MPI rank in filenames.
const QFieldMapType & GetQFieldMap() const
Get a const reference to the internal q-field map.
int GetCycle() const
Get time cycle (for time-dependent simulations)
MPI_Comm m_comm
Associated MPI communicator.
virtual void RegisterField(const std::string &field_name, GridFunction *gf)
Add a grid function to the collection.
int cycle
Time cycle; for time-dependent simulations cycle >= 0, otherwise = -1.
virtual void Load(int cycle=0)
Load the collection based blueprint data.
std::string to_padded_string(int i, int digits)
Convert an integer to a 0-padded string with the given number of digits.
Definition: text.hpp:54
double GetTimeStep() const
Get the simulation time step (for time-dependent simulations)
int DataCollectionToFmsDataCollection(DataCollection *mfem_dc, FmsDataCollection *dc)
int FmsDataCollectionToDataCollection(FmsDataCollection dc, DataCollection **mfem_dc)
void SetTime(double t)
Set physical time (for time-dependent simulations)
const FieldMapType & GetFieldMap() const
Get a const reference to the internal field map.
void SetOwnData(bool o)
Set the ownership of collection data.
int myid
MPI rank (in parallel)
std::string RootFileName()
Returns file name for the current cycle.
void SetTimeStep(double ts)
Set the simulation time step (for time-dependent simulations)
void SetProtocol(const std::string &protocol)
Set the FMS relay i/o protocol to use.
std::string prefix_path
A path where the directory with results is saved. If not empty, it has &#39;/&#39; at the end...
virtual void DeregisterField(const std::string &field_name)
Remove a grid function from the collection.
Mesh * GetMesh()
Get a pointer to the mesh in the collection.
FMSDataCollection(const std::string &collection_name, Mesh *mesh=NULL)
Constructor. The collection name is used when saving the data.
virtual void SetMesh(Mesh *new_mesh)
Set/change the mesh associated with the collection.
std::string name
Name of the collection, used as a directory name when saving.
int num_procs
Number of MPI ranks (in parallel)
const std::string & GetCollectionName() const
Get the name of the collection.
virtual ~FMSDataCollection()
We will delete the mesh and fields if we own them.
double GetTime() const
Get physical time (for time-dependent simulations)