67 : lambda(lambda_),
mu(mu_), u(NULL), si(0), sj(0) { }
70 void SetComponent(
int i,
int j) { si = i; sj = j; }
76 class VisMan :
public iostream
84 int win_x, win_y, win_w, win_h;
85 int win_stride_x, win_stride_y, win_nx;
90 void CloseConnection();
91 void PositionWindow();
99 ostream &
operator<<(ostream &v,
void (*
f)(VisMan&));
101 int main(
int argc,
char *argv[])
104 const char *mesh_file =
"../data/beam-tri.mesh";
109 bool visualization = 1;
112 args.
AddOption(&mesh_file,
"-m",
"--mesh",
113 "Mesh file to use.");
114 args.
AddOption(&ref_levels,
"-r",
"--refine",
115 "Number of times to refine the mesh uniformly, -1 for auto.");
117 "Finite element order (polynomial degree).");
119 "One of the two DG penalty parameters, typically +1/-1."
120 " See the documentation of class DGElasticityIntegrator.");
122 "One of the two DG penalty parameters, should be positive."
123 " Negative values are replaced with (order+1)^2.");
124 args.
AddOption(&visualization,
"-vis",
"--visualization",
"-no-vis",
125 "--no-visualization",
126 "Enable or disable GLVis visualization.");
135 kappa = (order+1)*(order+1);
140 Mesh mesh(mesh_file, 1, 1);
145 cerr <<
"\nInput mesh should have at least two materials and "
146 <<
"two boundary attributes! (See schematic in ex17.cpp)\n"
154 ref_levels = (int)floor(log(5000./mesh.
GetNE())/log(2.)/
dim);
156 for (
int l = 0; l < ref_levels; l++)
170 cout <<
"Number of finite element unknowns: " << fespace.
GetTrueVSize()
171 <<
"\nAssembling: " << flush;
211 cout <<
"r.h.s. ... " << flush;
214 init_x, lambda_c, mu_c, alpha, kappa), dir_bdr);
232 cout <<
"matrix ... " << flush;
238 cout <<
"done." << endl;
243 #ifndef MFEM_USE_SUITESPARSE
248 const double rtol = 1e-6;
251 PCG(A, M, B, X, 3, 5000, rtol*rtol, 0.0);
255 GMRES(A, M, B, X, 3, 5000, 100, rtol*rtol, 0.0);
260 umf_solver.
Control[UMFPACK_ORDERING] = UMFPACK_ORDERING_METIS;
262 umf_solver.
Mult(B, X);
273 if (visualization) { reference_nodes = *mesh.
GetNodes(); }
281 ofstream mesh_ofs(
"displaced.mesh");
282 mesh_ofs.precision(8);
283 mesh.
Print(mesh_ofs);
284 ofstream sol_ofs(
"sol.gf");
285 sol_ofs.precision(8);
294 VisMan vis(vishost, visport);
295 const char *glvis_keys = (dim < 3) ?
"Rjlc" :
"c";
299 <<
"solution\n" << mesh << x << flush
300 <<
"keys " << glvis_keys << endl
301 <<
"window_title 'Deformed configuration'" << endl
302 <<
"plot_caption 'Backward displacement'" << endl
306 const char *c =
"xyz";
309 StressCoefficient stress_c(lambda_c, mu_c);
312 stress_c.SetDisplacement(x);
313 for (
int si = 0; si <
dim; si++)
315 for (
int sj = si; sj <
dim; sj++)
317 stress_c.SetComponent(si, sj);
321 <<
"solution\n" << mesh << stress << flush
322 <<
"keys " << glvis_keys << endl
323 <<
"window_title |Stress " << c[si] << c[sj] <<
'|' << endl
336 u(u.
Size()-1) = -0.2*x(0);
343 MFEM_ASSERT(u != NULL,
"displacement field is not set");
345 double L = lambda.Eval(T, ip);
346 double M =
mu.Eval(T, ip);
347 u->GetVectorGradient(T, grad);
350 double div_u = grad.Trace();
351 return L*div_u + 2*M*grad(si,si);
355 return M*(grad(si,sj) + grad(sj,si));
362 host(vishost), port(visport), sid(0)
368 win_stride_x = win_w;
369 win_stride_y = win_h + 20;
373 void VisMan::NewWindow()
377 iostream::rdbuf(sock[sid]->rdbuf());
380 void VisMan::CloseConnection()
382 if (sid < sock.Size())
390 void VisMan::PositionWindow()
392 *
this <<
"window_geometry "
393 << win_x + win_stride_x*(sid%win_nx) <<
' '
394 << win_y + win_stride_y*(sid/win_nx) <<
' '
395 << win_w <<
' ' << win_h << endl;
400 for (
int i = sock.Size()-1; i >= 0; i--)
408 VisMan *vp =
dynamic_cast<VisMan*
>(&v);
409 if (vp) { (*f)(*vp); }
double Eval(ElementTransformation &T, const IntegrationPoint &ip, double t)
Evaluate the coefficient in the element described by T at the point ip at time t. ...
virtual void Print(std::ostream &out=mfem::out) const
Class for grid function - Vector with associated FE space.
Data type dense matrix using column-major storage.
int Size() const
Returns the size of the vector.
int GetNE() const
Returns number of elements.
void PrintInfo(std::ostream &out) const
Print various sparse matrix statistics.
int main(int argc, char *argv[])
Data type for Gauss-Seidel smoother of sparse matrix.
void position_window(VisMan &v)
Direct sparse solver using UMFPACK.
virtual void Save(std::ostream &out) const
Save the GridFunction to an output stream.
void Parse()
Parse the command-line options. Note that this function expects all the options provided through the ...
void new_window(VisMan &v)
void UniformRefinement(int i, const DSTable &, int *, int *, int *)
virtual void SetCurvature(int order, bool discont=false, int space_dim=-1, int ordering=1)
T Max() const
Find the maximal element in the array, using the comparison operator < for class T.
void PCG(const Operator &A, Solver &B, const Vector &b, Vector &x, int print_iter, int max_num_iter, double RTOLERANCE, double ATOLERANCE)
Preconditioned conjugate gradient method. (tolerances are squared)
virtual int GetTrueVSize() const
Return the number of vector true (conforming) dofs.
void SetNodalFESpace(FiniteElementSpace *nfes)
void PrintUsage(std::ostream &out) const
Print the usage message.
A general vector function coefficient.
void close_connection(VisMan &v)
Array< int > bdr_attributes
A list of all unique boundary attributes used by the Mesh.
double Control[UMFPACK_CONTROL]
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
Base class Coefficients that optionally depend on space and time. These are used by the BilinearFormI...
void AddOption(bool *var, const char *enable_short_name, const char *enable_long_name, const char *disable_short_name, const char *disable_long_name, const char *description, bool required=false)
Add a boolean option and set 'var' to receive the value. Enable/disable tags are used to set the bool...
NURBSExtension * NURBSext
Optional NURBS mesh extension.
A piecewise constant coefficient with the constants keyed off the element attribute numbers...
Class for integration point with weight.
void InitDisplacement(const Vector &x, Vector &u)
void PrintOptions(std::ostream &out) const
Print the options.
virtual void ProjectCoefficient(Coefficient &coeff)
std::ostream & operator<<(std::ostream &os, SparseMatrix const &mat)
void GetNodes(Vector &node_coord) const
Array< int > attributes
A list of all unique element attributes used by the Mesh.
virtual void Mult(const Vector &b, Vector &x) const
Operator application: y=A(x).
int GMRES(const Operator &A, Vector &x, const Vector &b, Solver &M, int &max_iter, int m, double &tol, double atol, int printit)
GMRES method. (tolerances are squared)
Arbitrary order "L2-conforming" discontinuous finite elements.
double f(const Vector &p)
void Neg()
(*this) = -(*this)
virtual void SetOperator(const Operator &op)
Factorize the given Operator op which must be a SparseMatrix.
bool Good() const
Return true if the command line options were parsed successfully.