53int main(
int argc,
char *argv[])
56 cout <<
"This miniapp is not supported in single precision.\n\n";
57 return MFEM_SKIP_RETURN_VALUE;
65 const char *mesh_file =
"../../data/ref-cube.mesh";
68 int num_parallel_refs = 3;
69 int number_of_particles = 3;
86 real_t level_set_threshold = 0.0;
87 bool paraview_export =
true;
88 bool glvis_export =
true;
89 bool uniform_rf =
false;
90 bool random_seed =
true;
91 bool compute_boundary_integrals =
false;
94 args.
AddOption(&mesh_file,
"-m",
"--mesh",
"Mesh file to use.");
96 "Finite element order (polynomial degree) or -1 for"
97 " isoparametric space.");
98 args.
AddOption(&num_refs,
"-r",
"--refs",
"Number of uniform refinements");
99 args.
AddOption(&num_parallel_refs,
"-rp",
"--refs-parallel",
100 "Number of uniform refinements");
101 args.
AddOption(&topological_support,
"-top",
"--topology",
102 "Topological support. 0 particles, 1 octet-truss");
103 args.
AddOption(&nu,
"-nu",
"--nu",
"Fractional exponent nu (smoothness)");
104 args.
AddOption(&tau,
"-t",
"--tau",
"Parameter for topology generation");
106 "First component of diagonal core of theta");
108 "Second component of diagonal core of theta");
110 "Third component of diagonal core of theta");
111 args.
AddOption(&e1,
"-e1",
"--e1",
"First euler angle for rotation of theta");
113 "Second euler angle for rotation of theta");
114 args.
AddOption(&e3,
"-e3",
"--e3",
"Third euler angle for rotation of theta");
115 args.
AddOption(&pl1,
"-pl1",
"--pl1",
"Length scale 1 of particles");
116 args.
AddOption(&pl2,
"-pl2",
"--pl2",
"Length scale 2 of particles");
117 args.
AddOption(&pl3,
"-pl3",
"--pl3",
"Length scale 3 of particles");
118 args.
AddOption(&uniform_min,
"-umin",
"--uniform-min",
119 "Minimum value of uniform distribution");
120 args.
AddOption(&uniform_max,
"-umax",
"--uniform-max",
121 "Maximum value of uniform distribution");
122 args.
AddOption(&offset,
"-off",
"--offset",
123 "Offset for random field u(x) -> u(x) + a");
125 "Scale for random field u(x) -> a * u(x)");
126 args.
AddOption(&level_set_threshold,
"-lst",
"--level-set-threshold",
127 "Level set threshold");
128 args.
AddOption(&number_of_particles,
"-n",
"--number-of-particles",
129 "Number of particles");
130 args.
AddOption(¶view_export,
"-pvis",
"--paraview-visualization",
131 "-no-pvis",
"--no-paraview-visualization",
132 "Enable or disable ParaView visualization.");
133 args.
AddOption(&glvis_export,
"-vis",
"--visualization",
"-no-vis",
134 "--no-visualization",
135 "Enable or disable GLVis visualization.");
136 args.
AddOption(&uniform_rf,
"-urf",
"--uniform-rf",
"-no-urf",
138 "Enable or disable the transformation of GRF to URF.");
139 args.
AddOption(&random_seed,
"-rs",
"--random-seed",
"-no-rs",
140 "--no-random-seed",
"Enable or disable random seed.");
141 args.
AddOption(&compute_boundary_integrals,
"-cbi",
142 "--compute-boundary-integrals",
"-no-cbi",
143 "--no-compute-boundary-integrals",
144 "Enable or disable computation of boundary integrals.");
157 Mesh mesh(mesh_file, 1, 1);
159 bool is_3d = (
dim == 3);
162 for (
int i = 0; i < num_refs; i++)
166 ParMesh pmesh(MPI_COMM_WORLD, mesh);
168 for (
int i = 0; i < num_parallel_refs; i++)
180 cout <<
"Number of finite element unknowns: " << size <<
"\n";
181 cout <<
"Boundary attributes: ";
182 boundary.
Print(cout, 6);
203 std::vector<real_t> random_positions(3 * number_of_particles);
204 std::vector<real_t> random_rotations(9 * number_of_particles);
215 MPI_Bcast(random_positions.data(), 3 * number_of_particles,
217 MPI_Bcast(random_rotations.data(), 9 * number_of_particles,
227 mfem::out <<
"Error: Selected topological support not valid."
234 auto topo = [&mdm, &tau](
const Vector &x)
263 const int seed = (random_seed) ? 0 : std::numeric_limits<int>::max();
268 if (compute_boundary_integrals)
void Print(std::ostream &out=mfem::out, int width=4) const
Prints array to stream with width elements per row.
A general function coefficient.
Arbitrary order H1-conforming (continuous) finite elements.
static void Init()
Initialize hypre by calling HYPRE_Init() and set default options. After calling Hypre::Init(),...
Virtual class to define the interface for defining the material topology.
virtual real_t ComputeMetric(const Vector &x)=0
Compute the metric rho describing the material topology.
Array< int > bdr_attributes
A list of all unique boundary attributes used by the Mesh.
void Clear()
Clear the contents of the Mesh.
int Dimension() const
Dimension of the reference space used within the elements.
void UniformRefinement(int i, const DSTable &, int *, int *, int *)
static bool Root()
Return true if the rank in MPI_COMM_WORLD is zero.
static void Init(int &argc, char **&argv, int required=default_thread_required, int *provided=nullptr)
Singleton creation with Mpi::Init(argc, argv).
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.
Abstract parallel finite element space.
HYPRE_BigInt GlobalTrueVSize() const
Class for parallel grid function.
void ProjectCoefficient(Coefficient &coeff) override
Project coeff Coefficient to this GridFunction. The projection computation depends on the choice of t...
Class for parallel meshes.
Class that implements the particle topology.
void GenerateRandomField(ParGridFunction &x)
void SetupRandomFieldGenerator(int seed=0)
void transformation(const Vector &p, Vector &v)
real_t u(const Vector &xvec)
void FillWithRandomRotations(std::vector< real_t > &x)
OutStream out(std::cout)
Global stream used by the library for standard output. Initially it uses the same std::streambuf as s...
void FillWithRandomNumbers(std::vector< real_t > &x, real_t a, real_t b)
Fills the vector x with random numbers between a and b.
Helper struct to convert a C++ type to an MPI type.
void PrintInfo(std::ostream &os=mfem::out) const
Print the information specifying the boundary conditions.
void ComputeBoundaryError(const ParGridFunction &solution)
void VerifyDefinedBoundaries(const Mesh &mesh) const