MFEM  v4.5.1
Finite element discretization library
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
hooke.cpp
Go to the documentation of this file.
1 // Copyright (c) 2010-2022, Lawrence Livermore National Security, LLC. Produced
2 // at the Lawrence Livermore National Laboratory. All Rights reserved. See files
3 // LICENSE and NOTICE for details. LLNL-CODE-806117.
4 //
5 // This file is part of the MFEM library. For more information and source code
6 // availability visit https://mfem.org.
7 //
8 // MFEM is free software; you can redistribute it and/or modify it under the
9 // terms of the BSD-3 license. We welcome feedback and contributions, see file
10 // CONTRIBUTING.md for details.
11 
12 // This miniapp solves a quasistatic solid mechanics problem assuming an elastic
13 // material and no body forces.
14 //
15 // The equation
16 // ∇⋅σ(∇u) = 0
17 //
18 // with stress σ is solved for displacement u.
19 //
20 // +----------+----------+
21 // fixed --->| |<--- constant displacement
22 // | |
23 // +----------+----------+
24 //
25 // This miniapp uses an elasticity operator that allows for a custom material.
26 // By default the NeoHookeanMaterial is used. A linear elastic material is also
27 // provided. Based on these examples, other materials could be implemented.
28 //
29 // The implementation of NeoHookeanMaterial also demonstrates the use of
30 // automatic differentiation using either a native dual number forward mode
31 // implementation or leveraging the Enzyme third party library.
32 
33 #include <mfem.hpp>
34 
36 #include "materials/neohookean.hpp"
40 
41 using namespace std;
42 using namespace mfem;
43 
44 /// This example only works in 3D. Kernels for 2D are not implemented.
45 constexpr int dimension = 3;
46 
47 void display_banner(ostream& os)
48 {
49  os << R"( ___ ___ ________ ________ ____ __.___________ / | \\_____ \ \_____ \ | |/ _|\_ _____/ / ~ \/ | \ / | \| < | __)_ \ Y / | \/ | \ | \ | \ \___|_ /\_______ /\_______ /____|__ \/_______ / \/ \/ \/ \/ \/ )"
50  << endl << flush;
51 }
52 
53 int main(int argc, char *argv[])
54 {
55  Mpi::Init(argc, argv);
56  int num_procs = Mpi::WorldSize();
57  int myid = Mpi::WorldRank();
58 
59  int order = 1;
60  const char *device_config = "cpu";
61  int diagpc_type = ElasticityDiagonalPreconditioner::Type::Diagonal;
62  int serial_refinement_levels = 0;
63  bool visualization = true;
64  bool paraview = false;
65 
66  if (Mpi::Root())
67  {
69  }
70 
71  OptionsParser args(argc, argv);
72  args.AddOption(&order, "-o", "--order",
73  "Finite element order (polynomial degree).");
74  args.AddOption(&device_config, "-d", "--device",
75  "Device configuration string, see Device::Configure().");
76  args.AddOption(&diagpc_type, "-pc", "--pctype",
77  "Select diagonal preconditioner type"
78  " (0:Diagonal, 1:BlockDiagonal).");
79  args.AddOption(&serial_refinement_levels, "-rs", "--ref-serial",
80  "Number of uniform refinements on the serial mesh.");
81  args.AddOption(&visualization, "-vis", "--visualization", "-no-vis",
82  "--no-visualization",
83  "Enable or disable GLVis visualization.");
84  args.AddOption(&paraview, "-pv", "--paraview", "-no-pv",
85  "--no-paraview",
86  "Enable or disable ParaView DataCollection output.");
87  args.ParseCheck();
88 
89  Device device(device_config);
90  if (Mpi::Root())
91  {
92  device.Print();
93  }
94 
95  auto mesh =
96  Mesh::MakeCartesian3D(8, 2, 2, Element::HEXAHEDRON, 8.0, 1.0, 1.0);
97  if (mesh.Dimension() != dimension)
98  {
99  MFEM_ABORT("This example only works in 3D.");
100  }
101  mesh.EnsureNodes();
102 
103  for (int l = 0; l < serial_refinement_levels; l++)
104  {
105  mesh.UniformRefinement();
106  }
107 
108  ParMesh pmesh(MPI_COMM_WORLD, mesh);
109  mesh.Clear();
110 
111  // Create the elasticity operator on the parallel mesh.
112  ElasticityOperator elasticity_op(pmesh, order);
113 
114  // Create and set the material type. We define its GradientType during
115  // instantiation.
116 
117  // As seen in materials/gradient_type.hpp there is a choice of the
118  // GradientType with either
119  // * Symbolic (Manually derived)
120  // * EnzymeFwd
121  // * EnzymeRev
122  // * FiniteDiff
123  // * InternalFwd
125  elasticity_op.SetMaterial(material);
126 
127  // Define all essential boundaries. In this specific example, this includes
128  // all fixed and statically displaced degrees of freedom on mesh entities in
129  // the defined attributes.
130  if (pmesh.bdr_attributes.Size())
131  {
132  Array<int> ess_attr(pmesh.bdr_attributes.Max());
133  ess_attr = 0;
134  ess_attr[4] = 1;
135  ess_attr[2] = 1;
136  elasticity_op.SetEssentialAttributes(ess_attr);
137  }
138 
139  // Define all statically displaced degrees of freedom on mesh entities in the
140  // defined attributes. On these degrees of freedom (determined from the mesh
141  // attributes), a fixed displacement is prescribed.
142  if (pmesh.bdr_attributes.Size())
143  {
144  Array<int> displaced_attr(pmesh.bdr_attributes.Max());
145  displaced_attr = 0;
146  displaced_attr[2] = 1;
147  elasticity_op.SetPrescribedDisplacement(displaced_attr);
148  }
149 
150  ParGridFunction U_gf(&elasticity_op.h1_fes_);
151  U_gf = 0.0;
152 
153  Vector U;
154  U_gf.GetTrueDofs(U);
155 
156  // Prescribe a fixed displacement to the displaced degrees of freedom.
157  U.SetSubVector(elasticity_op.GetPrescribedDisplacementTDofs(), 1.0e-2);
158 
159  // Define the type of preconditioner to use for the linear solver.
161  static_cast<ElasticityDiagonalPreconditioner::Type>(diagpc_type));
162 
163  CGSolver cg(MPI_COMM_WORLD);
164  cg.SetRelTol(1e-1);
165  cg.SetMaxIter(10000);
166  cg.SetPrintLevel(2);
167  cg.SetPreconditioner(diagonal_pc);
168 
169  NewtonSolver newton(MPI_COMM_WORLD);
170  newton.SetSolver(cg);
171  newton.SetOperator(elasticity_op);
172  newton.SetRelTol(1e-6);
173  newton.SetMaxIter(10);
174  newton.SetPrintLevel(1);
175 
176  Vector zero;
177  newton.Mult(zero, U);
178 
179  U_gf.Distribute(U);
180 
181  if (visualization)
182  {
183  char vishost[] = "localhost";
184  int visport = 19916;
185  socketstream sol_sock(vishost, visport);
186  sol_sock << "parallel " << num_procs << " " << myid << "\n";
187  sol_sock.precision(8);
188  sol_sock << "solution\n" << pmesh << U_gf << flush;
189  }
190 
191  if (paraview)
192  {
193  ParaViewDataCollection pd("elasticity_output", &pmesh);
194  pd.RegisterField("solution", &U_gf);
195  pd.SetLevelsOfDetail(order);
196  pd.SetDataFormat(VTKFormat::BINARY);
197  pd.SetHighOrderOutput(true);
198  pd.SetCycle(0);
199  pd.SetTime(0.0);
200  pd.Save();
201  }
202 
203  return 0;
204 }
205 
void SetSubVector(const Array< int > &dofs, const double value)
Set the entries listed in dofs to the given value.
Definition: vector.cpp:573
Conjugate gradient method.
Definition: solvers.hpp:465
Helper class for ParaView visualization data.
Neo-Hookean material.
Definition: neohookean.hpp:40
int material(Vector &x, Vector &xmin, Vector &xmax)
Definition: shaper.cpp:53
constexpr char vishost[]
constexpr int visport
Newton&#39;s method for solving F(x)=b for a given operator F.
Definition: solvers.hpp:613
constexpr int dimension
This example only works in 3D. Kernels for 2D are not implemented.
Definition: hooke.cpp:45
Vector data type.
Definition: vector.hpp:60
ElasticityDiagonalPreconditioner acts as a matrix-free preconditioner for ElasticityOperator.
Class for parallel grid function.
Definition: pgridfunc.hpp:32
OutStream out(std::cout)
Global stream used by the library for standard output. Initially it uses the same std::streambuf as s...
Definition: globals.hpp:66
The MFEM Device class abstracts hardware devices such as GPUs, as well as programming models such as ...
Definition: device.hpp:121
void display_banner(ostream &os)
Definition: joule.cpp:778
Class for parallel meshes.
Definition: pmesh.hpp:32
int main()