67 : lambda(lambda_),
mu(mu_),
u(NULL), si(0), sj(0) { }
70 void SetComponent(
int i,
int j) { si = i; sj = j; }
76class VisMan :
public iostream
84 int win_x, win_y, win_w, win_h;
85 int win_stride_x, win_stride_y, win_nx;
88 VisMan(
const char *
vishost,
const int visport);
90 void CloseConnection();
91 void PositionWindow();
101int 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;
212 b.AddBdrFaceIntegrator(
214 init_x, lambda_c, mu_c,
alpha,
kappa), dir_bdr);
226 a.AddInteriorFaceIntegrator(
228 a.AddBdrFaceIntegrator(
232 cout <<
"matrix ... " << flush;
237 a.FormLinearSystem(ess_tdof_list, x,
b, A, X, B);
238 cout <<
"done." << endl;
243#ifndef MFEM_USE_SUITESPARSE
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);
266 a.RecoverFEMSolution(X,
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);
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");
351 return L*div_u + 2*M*grad(si,si);
355 return M*(grad(si,sj) + grad(sj,si));
360VisMan::VisMan(
const char *
vishost,
const int visport)
362 host(
vishost), port(visport), sid(0)
368 win_stride_x = win_w;
369 win_stride_y = win_h + 20;
373void VisMan::NewWindow()
377 iostream::rdbuf(sock[sid]->rdbuf());
380void VisMan::CloseConnection()
382 if (sid < sock.Size())
390void 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); }
T Max() const
Find the maximal element in the array, using the comparison operator < for class T.
@ GaussLobatto
Closed type.
Base class Coefficients that optionally depend on space and time. These are used by the BilinearFormI...
virtual real_t Eval(ElementTransformation &T, const IntegrationPoint &ip)=0
Evaluate the coefficient in the element described by T at the point ip.
Data type dense matrix using column-major storage.
real_t Trace() const
Trace of a square matrix.
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
virtual int GetTrueVSize() const
Return the number of vector true (conforming) dofs.
Data type for Gauss-Seidel smoother of sparse matrix.
Class for grid function - Vector with associated FE space.
void GetVectorGradient(ElementTransformation &tr, DenseMatrix &grad) const
Compute the vector gradient with respect to the physical element variable.
virtual void Save(std::ostream &out) const
Save the GridFunction to an output stream.
virtual void ProjectCoefficient(Coefficient &coeff)
Project coeff Coefficient to this GridFunction. The projection computation depends on the choice of t...
Class for integration point with weight.
Arbitrary order "L2-conforming" discontinuous finite elements.
Array< int > bdr_attributes
A list of all unique boundary attributes used by the Mesh.
NURBSExtension * NURBSext
Optional NURBS mesh extension.
virtual void Print(std::ostream &os=mfem::out, const std::string &comments="") const
Print the mesh to the given stream using the default MFEM mesh format.
int GetNE() const
Returns number of elements.
int Dimension() const
Dimension of the reference space used within the elements.
void GetNodes(Vector &node_coord) const
virtual void SetNodalFESpace(FiniteElementSpace *nfes)
virtual void SetCurvature(int order, bool discont=false, int space_dim=-1, int ordering=1)
Set the curvature of the mesh nodes using the given polynomial degree.
void UniformRefinement(int i, const DSTable &, int *, int *, int *)
Array< int > attributes
A list of all unique element attributes used by the Mesh.
void Parse()
Parse the command-line options. Note that this function expects all the options provided through the ...
void PrintUsage(std::ostream &out) const
Print the usage message.
void PrintOptions(std::ostream &out) const
Print the options.
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...
bool Good() const
Return true if the command line options were parsed successfully.
A piecewise constant coefficient with the constants keyed off the element attribute numbers.
void PrintInfo(std::ostream &out) const
Print various sparse matrix statistics.
Direct sparse solver using UMFPACK.
real_t Control[UMFPACK_CONTROL]
void SetOperator(const Operator &op) override
Factorize the given Operator op which must be a SparseMatrix.
void Mult(const Vector &b, Vector &x) const override
Direct solution of the linear system using UMFPACK.
A general vector function coefficient.
void Neg()
(*this) = -(*this)
void close_connection(VisMan &v)
void InitDisplacement(const Vector &x, Vector &u)
void position_window(VisMan &v)
void new_window(VisMan &v)
std::ostream & operator<<(std::ostream &os, SparseMatrix const &mat)
real_t u(const Vector &xvec)
int GMRES(const Operator &A, Vector &x, const Vector &b, Solver &M, int &max_iter, int m, real_t &tol, real_t atol, int printit)
GMRES method. (tolerances are squared)
void PCG(const Operator &A, Solver &B, const Vector &b, Vector &x, int print_iter, int max_num_iter, real_t RTOLERANCE, real_t ATOLERANCE)
Preconditioned conjugate gradient method. (tolerances are squared)
std::function< real_t(const Vector &)> f(real_t mass_coeff)