G+Smo  24.08.0
Geometry + Simulation Modules
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
gsBaseUtils.h
Go to the documentation of this file.
1 
15 #pragma once
16 
17 #include <gsCore/gsConfig.h>
18 #include <gsCore/gsDebug.h>
19 #include <gsUtils/gsUtils.h>
20 #include <gsCore/gsBoundary.h>
21 
22 namespace gismo
23 {
24 
26 struct ale_method
27 {
28  enum method
29  {
30  HE = 0,
31  IHE = 1,
32  LE = 2,
33  ILE = 3,
34  TINE = 4,
35  TINE_StVK = 5,
36  BHE = 6,
37  IBHE = 7,
38  };
39 };
40 
43 {
44  enum type
45  {
46  ossen = 0,
49  };
50 };
51 
54 {
55  enum scheme
56  {
57  explicit_ = 0,
61  };
62 };
63 
66 {
67  enum solver
68  {
69  LU = 0,
70  LDLT = 1,
71  CGDiagonal = 2,
73  };
74 };
75 
76 
79 {
80  enum verbosity
81  {
82  none = 0,
83  some = 1,
84  all = 2
85  };
86 };
87 
90 {
91  enum type
92  {
93  update = 0,
94  next = 1
95  };
96 };
97 
99 enum solver_status { converged,
103 
110 {
112  {
113  von_mises = 0,
115  all_2D_matrix = 2,
118  shear_3D_vector = 4,
120  all_3D_matrix = 5
122  };
123 };
124 
127 {
128  enum law
129  {
130  mixed_hooke = -2,
131  hooke = -1,
136  muscle = 4
137  };
141 };
142 
143 struct GISMO_EXPORT gsBoundaryInterface
144 {
145  gsBoundaryInterface() {}
146 
148  std::vector<patchSide> sidesA;
149  std::vector<patchSide> sidesB;
151  std::vector<std::pair<index_t,index_t> > patches;
152 
153 
154  void addInterfaceSide(index_t patchA, boundary::side sideA,
155  index_t patchB, boundary::side sideB)
156  {
157  sidesA.push_back(patchSide(patchA,sideA));
158  sidesB.push_back(patchSide(patchB,sideB));
159  }
160 
161  void addPatches(index_t patchA, index_t patchB)
162  { patches.push_back(std::pair<index_t,index_t>(patchA,patchB)); }
163 };
164 
172 {
173 public:
175  gsProgressBar(index_t width = 25) : m_width(width) {}
176 
178  void display(double progress)
179  {
180  GISMO_ENSURE(progress >= 0. && progress <= 1.,"Invalid progress value! Must be between 0 and 1.");
181  index_t threshold = (index_t)(progress*m_width);
182  gsInfo << "[";
183  for(index_t i = 0; i < m_width; i++)
184  if(i < threshold)
185  gsInfo << "=";
186  else if(i == threshold)
187  gsInfo << ">";
188  else
189  gsInfo << " ";
190  gsInfo << "] " << (abs(progress - 1.) < 1e-12 ? 100 : (index_t)(progress*100)) << " %\r";
191  gsInfo.flush();
192 
193  if (abs(progress - 1.) < 1e-12)
194  gsInfo << std::endl;
195  }
196 
198  void display(index_t progress, index_t total)
199  {
200  GISMO_ENSURE(progress >= 0 && progress <= total && total >= 0,"Invalid progress value!");
201  index_t threshold = (index_t)(1.*progress*m_width/total);
202  gsInfo << "[";
203  for(index_t i = 0; i < m_width; i++)
204  if(i < threshold)
205  gsInfo << "=";
206  else if(i == threshold)
207  gsInfo << ">";
208  else
209  gsInfo << " ";
210  gsInfo << "] " << progress << "/" << total << " \r";
211  gsInfo.flush();
212 
213  if (progress == total)
214  gsInfo << std::endl;
215  }
216 
217 protected:
218  index_t m_width;
219 };
220 
221 template <class T>
222 std::string secToHMS(T sec)
223 {
224  if (sec < 10)
225  return util::to_string(sec) + "s";
226 
227  index_t days = (index_t)(sec)/(3600*24);
228  index_t residual = (index_t)(sec)- 3600*24*days;
229  index_t hours = residual/3600;
230  residual -= 3600*hours;
231  index_t minutes = residual/60;
232  residual -= 60*minutes;
233 
234  std::string result = util::to_string(residual) + "s";
235  if (minutes > 0)
236  result = util::to_string(minutes) + "m" + result;
237  if (hours > 0)
238  result = util::to_string(hours) + "h" + result;
239  if (days > 0)
240  result = util::to_string(days) + "d" + result;
241 
242  return result;
243 }
244 
245 }
gsProgressBar(index_t width=25)
Constructor. Width is a number of symbols the progress bar spans.
Definition: gsBaseUtils.h:175
Specifies the verbosity of the iterative solver.
Definition: gsBaseUtils.h:78
method
Definition: gsBaseUtils.h:28
method was interrupted because the current solution is invalid
Definition: gsBaseUtils.h:109
Specifies the material law to use.
Definition: gsBaseUtils.h:126
Specifies method used for mesh deformation in fluid-structure interaction.
Definition: gsBaseUtils.h:26
implicit scheme with linear problem (theta-scheme)
Definition: gsBaseUtils.h:60
Provides structs and classes related to interfaces and boundaries.
type
Definition: gsBaseUtils.h:91
std::string to_string(const C &value)
Converts value to string, assuming &quot;operator&lt;&lt;&quot; defined on C.
Definition: gsUtils.h:56
S = lambda*ln(J)*C^-1 + mu*(I-C^-1)
Definition: gsBaseUtils.h:134
Specifies linear solver to use if it is hidden within some other class (like Newton&#39;s method or time ...
Definition: gsBaseUtils.h:65
side
Identifiers for topological sides.
Definition: gsBoundary.h:58
type
Definition: gsBaseUtils.h:44
solver_status
Specifies the status of the iterative solver.
Definition: gsBaseUtils.h:99
#define index_t
Definition: gsConfig.h:32
#define GISMO_ENSURE(cond, message)
Definition: gsDebug.h:102
verbosity
Definition: gsBaseUtils.h:80
This file contains the debugging and messaging system of G+Smo.
sigma = 2*mu*eps + p*I
Definition: gsBaseUtils.h:131
return von Mises stress as a scala
Definition: gsBaseUtils.h:114
each iteration yields an update
Definition: gsBaseUtils.h:94
sigma = 2*mu*eps + lambda*tr(eps)*I
Definition: gsBaseUtils.h:132
Provides preprocessor directives configuration of G+Smo.
Specifies iteration type for an iterative solver
Definition: gsBaseUtils.h:89
tangential incremental nonlinear elasticity (with the St.Venant-Kirchhoff law)
Definition: gsBaseUtils.h:36
Conjugate gradient solver with diagonal (a.k.a. Jacobi) preconditioning: iterative(!), simmetric, Eigen only.
Definition: gsBaseUtils.h:72
void display(index_t progress, index_t total)
display the progress from 0 to 1
Definition: gsBaseUtils.h:198
no output
Definition: gsBaseUtils.h:83
harmonic extension
Definition: gsBaseUtils.h:31
only essential output
Definition: gsBaseUtils.h:84
S = p*C^-1 + mu*(I-C^-1)
Definition: gsBaseUtils.h:136
#define gsInfo
Definition: gsDebug.h:43
method successfully converged
Definition: gsBaseUtils.h:100
S = lambda/2*(J^2-1)C^-1 + mu*(I-C^-1)
Definition: gsBaseUtils.h:135
incremental linear elasticity
Definition: gsBaseUtils.h:34
void display(double progress)
display the progress from 0 to 1
Definition: gsBaseUtils.h:178
law
Definition: gsBaseUtils.h:128
tangential incremental nonlinear elasticity (with the neo-Hookean law)
Definition: gsBaseUtils.h:35
explicit scheme
Definition: gsBaseUtils.h:58
incremental harmonic extension
Definition: gsBaseUtils.h:32
scheme
Definition: gsBaseUtils.h:55
explicit scheme with lumped mass matrix
Definition: gsBaseUtils.h:59
bi-harmonic extension
Definition: gsBaseUtils.h:37
Simple progress bar class.
Definition: gsBaseUtils.h:171
Several utility functions for miscellaneous tasks.
return all components of the 2D stress tensor as a 2x2 matrix
Definition: gsBaseUtils.h:117
EIGEN_STRONG_INLINE abs_expr< E > abs(const E &u)
Absolute value.
Definition: gsExpressions.h:4486
LU decomposition: direct, no matrix requirements, robust but a bit slow, Eigen and Pardiso available...
Definition: gsBaseUtils.h:70
Definition: gsBaseUtils.h:119
Specifies the time integration scheme, see Wriggers, Nonlinear FEM, p. 205.
Definition: gsBaseUtils.h:53
S = 2*mu*E + lambda*tr(E)*I.
Definition: gsBaseUtils.h:133
linear elasticity
Definition: gsBaseUtils.h:33
Cholesky decomposition pivoting: direct, simmetric positive or negative semidefinite, rather fast, Eigen and Pardiso available.
Definition: gsBaseUtils.h:71
solver working
Definition: gsBaseUtils.h:102
Definition: gsBaseUtils.h:116
Specifies the iteration type used to solve nonlinear systems.
Definition: gsBaseUtils.h:42
components
Definition: gsBaseUtils.h:111
solver was interrupted after exceeding the limit of iterations
Definition: gsBaseUtils.h:101
stationary point iteration, 1st order, yields a new solution at each iteration
Definition: gsBaseUtils.h:47
solver
Definition: gsBaseUtils.h:67
newton&#39;s method, 2nd order, yields updates to the solution
Definition: gsBaseUtils.h:48
Definition: gsBaseUtils.h:121