MFEM v4.7.0
Finite element discretization library
Loading...
Searching...
No Matches
lsf_integral.cpp
Go to the documentation of this file.
1// Copyright (c) 2010-2024, 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// Miniapp: Integral of implicit domains defined by a level set function
14// ---------------------------------------------------------------------
15//
16// The miniapp demonstrates an interface to the Algoim library for computing
17// volumetric and surface integrals over domains and surfaces defined
18// implicitly by a level set function. The miniapp requires MFEM to be built
19// with Blitz and Algoim libraries (see INSTALL).
20//
21// Compile with: make lsf_integral
22//
23// Sample runs:
24//
25// Evaluates surface and volumetric integral for a circle with radius 1
26// lsf_integral -m ../../data/star-q3.mesh
27//
28// Evaluates surface and volumetric integral for a level set defined
29// by y=0.5-(0.1*sin(3*pi*x+pi/2))
30// lsf_integral -ls 2 -m ../../data/inline-quad.mesh -rs 2 -o 3 -ao 3
31
32#include "mfem.hpp"
33#include "integ_algoim.hpp"
34
35using namespace mfem;
36using namespace std;
37
38// Level set function for sphere in 3D and circle in 2D
40{
41 real_t r2= x*x;
42 return -sqrt(r2)+1.0;//the radius is 1.0
43}
44
45// Level set function for a sinusoidal wave.
46// Resulting zero isocontour is at y=0.5-(0.1*sin(3*pi*x+pi/2))
48{
49 real_t a1 = 20., a2 = 2., a3 = 3.;
50 return tanh(a1*(x(1)-0.5) + a2*sin(a3*(x(0)-0.5)*M_PI));
51}
52
53int main(int argc, char *argv[])
54{
55 // Parse command-line options
56 const char *mesh_file = "../../data/star-q3.mesh";
57 int ser_ref_levels = 1;
58 int order = 2;
59 int iorder = 2; // MFEM integration points
60 int aorder = 2; // Algoim integration points
61 bool visualization = true;
62 int print_level = 0;
63 int ls_type = 1;
64
65 OptionsParser args(argc, argv);
66 args.AddOption(&mesh_file, "-m", "--mesh", "Mesh file to use.");
67 args.AddOption(&ser_ref_levels,
68 "-rs",
69 "--refine-serial",
70 "Number of times to refine the mesh uniformly in serial.");
71 args.AddOption(&order,
72 "-o",
73 "--order",
74 "Order (degree) of the finite elements.");
75 args.AddOption(&iorder,
76 "-io",
77 "--iorder",
78 "MFEM Integration order.");
79 args.AddOption(&aorder,
80 "-ao",
81 "--aorder",
82 "Algoim Integration order.");
83 args.AddOption(&ls_type,
84 "-ls",
85 "--ls-type",
86 "Level set type: 1: circle, 2 sinusoidal wave");
87 args.AddOption(&visualization,
88 "-vis",
89 "--visualization",
90 "-no-vis",
91 "--no-visualization",
92 "Enable or disable GLVis visualization.");
93 args.AddOption((&print_level), "-prt", "--print-level", "Print level.");
94 args.Parse();
95 if (!args.Good())
96 {
97 args.PrintUsage(std::cout);
98 return 1;
99 }
100 args.PrintOptions(std::cout);
101
102 // Read the (serial) mesh from the given mesh file.
103 Mesh *mesh = new Mesh(mesh_file, 1, 1);
104 int dim = mesh->Dimension();
105 // Refine the mesh in serial to increase the resolution. In this example
106 // we do 'ser_ref_levels' of uniform refinement, where 'ser_ref_levels' is
107 // a command-line parameter.
108 for (int lev = 0; lev < ser_ref_levels; lev++)
109 {
110 mesh->UniformRefinement();
111 }
112
113 // Define the finite element space for the level-set function.
114 H1_FECollection fec(order, dim);
115 FiniteElementSpace fespace(mesh, &fec, 1, Ordering::byVDIM);
116 int glob_size = fespace.GetTrueVSize();
117 std::cout << "Number of finite element unknowns: " << glob_size << std::endl;
118
119 // Define the level set grid function
120 GridFunction x(&fespace);
121
122 // Define the level set coefficient
123 Coefficient *ls_coeff = nullptr;
124 if (ls_type == 1)
125 {
126 ls_coeff=new FunctionCoefficient(sphere_ls);
127 }
128 else if (ls_type == 2)
129 {
131 }
132 else
133 {
134 MFEM_ABORT("Level set coefficient not defined");
135 }
136
137 // Project the coefficient onto the LS grid function
138 x.ProjectCoefficient(*ls_coeff);
139
140 if (visualization)
141 {
142 osockstream sock(19916, "localhost");
143 sock << "solution\n";
144 mesh->Print(sock);
145 x.Save(sock);
146 sock.send();
147 sock << "window_title 'Level set'\n"
148 << "window_geometry "
149 << 1200 << " " << 0 << " " << 600 << " " << 600 << "\n"
150 << "keys jRmmclA" << endl;
151 }
152
153 // Exact volume and area
154 real_t exact_volume = -10, exact_area = -10;
155 if (ls_type == 1)
156 {
157 if (strncmp(mesh_file,"../../data/star-q3.mesh",100) == 0)
158 {
159 exact_volume = M_PI;
160 exact_area = M_PI*2;
161 }
162 else if (strncmp(mesh_file, "../../data/inline-quad.mesh",100) == 0)
163 {
164 exact_volume = M_PI/4;
165 exact_area = M_PI/2;
166 }
167 }
168 else if (ls_type == 2)
169 {
170 if (strncmp(mesh_file, "../../data/inline-quad.mesh",100) == 0)
171 {
172 exact_volume = 0.5;
173 exact_area = 1.194452300992437;
174 }
175 }
176 (void)(&exact_area); // suppress a warning
177
179 const IntegrationRule* ir=nullptr;
180
181 // Integration with Algoim
182 real_t vol=0.0;
183
184#ifdef MFEM_USE_ALGOIM
185 real_t area=0.0;
186 DenseMatrix bmat; // gradients of the shape functions in isoparametric space
187 DenseMatrix pmat; // gradients of the shape functions in physical space
188 Vector inormal; // normal to the level set in isoparametric space
189 Vector tnormal; // normal to the level set in physical space
190 Vector lsfun; // level set function restricted to an element
191 DofTransformation *doftrans;
192 Array<int> vdofs;
193 for (int i=0; i<fespace.GetNE(); i++)
194 {
195 const FiniteElement* el=fespace.GetFE(i);
196
197 // get the element transformation
198 trans = fespace.GetElementTransformation(i);
199
200 // extract the element vector from the level-set
201 doftrans = fespace.GetElementVDofs(i,vdofs);
202 x.GetSubVector(vdofs, lsfun);
203
204 // construct Algoim integration object
205 AlgoimIntegrationRule air(aorder, *el, *trans, lsfun);
206
207 // compute the volume contribution from the element
208 ir = air.GetVolumeIntegrationRule();
209 for (int j = 0; j < ir->GetNPoints(); j++)
210 {
211 const IntegrationPoint &ip = ir->IntPoint(j);
212 trans->SetIntPoint(&ip);
213 vol += ip.weight * trans->Weight();
214 }
215
216 // compute the perimeter/area contribution from the element
217 bmat.SetSize(el->GetDof(),el->GetDim());
218 pmat.SetSize(el->GetDof(),el->GetDim());
219 inormal.SetSize(el->GetDim());
220 tnormal.SetSize(el->GetDim());
221
222 ir = air.GetSurfaceIntegrationRule();
223 for (int j = 0; j < ir->GetNPoints(); j++)
224 {
225 const IntegrationPoint &ip = ir->IntPoint(j);
226 trans->SetIntPoint(&ip);
227
228 el->CalcDShape(ip,bmat);
229 Mult(bmat, trans->AdjugateJacobian(), pmat);
230 // compute the normal to the LS in isoparametric space
231 bmat.MultTranspose(lsfun,inormal);
232 // compute the normal to the LS in physical space
233 pmat.MultTranspose(lsfun,tnormal);
234 area += ip.weight * tnormal.Norml2() / inormal.Norml2();
235 }
236 }
237
238 if (exact_volume > 0)
239 {
240 std::cout<<"Algoim Volume="<<vol<<" Error="<<vol-exact_volume<<std::endl;
241 std::cout<<"Algoim Area="<<area<<" Error="<<area-exact_area<<std::endl;
242 }
243 else
244 {
245 std::cout<<"Algoim Volume="<<vol<<std::endl;
246 std::cout<<"Algoim Area="<<area<<std::endl;
247 }
248#endif
249
250 // Perform standard MFEM integration
251 vol=0.0;
252 for (int i=0; i<fespace.GetNE(); i++)
253 {
254 const FiniteElement* el=fespace.GetFE(i);
255 // get the element transformation
256 trans = fespace.GetElementTransformation(i);
257
258 // compute the volume contribution from the element
259 ir=&IntRules.Get(el->GetGeomType(), iorder);
260 for (int j = 0; j < ir->GetNPoints(); j++)
261 {
262 const IntegrationPoint &ip = ir->IntPoint(j);
263 trans->SetIntPoint(&ip);
264 real_t vlsf=x.GetValue(*trans,ip);
265 if (vlsf>=0.0)
266 {
267 vol += ip.weight * trans->Weight();
268 }
269 }
270 }
271
272 if (exact_volume > 0.0)
273 {
274 std::cout<<"MFEM Volume="<<vol<<" Error="<<vol-exact_volume<<std::endl;
275 }
276 else
277 {
278 std::cout<<"MFEM Volume="<<vol<<std::endl;
279 }
280
281 ParaViewDataCollection dacol("ParaViewLSF", mesh);
282 dacol.SetLevelsOfDetail(order);
283 dacol.RegisterField("LSF",&x);
284 dacol.SetTime(1.0);
285 dacol.SetCycle(1);
286 dacol.Save();
287
288 delete ls_coeff;
289 delete mesh;
290
291 return 0;
292}
const IntegrationRule * GetVolumeIntegrationRule()
const IntegrationRule * GetSurfaceIntegrationRule()
Base class Coefficients that optionally depend on space and time. These are used by the BilinearFormI...
virtual void RegisterField(const std::string &field_name, GridFunction *gf)
Add a grid function to the collection.
void SetCycle(int c)
Set time cycle (for time-dependent simulations)
void SetTime(real_t t)
Set physical time (for time-dependent simulations)
Data type dense matrix using column-major storage.
Definition densemat.hpp:24
void MultTranspose(const real_t *x, real_t *y) const
Multiply a vector with the transpose matrix.
Definition densemat.cpp:262
void SetSize(int s)
Change the size of the DenseMatrix to s x s.
Definition densemat.hpp:105
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
Definition fespace.hpp:220
virtual int GetTrueVSize() const
Return the number of vector true (conforming) dofs.
Definition fespace.hpp:716
ElementTransformation * GetElementTransformation(int i) const
Returns ElementTransformation for the i-th element.
Definition fespace.hpp:773
DofTransformation * GetElementVDofs(int i, Array< int > &vdofs) const
Returns indices of degrees of freedom for the i'th element. The returned indices are offsets into an ...
Definition fespace.cpp:287
virtual const FiniteElement * GetFE(int i) const
Returns pointer to the FiniteElement in the FiniteElementCollection associated with i'th element in t...
Definition fespace.cpp:3168
int GetNE() const
Returns number of elements in the mesh.
Definition fespace.hpp:740
Abstract class for all finite elements.
Definition fe_base.hpp:239
int GetDim() const
Returns the reference space dimension for the finite element.
Definition fe_base.hpp:316
virtual void CalcDShape(const IntegrationPoint &ip, DenseMatrix &dshape) const =0
Evaluate the gradients of all shape functions of a scalar finite element in reference space at the gi...
Geometry::Type GetGeomType() const
Returns the Geometry::Type of the reference element.
Definition fe_base.hpp:326
int GetDof() const
Returns the number of degrees of freedom in the finite element.
Definition fe_base.hpp:329
A general function coefficient.
Class for grid function - Vector with associated FE space.
Definition gridfunc.hpp:31
virtual real_t GetValue(int i, const IntegrationPoint &ip, int vdim=1) const
Definition gridfunc.cpp:449
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...
Arbitrary order H1-conforming (continuous) finite elements.
Definition fe_coll.hpp:260
Class for integration point with weight.
Definition intrules.hpp:35
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
const IntegrationRule & Get(int GeomType, int Order)
Returns an integration rule for given GeomType and Order.
Mesh data type.
Definition mesh.hpp:56
virtual void Print(std::ostream &os=mfem::out, const std::string &comments="") const
Definition mesh.hpp:2288
int Dimension() const
Dimension of the reference space used within the elements.
Definition mesh.hpp:1160
void UniformRefinement(int i, const DSTable &, int *, int *, int *)
Definition mesh.cpp:10970
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.
Helper class for ParaView visualization data.
void SetLevelsOfDetail(int levels_of_detail_)
virtual void Save() override
Vector data type.
Definition vector.hpp:80
real_t Norml2() const
Returns the l2 norm of the vector.
Definition vector.cpp:832
void SetSize(int s)
Resize the vector to size s.
Definition vector.hpp:538
void GetSubVector(const Array< int > &dofs, Vector &elemvect) const
Extract entries listed in dofs to the output Vector elemvect.
Definition vector.cpp:578
int dim
Definition ex24.cpp:53
void trans(const Vector &u, Vector &x)
Definition ex27.cpp:412
int main()
real_t sinusoidal_ls(const Vector &x)
real_t sphere_ls(const Vector &x)
void Mult(const Table &A, const Table &B, Table &C)
C = A * B (as boolean matrices)
Definition table.cpp:476
float real_t
Definition config.hpp:43
IntegrationRules IntRules(0, Quadrature1D::GaussLegendre)
A global object with all integration rules (defined in intrules.cpp)
Definition intrules.hpp:486