G+Smo  25.01.0
Geometry + Simulation Modules
 
Loading...
Searching...
No Matches
gsDynamicXBraid.h
1
17#pragma once
19#include <gsXBraid/gsXBraid.h>
20#include <gsIO/gsOptionList.h>
21
22namespace gismo
23{
24
32template <class T>
33class gsDynamicXBraid : public gsXBraid< gsMatrix<T> >
34{
35
36public:
37 typedef typename std::function<void(const index_t&,const T&,const gsVector<T>&,const gsVector<T>&,const gsVector<T>&)> callback_type;
38
39 virtual ~gsDynamicXBraid() {};
40
43 const gsDynamicBase<T> * solver,
44 const gsMpiComm& comm,
45 const T& tstart,
46 const T& tend,
47 const index_t numDofs,
48 index_t numSteps = 10)
49 :
50 gsXBraid< gsMatrix<T> >::gsXBraid(comm, tstart, tend, (int)numSteps),
51 m_solver(solver),
52 m_numDofs(numDofs)
53 {
54 defaultOptions();
55 m_callback = [](const index_t&,const T&,const gsVector<T>&,const gsVector<T>&,const gsVector<T>&){return;}; // set empty callback
56 }
57
58// Member functions
59public:
60 void defaultOptions()
61 {
62 // see gsXBraid/filedata/pde/heat2d_square_ibvp1.xml
63 m_options.addInt("CFactor", "Coarsening factor of the parallel-in-time multigrid solver", 2);
64 m_options.addInt("access", "Access level (never [=0], =after finished [=1(default)], each iteration [=2]", 1);
65 m_options.addInt("maxIter", "Maximum iteration numbers of the parallel-in-time multigrid solver", 100);
66 m_options.addInt("maxLevel", "Maximum numbers of parallel-in-time multigrid levels", 30);
67 m_options.addInt("minCLevel", "Minimum level of the parallel-in-time multigrid solver", 2);
68 m_options.addInt("norm", "Temporal norm of the parallel-in-time multigrid solver (1-norm [=1], 2-norm [=2(default)], inf-norm [=3])", 2);
69 m_options.addInt("numFMG", "Number of full multigrid steps of the parallel-in-time multigrid solver", 1);
70 m_options.addInt("numFMGVcyc", "Number of full multigrid V-cycles of the parallel-in-time multigrid solver", 1);
71 m_options.addInt("numMaxRef", "Maximum number of refinements of the parallel-in-time multigrid solver", 1);
72 m_options.addInt("numRelax", "Number of relaxation steps of the parallel-in-time multigrid solver", 1);
73 m_options.addInt("numStorage", "Number of storage of the parallel-in-time multigrid solver", -1);
74 m_options.addInt("print", "Print level (no output [=0], runtime inforation [=1], run statistics [=2(default)], debug [=3])", 2);
75 m_options.addReal("absTol", "Absolute tolerance of the parallel-in-time multigrid solver", 1e-6);
76 m_options.addReal("relTol", "Relative tolerance of the parallel-in-time multigrid solver", 1e-3);
77 m_options.addSwitch("fmg", "Perform full multigrid (default is off)", 0);
78 m_options.addSwitch("incrMaxLevels", "Increase the maximum number of parallel-in-time multigrid levels after performing a refinement (default is off)", 0);
79 m_options.addSwitch("periodic", "Periodic time grid (default is off)", 0);
80 m_options.addSwitch("refine", "Perform refinement in time (default off)", 0);
81 m_options.addSwitch("sequential", "Set the initial guess of the parallel-in-time multigrid solver as the sequential time stepping solution (default is off)", 0);
82 m_options.addSwitch("skip", "Skip all work on the first down cycle of the parallel-in-time multigrid solver (default on)", 1);
83 m_options.addSwitch("spatial", "Perform spatial coarsening and refinement (default is off)", 0);
84 m_options.addSwitch("tol", "Tolerance type (absolute [=true], relative [=false(default)]", 0);
85
86
87 m_options.addSwitch("extraVerbose", "Extra verbosity", 0);
88 }
89
90 void initialize()
91 {
92 // see gsXBraid/examples/xbraid_heatEquation_example.cpp
93 this->SetCFactor(m_options.getInt("CFactor"));
94 this->SetMaxIter(m_options.getInt("maxIter"));
95 this->SetMaxLevels(m_options.getInt("maxLevel"));
96 this->SetMaxRefinements(m_options.getInt("numMaxRef"));
97 this->SetMinCoarse(m_options.getInt("minCLevel"));
98 this->SetNFMG(m_options.getInt("numFMG"));
99 this->SetNFMGVcyc(m_options.getInt("numFMGVcyc"));
100 this->SetNRelax(m_options.getInt("numRelax"));
101 this->SetAccessLevel(m_options.getInt("access"));
102 this->SetPrintLevel(m_options.getInt("print"));
103 this->SetStorage(m_options.getInt("numStorage"));
104 this->SetTemporalNorm(m_options.getInt("norm"));
105
106 if (m_options.getSwitch("fmg")) this->SetFMG();
107 if (m_options.getSwitch("incrMaxLevels")) this->SetIncrMaxLevels();
108 if (m_options.getSwitch("periodic")) this->SetPeriodic(1); else this->SetPeriodic(0);
109 if (m_options.getSwitch("refine")) this->SetRefine(1); else this->SetRefine(0);
110 if (m_options.getSwitch("sequential")) this->SetSeqSoln(1); else this->SetSeqSoln(0);
111 if (m_options.getSwitch("skip")) this->SetSkip(1); else this->SetSkip(0);
112 if (m_options.getSwitch("spatial")) this->SetSpatialCoarsenAndRefine();
113 if (m_options.getSwitch("tol")) this->SetAbsTol(m_options.getReal("absTol"));
114 else this->SetRelTol(m_options.getReal("relTol"));
115 }
116
118 braid_Int Init(braid_Real t, braid_Vector *u_ptr) override
119 {
120 gsMatrix<T>* u = new gsMatrix<T>(3*m_numDofs, 1);
121
122 // Does this mean zero displacements?
123 u->setZero();
124
125 if (m_solver->solutionU().rows()==m_numDofs)
126 u->col(0).segment(0 ,m_numDofs) = m_solver->solutionU();
127 if (m_solver->solutionV().rows()==m_numDofs)
128 u->col(0).segment(m_numDofs ,m_numDofs) = m_solver->solutionV();
129 if (m_solver->solutionA().rows()==m_numDofs)
130 u->col(0).segment(2*m_numDofs,m_numDofs) = m_solver->solutionA();
131
132 *u_ptr = (braid_Vector) u;
133 return braid_Int(0);
134 }
135
136 gsOptionList & options() { return m_options; }
137 void setOptions(gsOptionList options) {m_options.update(options,gsOptionList::addIfUnknown); };
138
140 braid_Int Step(braid_Vector u, braid_Vector ustop, braid_Vector fstop, BraidStepStatus &status) override
141 {
142 gsVector<T>* u_ptr = (gsVector<T>*) u;
143 // gsMatrix<T>* ustop_ptr = (gsMatrix<T>*) ustop; // the guess is not used
144
145 // XBraid forcing
146 if (fstop != NULL)
147 {
148 gsVector<T>* fstop_ptr = (gsVector<T>*) fstop;
149 *u_ptr += *fstop_ptr;
150 }
151
152 // Get time step information
153 std::pair<braid_Real, braid_Real> time = static_cast<gsXBraidStepStatus&>(status).timeInterval();
154 if (m_options.getSwitch("extraVerbose")) gsInfo<<"Solving interval ["<<time.first<<" , "<<time.second<<"] (level "<<static_cast<gsXBraidStepStatus&>(status).level()<<")\n";
155 T t = time.first;
156 T dt = time.second - time.first;
157
158 // Solve time step
159 gsVector<T> U = (*u_ptr).segment(0 ,m_numDofs);
160 gsVector<T> V = (*u_ptr).segment(m_numDofs ,m_numDofs);
161 gsVector<T> A = (*u_ptr).segment(2*m_numDofs,m_numDofs);
162
163 gsStatus stepStatus = m_solver->step(t,dt,U,V,A);
164
165 u_ptr->segment(0 ,m_numDofs) = U;
166 u_ptr->segment(m_numDofs ,m_numDofs) = V;
167 u_ptr->segment(2*m_numDofs,m_numDofs) = A;
168
169 // Carry out adaptive refinement in time
170 if (static_cast<gsXBraidStepStatus&>(status).level() == 0)
171 {
172 if (stepStatus==gsStatus::Success)
173 {
174 braid_Real error = static_cast<gsXBraidStepStatus&>(status).error();
175 if (error != braid_Real(-1.0))
176 {
177 braid_Int rfactor = (braid_Int) std::ceil( std::sqrt( error / 1e-3) );
178 status.SetRFactor(rfactor);
179 }
180 else
181 status.SetRFactor(1);
182 }
183 // Refine if solution interval failed
184 else
185 {
186 if (m_options.getSwitch("extraVerbose")) gsInfo<<"Step "<<(static_cast<gsXBraidStepStatus&>(status)).timeIndex()<<" did not converge";
187 status.SetRFactor((braid_Int)2);
188 }
189 }
190 return braid_Int(0);
191 }
192
194 braid_Int SpatialNorm( braid_Vector u,
195 braid_Real *norm_ptr) override
196 {
197 gsVector<T>* u_ptr = (gsVector<T>*) u;
198 *norm_ptr = u_ptr->segment(0,m_numDofs).norm(); // Displacement-based norm
199 // *norm_ptr = u_ptr->norm();
200 return braid_Int(0);
201 }
202
204 braid_Int BufSize(braid_Int *size_ptr, BraidBufferStatus &status) override
205 {
206 *size_ptr = sizeof(T)*(m_numDofs*3+2); // +2 comes from rows, cols of the solution vector u.
207 return braid_Int(0);
208 }
209
210 void setCallback(callback_type callback) const {m_callback = callback;}
211
213 braid_Int Access(braid_Vector u, BraidAccessStatus &status) override
214 {
215 gsVector<T>* u_ptr = (gsVector<T>*) u;
216 m_callback((index_t) static_cast<gsXBraidAccessStatus&>(status).timeIndex(),
217 (T) static_cast<gsXBraidAccessStatus&>(status).time(),
218 (gsVector<T>)(*u_ptr).segment(0 ,m_numDofs),
219 (gsVector<T>)(*u_ptr).segment(m_numDofs ,m_numDofs),
220 (gsVector<T>)(*u_ptr).segment(2*m_numDofs,m_numDofs)
221 );
222 return braid_Int(0);
223 }
224
226 /*
227 NOTE: This routine is not implemented. How to do it:
228 1. Make the Coarsen and Refine routines virtual
229 2. Within the example, overload this class, and define Coarsen and Refine for the problem
230 3. Update the m_numDoFs within!!
231 */
232 braid_Int Coarsen(braid_Vector fu, braid_Vector *cu_ptr, BraidCoarsenRefStatus &status) override
233 {
234 gsMatrix<T> *fu_ptr = (gsMatrix<T>*) fu;
235 gsMatrix<T>* cu = new gsMatrix<T>();
236 *cu = *fu_ptr;
237 *cu_ptr = (braid_Vector) cu;
238 return braid_Int(0);
239 }
240
241 // Performs spatial refinement
242 braid_Int Refine(braid_Vector cu, braid_Vector *fu_ptr, BraidCoarsenRefStatus &status) override
243 {
244 gsMatrix<T> *cu_ptr = (gsMatrix<T>*) cu;
245 gsMatrix<T>* fu = new gsMatrix<T>();
246 *fu = *cu_ptr;
247 *fu_ptr = (braid_Vector) fu;
248 return braid_Int(0);
249 }
250
251
252// Class members
253protected:
254
255 const gsDynamicBase<T> * m_solver;
256 index_t m_numDofs;
257 gsOptionList m_options;
258 mutable callback_type m_callback;
259
260};
261
262} // namespace gismo
263
264#ifndef GISMO_BUILD_LIB
265#include GISMO_HPP_HEADER(gsDynamicImplicitEuler.hpp)
266#endif
Performs the arc length method to solve a nonlinear system of equations.
Definition gsDynamicBase.h:35
Performs the arc length method to solve a nonlinear system of equations.
Definition gsDynamicXBraid.h:34
braid_Int Step(braid_Vector u, braid_Vector ustop, braid_Vector fstop, BraidStepStatus &status) override
Performs a single step of the parallel-in-time multigrid.
Definition gsDynamicXBraid.h:140
braid_Int Coarsen(braid_Vector fu, braid_Vector *cu_ptr, BraidCoarsenRefStatus &status) override
Performs spatial coarsening.
Definition gsDynamicXBraid.h:232
braid_Int BufSize(braid_Int *size_ptr, BraidBufferStatus &status) override
Sets the size of the MPI communication buffer.
Definition gsDynamicXBraid.h:204
gsDynamicXBraid(const gsDynamicBase< T > *solver, const gsMpiComm &comm, const T &tstart, const T &tend, const index_t numDofs, index_t numSteps=10)
Constructor.
Definition gsDynamicXBraid.h:42
braid_Int SpatialNorm(braid_Vector u, braid_Real *norm_ptr) override
Computes the spatial norm of the given vector.
Definition gsDynamicXBraid.h:194
braid_Int Init(braid_Real t, braid_Vector *u_ptr) override
Initializes a vector.
Definition gsDynamicXBraid.h:118
braid_Int Access(braid_Vector u, BraidAccessStatus &status) override
Handles access for input/output.
Definition gsDynamicXBraid.h:213
A matrix with arbitrary coefficient type and fixed or dynamic size.
Definition gsMatrix.h:41
Class which holds a list of parameters/options, and provides easy access to them.
Definition gsOptionList.h:33
void addInt(const std::string &label, const std::string &desc, const index_t &value)
Adds a option named label, with description desc and value value.
Definition gsOptionList.cpp:201
void update(const gsOptionList &other, updateType type=ignoreIfUnknown)
Updates the object using the data from other.
Definition gsOptionList.cpp:253
bool getSwitch(const std::string &label) const
Reads value for option label from options.
Definition gsOptionList.cpp:51
const index_t & getInt(const std::string &label) const
Reads value for option label from options.
Definition gsOptionList.cpp:37
Real getReal(const std::string &label) const
Reads value for option label from options.
Definition gsOptionList.cpp:44
void addReal(const std::string &label, const std::string &desc, const Real &value)
Adds a option named label, with description desc and value value.
Definition gsOptionList.cpp:211
void addSwitch(const std::string &label, const std::string &desc, const bool &value)
Adds a option named label, with description desc and value value.
Definition gsOptionList.cpp:235
A serial communication class.
Definition gsMpiComm.h:289
A vector with arbitrary coefficient type and fixed or dynamic size.
Definition gsVector.h:37
Class defining the XBraid access status wrapper.
Definition gsXBraid.h:474
Class defining the XBraid step status wrapper.
Definition gsXBraid.h:641
braid_Int level()
Returns the current multigrid level.
Definition gsXBraid.h:651
Class defining the XBraid wrapper.
Definition gsXBraid.h:74
void SetStorage(braid_Int storage)
Definition gsXBraid.h:179
void SetSpatialCoarsenAndRefine()
Sets user-defined coarsening and refinement routine.
Definition gsXBraid.h:211
void SetMaxRefinements(braid_Int max_refinements)
Sets the max number of time grid refinement levels allowed.
Definition gsXBraid.h:185
void SetFMG()
Sets FMG (F-cycle)
Definition gsXBraid.h:167
void SetPeriodic(braid_Int periodic)
Sets periodic time grid (default is 0)
Definition gsXBraid.h:136
void SetMinCoarse(braid_Int min_coarse)
Definition gsXBraid.h:109
void SetMaxLevels(braid_Int max_levels)
Sets the maximum number of multigrid levels.
Definition gsXBraid.h:97
void SetAbsTol(braid_Real tol)
Sets absolute stopping tolerance.
Definition gsXBraid.h:121
void SetCFactor(braid_Int level, braid_Int cfactor)
Sets the coarsening factor cfactor on grid level (default is 2)
Definition gsXBraid.h:130
void SetNFMGVcyc(braid_Int nfmg_Vcyc)
Sets the number of V-cycles to do at each FMG level (default is 1)
Definition gsXBraid.h:173
void SetRefine(braid_Int refine)
Turns time refinement on (refine = 1) or off (refine = 0).
Definition gsXBraid.h:182
void SetSkip(braid_Int skip)
Sets whether to skip all work on the first down cycle (skip = 1). On by default.
Definition gsXBraid.h:103
void SetSeqSoln(braid_Int use_seq_soln)
Definition gsXBraid.h:155
void SetIncrMaxLevels()
Increases the max number of multigrid levels after performing a refinement.
Definition gsXBraid.h:100
void SetRelTol(braid_Real tol)
Sets relative stopping tolerance.
Definition gsXBraid.h:124
void SetMaxIter(braid_Int max_iter)
Sets max number of multigrid iterations.
Definition gsXBraid.h:139
void SetTemporalNorm(braid_Int tnorm)
Sets the temporal norm: 1-norm (1), 2-norm (2:default), inf-norm (3)
Definition gsXBraid.h:127
void SetPrintLevel(braid_Int print_level)
Definition gsXBraid.h:146
void SetNFMG(braid_Int k)
Sets the number of initial F-cycles to do before switching to V-cycles.
Definition gsXBraid.h:170
void SetNRelax(braid_Int level, braid_Int nrelax)
Definition gsXBraid.h:114
void SetAccessLevel(braid_Int access_level)
Definition gsXBraid.h:164
#define index_t
Definition gsConfig.h:32
#define gsInfo
Definition gsDebug.h:43
Base class to perform time integration of second-order structural dynamics systems.
Provides a list of labeled parameters/options that can be set and accessed easily.
Provides declarations of the XBraid wrapper.
The G+Smo namespace, containing all definitions for the library.
gsStatus
Definition gsStructuralAnalysisTypes.h:21
@ Success
Successful.