MFEM  v4.6.0
Finite element discretization library
generate_random_field.cpp
Go to the documentation of this file.
1 // Copyright (c) 2010-2023, 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 // -------------------------------------------------------------------
13 // Gaussian Random Fields of Matern Covariance for Imperfect Materials
14 // -------------------------------------------------------------------
15 //
16 // See README.md for detailed description.
17 //
18 // Compile with: make generate_random_field
19 //
20 // Sample runs:
21 // (Basic usage)
22 // mpirun -np 4 generate_random_field
23 //
24 // (Generate 5 particles with random imperfections)
25 // mpirun -np 4 generate_random_field -o 1 -r 3 -rp 3 -nu 2 -l1 0.015 -l2 0.015 -l3 0.015 -s 0.01 -t 0.08 -n 5 -pl2 3 -top 0 -rs
26 //
27 // (Generate an Octet-Truss with random imperfections)
28 // mpirun -np 4 generate_random_field -o 1 -r 3 -rp 3 -nu 2 -l1 0.02 -l2 0.02 -l3 0.02 -s 0.01 -t 0.08 -top 1 -rs
29 //
30 // (Generate an Octet-Truss with random imperfections following a uniform distribution)
31 // mpirun -np 4 generate_random_field -o 1 -r 3 -rp 3 -nu 2 -l1 0.02 -l2 0.02 -l3 0.02 -umin 0.01 -umax 0.05 -t 0.08 -top 1 -urf -rs
32 //
33 // (2D random field with anisotropy)
34 // mpirun -np 4 generate_random_field -o 1 -r 3 -rp 3 -nu 4 -l1 0.09 -l2 0.03 -l3 0.05 -s 0.01 -t 0.08 -top 1 -no-rs -m ../../data/ref-square.mesh
35 
36 #include <math.h>
37 #include <fstream>
38 #include <iostream>
39 #include <string>
40 #include "mfem.hpp"
41 
42 #include "material_metrics.hpp"
43 #include "spde_solver.hpp"
44 #include "transformation.hpp"
45 #include "util.hpp"
46 #include "visualizer.hpp"
47 
48 using namespace std;
49 using namespace mfem;
50 
52 
53 int main(int argc, char *argv[])
54 {
55  // 0. Initialize MPI.
56  Mpi::Init(argc, argv);
57  Hypre::Init();
58 
59  // 1. Parse command-line options.
60  const char *mesh_file = "../../data/ref-cube.mesh";
61  int order = 1;
62  int num_refs = 3;
63  int num_parallel_refs = 3;
64  int number_of_particles = 3;
65  int topological_support = TopologicalSupport::kOctetTruss;
66  double nu = 2.0;
67  double tau = 0.08;
68  double l1 = 0.02;
69  double l2 = 0.02;
70  double l3 = 0.02;
71  double e1 = 0;
72  double e2 = 0;
73  double e3 = 0;
74  double pl1 = 1.0;
75  double pl2 = 1.0;
76  double pl3 = 1.0;
77  double uniform_min = 0.0;
78  double uniform_max = 1.0;
79  double offset = 0.0;
80  double scale = 0.01;
81  double level_set_threshold = 0.0;
82  bool paraview_export = true;
83  bool glvis_export = true;
84  bool uniform_rf = false;
85  bool random_seed = true;
86  bool compute_boundary_integrals = false;
87 
88  OptionsParser args(argc, argv);
89  args.AddOption(&mesh_file, "-m", "--mesh", "Mesh file to use.");
90  args.AddOption(&order, "-o", "--order",
91  "Finite element order (polynomial degree) or -1 for"
92  " isoparametric space.");
93  args.AddOption(&num_refs, "-r", "--refs", "Number of uniform refinements");
94  args.AddOption(&num_parallel_refs, "-rp", "--refs-parallel",
95  "Number of uniform refinements");
96  args.AddOption(&topological_support, "-top", "--topology",
97  "Topological support. 0 particles, 1 octet-truss");
98  args.AddOption(&nu, "-nu", "--nu", "Fractional exponent nu (smoothness)");
99  args.AddOption(&tau, "-t", "--tau", "Parameter for topology generation");
100  args.AddOption(&l1, "-l1", "--l1",
101  "First component of diagonal core of theta");
102  args.AddOption(&l2, "-l2", "--l2",
103  "Second component of diagonal core of theta");
104  args.AddOption(&l3, "-l3", "--l3",
105  "Third component of diagonal core of theta");
106  args.AddOption(&e1, "-e1", "--e1", "First euler angle for rotation of theta");
107  args.AddOption(&e2, "-e2", "--e2",
108  "Second euler angle for rotation of theta");
109  args.AddOption(&e3, "-e3", "--e3", "Third euler angle for rotation of theta");
110  args.AddOption(&pl1, "-pl1", "--pl1", "Length scale 1 of particles");
111  args.AddOption(&pl2, "-pl2", "--pl2", "Length scale 2 of particles");
112  args.AddOption(&pl3, "-pl3", "--pl3", "Length scale 3 of particles");
113  args.AddOption(&uniform_min, "-umin", "--uniform-min",
114  "Minimum value of uniform distribution");
115  args.AddOption(&uniform_max, "-umax", "--uniform-max",
116  "Maximum value of uniform distribution");
117  args.AddOption(&offset, "-off", "--offset",
118  "Offset for random field u(x) -> u(x) + a");
119  args.AddOption(&scale, "-s", "--scale",
120  "Scale for random field u(x) -> a * u(x)");
121  args.AddOption(&level_set_threshold, "-lst", "--level-set-threshold",
122  "Level set threshold");
123  args.AddOption(&number_of_particles, "-n", "--number-of-particles",
124  "Number of particles");
125  args.AddOption(&paraview_export, "-pvis", "--paraview-visualization",
126  "-no-pvis", "--no-paraview-visualization",
127  "Enable or disable ParaView visualization.");
128  args.AddOption(&glvis_export, "-vis", "--visualization", "-no-vis",
129  "--no-visualization",
130  "Enable or disable GLVis visualization.");
131  args.AddOption(&uniform_rf, "-urf", "--uniform-rf", "-no-urf",
132  "--no-uniform-rf",
133  "Enable or disable the transformation of GRF to URF.");
134  args.AddOption(&random_seed, "-rs", "--random-seed", "-no-rs",
135  "--no-random-seed", "Enable or disable random seed.");
136  args.AddOption(&compute_boundary_integrals, "-cbi",
137  "--compute-boundary-integrals", "-no-cbi",
138  "--no-compute-boundary-integrals",
139  "Enable or disable computation of boundary integrals.");
140  args.Parse();
141  if (!args.Good())
142  {
143  args.PrintUsage(cout);
144  return 1;
145  }
146  if (Mpi::Root())
147  {
148  args.PrintOptions(cout);
149  }
150 
151  // 2. Read the mesh from the given mesh file.
152  Mesh mesh(mesh_file, 1, 1);
153  int dim = mesh.Dimension();
154  bool is_3d = (dim == 3);
155 
156  // 3. Refine the mesh to increase the resolution.
157  for (int i = 0; i < num_refs; i++)
158  {
159  mesh.UniformRefinement();
160  }
161  ParMesh pmesh(MPI_COMM_WORLD, mesh);
162  mesh.Clear();
163  for (int i = 0; i < num_parallel_refs; i++)
164  {
165  pmesh.UniformRefinement();
166  }
167 
168  // 4. Define a finite element space on the mesh.
169  H1_FECollection fec(order, dim);
170  ParFiniteElementSpace fespace(&pmesh, &fec);
171  HYPRE_BigInt size = fespace.GlobalTrueVSize();
172  if (Mpi::Root())
173  {
174  const Array<int> boundary(pmesh.bdr_attributes);
175  cout << "Number of finite element unknowns: " << size << "\n";
176  cout << "Boundary attributes: ";
177  boundary.Print(cout, 6);
178  }
179 
180  // ========================================================================
181  // II. Generate topological support
182  // ========================================================================
183 
184  ParGridFunction v(&fespace);
185  v = 0.0;
186  MaterialTopology *mdm = nullptr;
187 
188  // II.1 Define the metric for the topological support.
189  if (is_3d)
190  {
191  if (topological_support == TopologicalSupport::kOctetTruss)
192  {
193  mdm = new OctetTrussTopology();
194  }
195  else if (topological_support == TopologicalSupport::kParticles)
196  {
197  // Create the same random particles on all processors.
198  std::vector<double> random_positions(3 * number_of_particles);
199  std::vector<double> random_rotations(9 * number_of_particles);
200  if (Mpi::Root())
201  {
202  // Generate random positions and rotations. We generate them on the root
203  // process and then broadcast them to all processes because we need the
204  // same random positions and rotations on all processes.
205  FillWithRandomNumbers(random_positions, 0.2, 0.8);
206  FillWithRandomRotations(random_rotations);
207  }
208 
209  // Broadcast the random positions and rotations to all processes.
210  MPI_Bcast(random_positions.data(), 3 * number_of_particles, MPI_DOUBLE, 0,
211  MPI_COMM_WORLD);
212  MPI_Bcast(random_rotations.data(), 9 * number_of_particles, MPI_DOUBLE, 0,
213  MPI_COMM_WORLD);
214 
215  mdm = new ParticleTopology(pl1, pl2, pl3, random_positions,
216  random_rotations);
217  }
218  else
219  {
220  if (Mpi::Root())
221  {
222  mfem::out << "Error: Selected topological support not valid."
223  << std::endl;
224  }
225  return 1;
226  }
227 
228  // II.2 Define lambda to wrap the call to the distance metric.
229  auto topo = [&mdm, &tau](const Vector &x)
230  {
231  return (tau - mdm->ComputeMetric(x));
232  };
233 
234  // II.3 Create a GridFunction for the topological support.
235  FunctionCoefficient topo_coeff(topo);
236  v.ProjectCoefficient(topo_coeff);
237  }
238 
239  // ========================================================================
240  // III. Generate random imperfections via fractional PDE
241  // ========================================================================
242 
243  /// III.1 Define the fractional PDE solution
244  ParGridFunction u(&fespace);
245  u = 0.0;
246 
247  // III.2 Define the boundary conditions.
248  spde::Boundary bc;
249  if (Mpi::Root())
250  {
251  bc.PrintInfo();
252  bc.VerifyDefinedBoundaries(pmesh);
253  }
254 
255  // III.3 Solve the SPDE problem
256  spde::SPDESolver solver(nu, bc, &fespace, l1, l2, l3, e1, e2,
257  e3);
258  const int seed = (random_seed) ? 0 : std::numeric_limits<int>::max();
259  solver.SetupRandomFieldGenerator(seed);
260  solver.GenerateRandomField(u);
261 
262  /// III.4 Verify boundary conditions
263  if (compute_boundary_integrals)
264  {
266  }
267 
268  // ========================================================================
269  // III. Combine topological support and random field
270  // ========================================================================
271 
272  if (uniform_rf)
273  {
274  /// Transform the random field to a uniform random field.
275  spde::UniformGRFTransformer transformation(uniform_min, uniform_max);
276  transformation.Transform(u);
277  }
278  if (scale != 1.0)
279  {
280  /// Scale the random field.
282  transformation.Transform(u);
283  }
284  if (offset != 0.0)
285  {
286  /// Add an offset to the random field.
288  transformation.Transform(u);
289  }
290  ParGridFunction w(&fespace); // Noisy material field.
291  w = 0.0;
292  w += u;
293  w += v;
294  ParGridFunction level_set(w); // Level set field.
295  {
296  spde::LevelSetTransformer transformation(level_set_threshold);
297  transformation.Transform(level_set);
298  }
299 
300  // ========================================================================
301  // IV. Export visualization to ParaView and GLVis
302  // ========================================================================
303 
304  spde::Visualizer vis(pmesh, order, u, v, w, level_set, is_3d);
305  if (paraview_export)
306  {
307  vis.ExportToParaView();
308  }
309  if (glvis_export)
310  {
311  vis.SendToGLVis();
312  }
313 
314  delete mdm;
315  return 0;
316 }
void ComputeBoundaryError(const ParGridFunction &solution)
Level Set Transformer, 1 for u(x) >= threshold, 0 otherwise.
Virtual class to define the interface for defining the material topology.
Transforms a grid function by scaling it by a constant factor.
void PrintOptions(std::ostream &out) const
Print the options.
Definition: optparser.cpp:331
int Dimension() const
Dimension of the reference space used within the elements.
Definition: mesh.hpp:1020
void SendToGLVis() const
Definition: visualizer.cpp:39
void PrintUsage(std::ostream &out) const
Print the usage message.
Definition: optparser.cpp:462
bool Good() const
Return true if the command line options were parsed successfully.
Definition: optparser.hpp:159
Abstract parallel finite element space.
Definition: pfespace.hpp:28
virtual void ProjectCoefficient(Coefficient &coeff)
Project coeff Coefficient to this GridFunction. The projection computation depends on the choice of t...
Definition: pgridfunc.cpp:561
STL namespace.
void PrintInfo(std::ostream &os=mfem::out) const
Print the information specifying the boundary conditions.
Definition: spde_solver.cpp:31
void GenerateRandomField(ParGridFunction &x)
void Parse()
Parse the command-line options. Note that this function expects all the options provided through the ...
Definition: optparser.cpp:151
Class that implements the particle topology.
virtual double ComputeMetric(const Vector &x)=0
Compute the metric rho describing the material topology.
void transformation(const Vector &p, Vector &v)
void UniformRefinement(int i, const DSTable &, int *, int *, int *)
Definition: mesh.cpp:10232
void SetupRandomFieldGenerator(int seed=0)
Adds an constant offset to a grid function, i.e. u(x) = u(x) + offset.
void Print(std::ostream &out=mfem::out, int width=4) const
Prints array to stream with width elements per row.
Definition: array.cpp:23
void FillWithRandomRotations(std::vector< double > &x)
Definition: util.cpp:26
HYPRE_BigInt GlobalTrueVSize() const
Definition: pfespace.hpp:281
Array< int > bdr_attributes
A list of all unique boundary attributes used by the Mesh.
Definition: mesh.hpp:275
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
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 &#39;var&#39; to receive the value. Enable/disable tags are used to set the bool...
Definition: optparser.hpp:82
HYPRE_Int HYPRE_BigInt
void FillWithRandomNumbers(std::vector< double > &x, double a, double b)
Fills the vector x with random numbers between a and b.
Definition: util.cpp:18
int dim
Definition: ex24.cpp:53
void Clear()
Clear the contents of the Mesh.
Definition: mesh.hpp:678
A general function coefficient.
Vector data type.
Definition: vector.hpp:58
Arbitrary order H1-conforming (continuous) finite elements.
Definition: fe_coll.hpp:259
double u(const Vector &xvec)
Definition: lor_mms.hpp:22
Class for parallel grid function.
Definition: pgridfunc.hpp:32
void VerifyDefinedBoundaries(const Mesh &mesh) const
Definition: spde_solver.cpp:80
Class for parallel meshes.
Definition: pmesh.hpp:32
int main(int argc, char *argv[])