IFEM  90A354
SIMExplicitLMM.h
Go to the documentation of this file.
1 //==============================================================================
11 //==============================================================================
12 
13 #ifndef SIM_EXPLICIT_LMM_H_
14 #define SIM_EXPLICIT_LMM_H_
15 
16 #include "SystemMatrix.h"
17 #include "SIMenums.h"
18 #include "TimeIntUtils.h"
19 #include "TimeStep.h"
20 
21 class DataExporter;
22 
23 
24 namespace TimeIntegration {
25 
28  // and which derive from SIMbase.
29  template<class Solver>
31 {
32 public:
38  SIMExplicitLMM(Solver& solv, Method type, bool standalone = true,
39  const std::string& solField = "") :
40  solver(solv), alone(standalone), fieldName(solField)
41  {
42  if (type == AB2)
43  order = 2;
44  else if (type == AB3)
45  order = 3;
46  else if (type == AB4)
47  order = 4;
48  else if (type == AB5)
49  order = 5;
50  else
51  order = 1;
52 
53  loads.resize(order, nullptr);
54 
55  Solver::msgLevel = 1; // prints primary solution summary only
56  }
57 
60  {
61  for (SystemVector*& v : loads)
62  delete v;
63  }
64 
66  const ProcessAdm& getProcessAdm() const { return solver.getProcessAdm(); }
67 
69  bool solveStep(TimeStep& tp)
70  {
71  if (alone)
72  solver.getProcessAdm().cout <<"\n step = "<< tp.step <<" time = "<< tp.time.t << std::endl;
73 
74  TimeDomain time(tp.time);
75  time.t = tp.time.t - tp.time.dt;
76  if (!solver.assembleSystem(time, Vectors(1, solver.getSolution(1)),
77  !linear || (tp.step == 1)))
78  return false;
79 
80  loads[0] = solver.getRHSvector(0, true);
81 
82  const std::vector<std::vector<double>> AB_coefs =
83  {{1.0},
84  {-0.5, 1.5},
85  {5.0/12.0, -16.0/12.0, 23.0/12.0},
86  {-9.0/24.0, 37.0/24, -59.0/24.0, 55.0/24.0},
87  {251.0/720.0, -1274.0/720.0, 2616.0/720.0, -2774.0/720.0, 1901.0/720.0}};
88 
89  const int c_order = hasICs ? order-1 : std::min(order-1, tp.step-1);
90  const std::vector<double>& AB_coef = AB_coefs[c_order];
91  solver.getRHSvector(0, false)->mult(AB_coef.back() * tp.time.dt);
92  for (size_t j = 0; j < AB_coef.size()-1; ++j)
93  solver.addToRHSvector(0, *loads[c_order-j], AB_coef[j]*tp.time.dt);
94 
95  if (!solver.solveSystem(solver.getSolution()))
96  return false;
97 
98  if (linear)
99  solver.setMode(SIM::RHS_ONLY);
100 
101  solver.getSolution() += solver.getSolution(1);
102 
103  Vector dum;
104  solver.updateDirichlet(tp.time.t, &dum);
105  solver.applyDirichlet(solver.getSolution());
106 
107  if (alone)
108  solver.printSolutionSummary(solver.getSolution(), 0,
109  solver.getProblem()->getField1Name(1).c_str());
110 
111  return true;
112  }
113 
116  {
117  // Evaluate fluxes for initial conditions.
118  if (tp.step == 1) {
119  hasICs = true;
120  for (int j = 2; j <= order; ++j) {
121  std::stringstream str;
122  str << fieldName << j;
123  if (solver.hasIC(str.str())) {
124  TimeDomain time(tp.time);
125  time.t = tp.time.t - j*tp.time.dt;
126  if (!solver.assembleSystem(time, Vectors(1, solver.getSolution(j-1))))
127  return false;
128 
129  loads[j-2] = solver.getRHSvector(0, true);
130  } else {
131  hasICs = false;
132  break;
133  }
134  }
135  }
136 
137  delete loads.back();
138  for (int j = order-2; j >= 0; --j)
139  loads[j+1] = loads[j];
140 
141  return solver.advanceStep(tp);
142  }
143 
145  bool saveModel(char* fileName, int& geoBlk, int& nBlock)
146  {
147  return solver.saveModel(fileName, geoBlk, nBlock);
148  }
149 
151  bool saveStep(const TimeStep& tp, int& nBlock)
152  {
153  return solver.saveStep(tp, nBlock);
154  }
155 
157  void registerFields(DataExporter& exporter)
158  {
159  solver.registerFields(exporter);
160  }
161 
164  bool serialize(std::map<std::string,std::string>& data)
165  {
166  return solver.serialize(data);
167  }
168 
171  bool deSerialize(const std::map<std::string,std::string>& data)
172  {
173  return solver.deSerialize(data);
174  }
175 
177  void setLinear(bool enable) { linear = enable; }
178 
179 protected:
180  Solver& solver;
181  std::vector<SystemVector*> loads;
182  int order;
183  bool alone;
184  const std::string fieldName;
185  bool hasICs = false;
186  bool linear = false;
187 };
188 
189 }
190 
191 #endif
std::vector< Vector > Vectors
An array of real-valued vectors with algebraic operations.
Definition: MatVec.h:37
Various enums for simulation scope.
General representation of system matrices and vectors.
Various helpers for time integration.
Class for encapsulation of general time stepping parameters.
Administer and write data using DataWriters.
Definition: DataExporter.h:38
Class for administration of MPI processes in IFEM library.
Definition: ProcessAdm.h:33
Base class for representing a system vector on different formats.
Definition: SystemMatrix.h:32
Explicit linear multistep time stepping for SIM classes.
Definition: SIMExplicitLMM.h:31
bool serialize(std::map< std::string, std::string > &data)
Serialize internal state for restarting purposes.
Definition: SIMExplicitLMM.h:164
bool hasICs
If true, start with full order.
Definition: SIMExplicitLMM.h:185
const ProcessAdm & getProcessAdm() const
Returns the parallel process administrator.
Definition: SIMExplicitLMM.h:66
bool saveStep(const TimeStep &tp, int &nBlock)
Saves the converged results of a given time step to VTF file.
Definition: SIMExplicitLMM.h:151
bool linear
If true, mass matrix is constant.
Definition: SIMExplicitLMM.h:186
bool solveStep(TimeStep &tp)
Computes the solution for the current time step.
Definition: SIMExplicitLMM.h:69
Solver & solver
Reference to simulator.
Definition: SIMExplicitLMM.h:180
std::vector< SystemVector * > loads
Unscaled load vectors.
Definition: SIMExplicitLMM.h:181
bool advanceStep(TimeStep &tp)
Advances the time step one step forward.
Definition: SIMExplicitLMM.h:115
const std::string fieldName
Name of primary solution fields (for ICs)
Definition: SIMExplicitLMM.h:184
bool alone
If true, this is a standalone solver.
Definition: SIMExplicitLMM.h:183
bool deSerialize(const std::map< std::string, std::string > &data)
Set internal state from a serialized state.
Definition: SIMExplicitLMM.h:171
bool saveModel(char *fileName, int &geoBlk, int &nBlock)
Opens a new VTF-file and writes the model geometry to it.
Definition: SIMExplicitLMM.h:145
~SIMExplicitLMM()
Destructor frees up the load vectors.
Definition: SIMExplicitLMM.h:59
int order
Order of method.
Definition: SIMExplicitLMM.h:182
SIMExplicitLMM(Solver &solv, Method type, bool standalone=true, const std::string &solField="")
Constructor.
Definition: SIMExplicitLMM.h:38
void setLinear(bool enable)
Mark operator as linear to avoid repeated assembly and factorization.
Definition: SIMExplicitLMM.h:177
void registerFields(DataExporter &exporter)
Registers fields for output to a data exporter.
Definition: SIMExplicitLMM.h:157
Class for encapsulation of general time stepping parameters.
Definition: TimeStep.h:31
int step
Time step counter.
Definition: TimeStep.h:72
TimeDomain time
Time domain data.
Definition: TimeStep.h:74
A vector class with some added algebraic operations.
Definition: matrix.h:64
Utilities for time integration.
Definition: BDF.h:21
Method
Enum defining various solution methods.
Definition: TimeIntUtils.h:28
@ AB2
Second order Adams-Bashforth, explicit.
Definition: TimeIntUtils.h:37
@ AB3
Third order Adams-Bashforth, explicit.
Definition: TimeIntUtils.h:38
@ AB4
Fourth order Adams-Bashforth, explicit.
Definition: TimeIntUtils.h:39
@ AB5
Fifth order Adams-Bashforth, explicit.
Definition: TimeIntUtils.h:40
Struct representing the time domain.
Definition: TimeDomain.h:23
double dt
Current timestep (or load parameter) increment.
Definition: TimeDomain.h:25
double t
Current time (or pseudo time, load parameter)
Definition: TimeDomain.h:24