MFEM v4.8.0
Finite element discretization library
Loading...
Searching...
No Matches
shaper.cpp
Go to the documentation of this file.
1// Copyright (c) 2010-2025, 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// Shaper Miniapp: Resolve material interfaces by mesh refinement
14// --------------------------------------------------------------
15//
16// This miniapp performs multiple levels of adaptive mesh refinement to resolve
17// the interfaces between different "materials" in the mesh, as specified by the
18// given material() function. It can be used as a simple initial mesh generator,
19// for example in the case when the interface is too complex to describe without
20// local refinement. Both conforming and non-conforming refinements are supported.
21//
22// Two additional versions of this miniapp can be found in the miniapps/toys
23// directory: Mandel uses the Shaper algorithm for fractal visualization, while
24// Mondrian convert an image to an AMR mesh suitable for MFEM computations.
25//
26// Compile with: make shaper
27//
28// Sample runs: shaper
29// shaper -m ../../data/inline-tri.mesh
30// shaper -m ../../data/inline-hex.mesh
31// shaper -m ../../data/inline-tet.mesh
32// shaper -m ../../data/amr-quad.mesh
33// shaper -m ../../data/beam-quad.mesh -a -ncl -1 -sd 4
34// shaper -m ../../data/ball-nurbs.mesh
35// shaper -m ../../data/mobius-strip.mesh
36// shaper -m ../../data/square-disc-surf.mesh
37// shaper -m ../../data/star-q3.mesh -sd 2 -ncl -1
38// shaper -m ../../data/fichera-amr.mesh -a -ncl -1
39
40#include "mfem.hpp"
41#include <fstream>
42#include <iostream>
43
44using namespace mfem;
45using namespace std;
46
47// Given a point x, return its material id as an integer. The ids should be
48// positive. If the point is exactly on the interface, return 0.
49//
50// This particular implementation, rescales the mesh to [-1,1]^sdim given the
51// xmin/xmax bounding box, and shapes in a simple annulus/shell with respect to
52// the rescaled coordinates.
53int material(Vector &x, Vector &xmin, Vector &xmax)
54{
55 static real_t p = 2.0;
56
57 // Rescaling to [-1,1]^sdim
58 for (int i = 0; i < x.Size(); i++)
59 {
60 x(i) = (2*x(i)-xmin(i)-xmax(i))/(xmax(i)-xmin(i));
61 }
62
63 // A simple annulus/shell
64 if (x.Normlp(p) > 0.4 && x.Normlp(p) < 0.6) { return 1; }
65 if (x.Normlp(p) < 0.4 || x.Normlp(p) > 0.6) { return 2; }
66 return 0;
67}
68
69int main(int argc, char *argv[])
70{
71 int sd = 2;
72 int nclimit = 1;
73 const char *mesh_file = "../../data/inline-quad.mesh";
74 bool aniso = false;
75 int visport = 19916;
76
77 // Parse command line
78 OptionsParser args(argc, argv);
79 args.AddOption(&mesh_file, "-m", "--mesh",
80 "Input mesh file to shape materials in.");
81 args.AddOption(&sd, "-sd", "--sub-divisions",
82 "Number of element subdivisions for interface detection.");
83 args.AddOption(&nclimit, "-ncl", "--nc-limit",
84 "Level of hanging nodes allowed (-1 = unlimited).");
85 args.AddOption(&aniso, "-a", "--aniso", "-i", "--iso",
86 "Enable anisotropic refinement of quads and hexes.");
87 args.AddOption(&visport, "-p", "--send-port", "Socket for GLVis.");
88 args.Parse();
89 if (!args.Good()) { args.PrintUsage(cout); return 1; }
90 args.PrintOptions(cout);
91
92 // Read initial mesh, get dimensions and bounding box
93 Mesh mesh(mesh_file, 1, 1);
94 int dim = mesh.Dimension();
95 int sdim = mesh.SpaceDimension();
96 Vector xmin, xmax;
97 mesh.GetBoundingBox(xmin, xmax);
98
99 // NURBS meshes don't support non-conforming refinement for now
100 if (mesh.NURBSext) { mesh.SetCurvature(2); }
101
102 // Anisotropic refinement not supported for simplex meshes.
103 if (mesh.MeshGenerator() & 1) { aniso = false; }
104
105 // Mesh attributes will be visualized as piece-wise constants
106 L2_FECollection attr_fec(0, dim);
107 FiniteElementSpace attr_fespace(&mesh, &attr_fec);
108 GridFunction attr(&attr_fespace);
109
110 // GLVis server to visualize to
111 char vishost[] = "localhost";
112 socketstream sol_sock(vishost, visport);
113 sol_sock.precision(8);
114
115 // Shaping loop
116 for (int iter = 0; 1; iter++)
117 {
119 for (int i = 0; i < mesh.GetNE(); i++)
120 {
121 bool refine = false;
122
123 // Sample materials in each element using "sd" sub-divisions
124 Vector pt;
127 RefinedGeometry *RefG = GlobGeometryRefiner.Refine(geom, sd, 1);
128 IntegrationRule &ir = RefG->RefPts;
129
130 // Refine any element where different materials are detected. A more
131 // sophisticated logic can be implemented here -- e.g. don't refine
132 // the interfaces between certain materials.
133 Array<int> mat(ir.GetNPoints());
134 real_t matsum = 0.0;
135 for (int j = 0; j < ir.GetNPoints(); j++)
136 {
137 T->Transform(ir.IntPoint(j), pt);
138 int m = material(pt, xmin, xmax);
139 mat[j] = m;
140 matsum += m;
141 if ((int)matsum != m*(j+1))
142 {
143 refine = true;
144 }
145 }
146
147 // Set the element attribute as the "average". Other choices are
148 // possible here too, e.g. attr(i) = mat;
149 attr(i) = round(matsum/ir.GetNPoints());
150
151 // Mark the element for refinement
152 if (refine)
153 {
154 int type = 7;
155 if (aniso)
156 {
157 // Determine the XYZ bitmask for anisotropic refinement.
158 int dx = 0, dy = 0, dz = 0;
159 const int s = sd+1;
160 if (dim == 2)
161 {
162 for (int jj = 0; jj <= sd; jj++)
163 for (int ii = 0; ii < sd; ii++)
164 {
165 dx += abs(mat[jj*s + ii+1] - mat[jj*s + ii]);
166 dy += abs(mat[(ii+1)*s + jj] - mat[ii*s + jj]);
167 }
168 }
169 else if (dim == 3)
170 {
171 for (int kk = 0; kk <= sd; kk++)
172 for (int jj = 0; jj <= sd; jj++)
173 for (int ii = 0; ii < sd; ii++)
174 {
175 dx += abs(mat[(kk*s + jj)*s + ii+1] - mat[(kk*s + jj)*s + ii]);
176 dy += abs(mat[(kk*s + ii+1)*s + jj] - mat[(kk*s + ii)*s + jj]);
177 dz += abs(mat[((ii+1)*s + jj)*s + kk] - mat[(ii*s + jj)*s + kk]);
178 }
179 }
180 type = 0;
181 const int tol = mat.Size() / 10;
182 if (dx > tol) { type |= 1; }
183 if (dy > tol) { type |= 2; }
184 if (dz > tol) { type |= 4; }
185 if (!type) { type = 7; } // because of tol
186 }
187
188 refs.Append(Refinement(i, type));
189 }
190 }
191
192 // Visualization
193 sol_sock << "solution\n" << mesh << attr;
194 if (iter == 0 && sdim == 2)
195 {
196 sol_sock << "keys 'RjlmpppppppppppppA*************'\n";
197 }
198 if (iter == 0 && sdim == 3)
199 {
200 sol_sock << "keys 'YYYYYYYYYXXXXXXXmA********8888888pppttt";
201 if (dim == 3) { sol_sock << "iiM"; }
202 sol_sock << "'\n";
203 }
204 sol_sock << flush;
205
206 // Ask the user if we should continue refining
207 char yn;
208 cout << "Mesh has " << mesh.GetNE() << " elements. \n"
209 << "Continue shaping? --> ";
210 cin >> yn;
211 if (yn == 'n' || yn == 'q') { break; }
212
213 // Perform refinement, update spaces and grid functions
214 mesh.GeneralRefinement(refs, -1, nclimit);
215 attr_fespace.Update();
216 attr.Update();
217 }
218
219 // Set element attributes in the mesh object before saving
220 for (int i = 0; i < mesh.GetNE(); i++)
221 {
222 mesh.SetAttribute(i, static_cast<int>(attr(i)));
223 }
224 mesh.SetAttributes();
225
226 // Save the final mesh
227 ofstream mesh_ofs("shaper.mesh");
228 mesh_ofs.precision(8);
229 mesh.Print(mesh_ofs);
230}
int Size() const
Return the logical size of the array.
Definition array.hpp:147
int Append(const T &el)
Append element 'el' to array, resize if necessary.
Definition array.hpp:830
virtual void Transform(const IntegrationPoint &, Vector &)=0
Transform integration point from reference coordinates to physical coordinates and store them in the ...
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
Definition fespace.hpp:244
virtual void Update(bool want_transform=true)
Reflect changes in the mesh: update number of DOFs, etc. Also, calculate GridFunction transformation ...
Definition fespace.cpp:4157
RefinedGeometry * Refine(Geometry::Type Geom, int Times, int ETimes=1)
Definition geom.cpp:1136
Class for grid function - Vector with associated FE space.
Definition gridfunc.hpp:31
Class for an integration rule - an Array of IntegrationPoint.
Definition intrules.hpp:100
int GetNPoints() const
Returns the number of the points in the integration rule.
Definition intrules.hpp:256
IntegrationPoint & IntPoint(int i)
Returns a reference to the i-th integration point.
Definition intrules.hpp:259
Arbitrary order "L2-conforming" discontinuous finite elements.
Definition fe_coll.hpp:346
Mesh data type.
Definition mesh.hpp:64
NURBSExtension * NURBSext
Optional NURBS mesh extension.
Definition mesh.hpp:298
void GeneralRefinement(const Array< Refinement > &refinements, int nonconforming=-1, int nc_limit=0)
Definition mesh.cpp:10883
void SetAttribute(int i, int attr)
Set the attribute of element i.
Definition mesh.cpp:7629
virtual void Print(std::ostream &os=mfem::out, const std::string &comments="") const
Definition mesh.hpp:2433
int GetNE() const
Returns number of elements.
Definition mesh.hpp:1282
void GetBoundingBox(Vector &min, Vector &max, int ref=2)
Returns the minimum and maximum corners of the mesh bounding box.
Definition mesh.cpp:138
int Dimension() const
Dimension of the reference space used within the elements.
Definition mesh.hpp:1216
void GetElementTransformation(int i, IsoparametricTransformation *ElTr) const
Builds the transformation defining the i-th element in ElTr. ElTr must be allocated in advance and wi...
Definition mesh.cpp:357
int SpaceDimension() const
Dimension of the physical space containing the mesh.
Definition mesh.hpp:1219
int MeshGenerator() const
Get the mesh generator/type.
Definition mesh.hpp:1243
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.
Definition mesh.cpp:6484
Geometry::Type GetElementBaseGeometry(int i) const
Definition mesh.hpp:1455
virtual void SetAttributes()
Determine the sets of unique attribute values in domain and boundary elements.
Definition mesh.cpp:1819
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...
Definition optparser.hpp:82
bool Good() const
Return true if the command line options were parsed successfully.
IntegrationRule RefPts
Definition geom.hpp:321
Vector data type.
Definition vector.hpp:82
int Size() const
Returns the size of the vector.
Definition vector.hpp:226
real_t Normlp(real_t p) const
Returns the l_p norm of the vector.
Definition vector.cpp:998
int dim
Definition ex24.cpp:53
int main()
GeometryRefiner GlobGeometryRefiner
Definition geom.cpp:2014
float real_t
Definition config.hpp:43
const char vishost[]
real_t p(const Vector &x, real_t t)
int material(Vector &x, Vector &xmin, Vector &xmax)
Definition shaper.cpp:53