G+Smo  25.01.0
Geometry + Simulation Modules
 
Loading...
Searching...
No Matches
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
22namespace gismo
23{
24
27{
28 enum method
29 {
30 HE = 0,
31 IHE = 1,
32 LE = 2,
33 ILE = 3,
34 TINE = 4,
36 BHE = 6,
37 IBHE = 7,
38 };
39};
40
43{
44 enum type
45 {
46 ossen = 0,
48 newton_next = 2
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,
73 };
74};
75
76
79{
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
103
110{
123};
124
127{
128 enum law
129 {
130 mixed_hooke = -2,
131 hooke = -1,
136 muscle = 4
140 };
141};
142
143struct 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{
173public:
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
217protected:
218 index_t m_width;
219};
220
221template <class T>
222std::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}
Simple progress bar class.
Definition gsBaseUtils.h:172
void display(double progress)
display the progress from 0 to 1
Definition gsBaseUtils.h:178
gsProgressBar(index_t width=25)
Constructor. Width is a number of symbols the progress bar spans.
Definition gsBaseUtils.h:175
void display(index_t progress, index_t total)
display the progress from 0 to 1
Definition gsBaseUtils.h:198
std::string to_string(const C &value)
Converts value to string, assuming "operator<<" defined on C.
Definition gsUtils.h:56
Provides structs and classes related to interfaces and boundaries.
Provides preprocessor directives configuration of G+Smo.
#define index_t
Definition gsConfig.h:32
This file contains the debugging and messaging system of G+Smo.
#define GISMO_ENSURE(cond, message)
Definition gsDebug.h:102
#define gsInfo
Definition gsDebug.h:43
Several utility functions for miscellaneous tasks.
The G+Smo namespace, containing all definitions for the library.
solver_status
Specifies the status of the iterative solver.
Definition gsBaseUtils.h:99
@ interrupted
method successfully converged
Definition gsBaseUtils.h:100
@ working
solver was interrupted after exceeding the limit of iterations
Definition gsBaseUtils.h:101
@ bad_solution
solver working
Definition gsBaseUtils.h:102
Specifies method used for mesh deformation in fluid-structure interaction.
Definition gsBaseUtils.h:27
method
Definition gsBaseUtils.h:29
@ IBHE
bi-harmonic extension
Definition gsBaseUtils.h:37
@ LE
incremental harmonic extension
Definition gsBaseUtils.h:32
@ TINE
incremental linear elasticity
Definition gsBaseUtils.h:34
@ TINE_StVK
tangential incremental nonlinear elasticity (with the neo-Hookean law)
Definition gsBaseUtils.h:35
@ BHE
tangential incremental nonlinear elasticity (with the St.Venant-Kirchhoff law)
Definition gsBaseUtils.h:36
@ ILE
linear elasticity
Definition gsBaseUtils.h:33
@ IHE
harmonic extension
Definition gsBaseUtils.h:31
@briefly Specifies iteration type for an iterative solver
Definition gsBaseUtils.h:90
type
Definition gsBaseUtils.h:92
@ next
each iteration yields an update
Definition gsBaseUtils.h:94
Specifies linear solver to use if it is hidden within some other class (like Newton's method or time ...
Definition gsBaseUtils.h:66
solver
Definition gsBaseUtils.h:68
@ BiCGSTABDiagonal
Conjugate gradient solver with diagonal (a.k.a. Jacobi) preconditioning: iterative(!...
Definition gsBaseUtils.h:72
@ LDLT
LU decomposition: direct, no matrix requirements, robust but a bit slow, Eigen and Pardiso available.
Definition gsBaseUtils.h:70
@ CGDiagonal
Cholesky decomposition pivoting: direct, simmetric positive or negative semidefinite,...
Definition gsBaseUtils.h:71
Specifies the material law to use.
Definition gsBaseUtils.h:127
law
Definition gsBaseUtils.h:129
@ mixed_neo_hooke_ln
S = lambda/2*(J^2-1)C^-1 + mu*(I-C^-1)
Definition gsBaseUtils.h:135
@ neo_hooke_quad
S = lambda*ln(J)*C^-1 + mu*(I-C^-1)
Definition gsBaseUtils.h:134
@ neo_hooke_ln
S = 2*mu*E + lambda*tr(E)*I.
Definition gsBaseUtils.h:133
@ saint_venant_kirchhoff
sigma = 2*mu*eps + lambda*tr(eps)*I
Definition gsBaseUtils.h:132
@ muscle
S = p*C^-1 + mu*(I-C^-1)
Definition gsBaseUtils.h:136
@ hooke
sigma = 2*mu*eps + p*I
Definition gsBaseUtils.h:131
Specifies the iteration type used to solve nonlinear systems.
Definition gsBaseUtils.h:43
type
Definition gsBaseUtils.h:45
@ newton_next
newton's method, 2nd order, yields updates to the solution
Definition gsBaseUtils.h:48
@ newton_update
stationary point iteration, 1st order, yields a new solution at each iteration
Definition gsBaseUtils.h:47
Specifies the verbosity of the iterative solver.
Definition gsBaseUtils.h:79
verbosity
Definition gsBaseUtils.h:81
@ some
no output
Definition gsBaseUtils.h:83
@ all
only essential output
Definition gsBaseUtils.h:84
method was interrupted because the current solution is invalid
Definition gsBaseUtils.h:110
components
Definition gsBaseUtils.h:112
@ normal_3D_vector
return all components of the 2D stress tensor as a 2x2 matrix
Definition gsBaseUtils.h:117
@ shear_3D_vector
Definition gsBaseUtils.h:119
@ all_3D_matrix
Definition gsBaseUtils.h:121
@ all_2D_vector
return von Mises stress as a scala
Definition gsBaseUtils.h:114
@ all_2D_matrix
Definition gsBaseUtils.h:116
Specifies the time integration scheme, see Wriggers, Nonlinear FEM, p. 205.
Definition gsBaseUtils.h:54
scheme
Definition gsBaseUtils.h:56
@ implicit_nonlinear
implicit scheme with linear problem (theta-scheme)
Definition gsBaseUtils.h:60
@ explicit_lumped
explicit scheme
Definition gsBaseUtils.h:58
@ implicit_linear
explicit scheme with lumped mass matrix
Definition gsBaseUtils.h:59