MFEM  v4.2.0
Finite element discretization library
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
mesh-optimizer.hpp
Go to the documentation of this file.
1 // Copyright (c) 2010-2020, 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 // MFEM Mesh Optimizer Miniapp - Serial/Parallel Shared Code
13 
14 #include "mfem.hpp"
15 #include <fstream>
16 #include <iostream>
17 
18 using namespace mfem;
19 using namespace std;
20 
21 double discrete_size_2d(const Vector &x)
22 {
23  int opt = 2;
24  const double small = 0.001, big = 0.01;
25  double val = 0.;
26 
27  if (opt == 1) // sine wave.
28  {
29  const double X = x(0), Y = x(1);
30  val = std::tanh((10*(Y-0.5) + std::sin(4.0*M_PI*X)) + 1) -
31  std::tanh((10*(Y-0.5) + std::sin(4.0*M_PI*X)) - 1);
32  }
33  else if (opt == 2) // semi-circle
34  {
35  const double xc = x(0) - 0.0, yc = x(1) - 0.5;
36  const double r = sqrt(xc*xc + yc*yc);
37  double r1 = 0.45; double r2 = 0.55; double sf=30.0;
38  val = 0.5*(1+std::tanh(sf*(r-r1))) - 0.5*(1+std::tanh(sf*(r-r2)));
39  }
40 
41  val = std::max(0.,val);
42  val = std::min(1.,val);
43 
44  return val * small + (1.0 - val) * big;
45 }
46 
47 double material_indicator_2d(const Vector &x)
48 {
49  double xc = x(0)-0.5, yc = x(1)-0.5;
50  double th = 22.5*M_PI/180.;
51  double xn = cos(th)*xc + sin(th)*yc;
52  double yn = -sin(th)*xc + cos(th)*yc;
53  double th2 = (th > 45.*M_PI/180) ? M_PI/2 - th : th;
54  double stretch = 1/cos(th2);
55  xc = xn/stretch; yc = yn/stretch;
56  double tfac = 20;
57  double s1 = 3;
58  double s2 = 3;
59  double wgt = std::tanh((tfac*(yc) + s2*std::sin(s1*M_PI*xc)) + 1);
60  if (wgt > 1) { wgt = 1; }
61  if (wgt < 0) { wgt = 0; }
62  return wgt;
63 }
64 
65 double discrete_ori_2d(const Vector &x)
66 {
67  return M_PI * x(1) * (1.0 - x(1)) * cos(2 * M_PI * x(0));
68 }
69 
70 double discrete_aspr_2d(const Vector &x)
71 {
72  double xc = x(0)-0.5, yc = x(1)-0.5;
73  double th = 22.5*M_PI/180.;
74  double xn = cos(th)*xc + sin(th)*yc;
75  double yn = -sin(th)*xc + cos(th)*yc;
76  //double th2 = (th > 45.*M_PI/180) ? M_PI/2 - th : th;
77  //double stretch = 1/cos(th2);
78  xc = xn; yc = yn;
79 
80  double tfac = 20;
81  double s1 = 3;
82  double s2 = 2;
83  double wgt = std::tanh((tfac*(yc) + s2*std::sin(s1*M_PI*xc)) + 1)
84  - std::tanh((tfac*(yc) + s2*std::sin(s1*M_PI*xc)) - 1);
85  if (wgt > 1) { wgt = 1; }
86  if (wgt < 0) { wgt = 0; }
87  return 0.1 + 1*(1-wgt)*(1-wgt);
88 }
89 
90 void discrete_aspr_3d(const Vector &x, Vector &v)
91 {
92  int dim = x.Size();
93  v.SetSize(dim);
94  double l1, l2, l3;
95  l1 = 1.;
96  l2 = 1. + 5*x(1);
97  l3 = 1. + 10*x(2);
98  v[0] = l1/pow(l2*l3,0.5);
99  v[1] = l2/pow(l1*l3,0.5);
100  v[2] = l3/pow(l2*l1,0.5);
101 }
102 
104 {
105 private:
106  int metric;
107 
108 public:
109  HessianCoefficient(int dim, int metric_id)
110  : TMOPMatrixCoefficient(dim), metric(metric_id) { }
111 
113  const IntegrationPoint &ip)
114  {
115  Vector pos(3);
116  T.Transform(ip, pos);
117  if (metric != 14 && metric != 85)
118  {
119  const double xc = pos(0) - 0.5, yc = pos(1) - 0.5;
120  const double r = sqrt(xc*xc + yc*yc);
121  double r1 = 0.15; double r2 = 0.35; double sf=30.0;
122  const double eps = 0.5;
123 
124  const double tan1 = std::tanh(sf*(r-r1)),
125  tan2 = std::tanh(sf*(r-r2));
126 
127  K(0, 0) = eps + 1.0 * (tan1 - tan2);
128  K(0, 1) = 0.0;
129  K(1, 0) = 0.0;
130  K(1, 1) = 1.0;
131  }
132  else if (metric == 14) // Size + Alignment
133  {
134  const double xc = pos(0), yc = pos(1);
135  double theta = M_PI * yc * (1.0 - yc) * cos(2 * M_PI * xc);
136  double alpha_bar = 0.1;
137 
138  K(0, 0) = cos(theta);
139  K(1, 0) = sin(theta);
140  K(0, 1) = -sin(theta);
141  K(1, 1) = cos(theta);
142 
143  K *= alpha_bar;
144  }
145  else if (metric == 85) // Shape + Alignment
146  {
147  Vector x = pos;
148  double xc = x(0)-0.5, yc = x(1)-0.5;
149  double th = 22.5*M_PI/180.;
150  double xn = cos(th)*xc + sin(th)*yc;
151  double yn = -sin(th)*xc + cos(th)*yc;
152  xc = xn; yc=yn;
153 
154  double tfac = 20;
155  double s1 = 3;
156  double s2 = 2;
157  double wgt = std::tanh((tfac*(yc) + s2*std::sin(s1*M_PI*xc)) + 1)
158  - std::tanh((tfac*(yc) + s2*std::sin(s1*M_PI*xc)) - 1);
159  if (wgt > 1) { wgt = 1; }
160  if (wgt < 0) { wgt = 0; }
161 
162  xc = pos(0), yc = pos(1);
163  double theta = M_PI * (yc) * (1.0 - yc) * cos(2 * M_PI * xc);
164 
165  K(0, 0) = cos(theta);
166  K(1, 0) = sin(theta);
167  K(0, 1) = -sin(theta);
168  K(1, 1) = cos(theta);
169 
170  double asp_ratio_tar = 0.1 + 1*(1-wgt)*(1-wgt);
171 
172  K(0, 0) *= 1/pow(asp_ratio_tar,0.5);
173  K(1, 0) *= 1/pow(asp_ratio_tar,0.5);
174  K(0, 1) *= pow(asp_ratio_tar,0.5);
175  K(1, 1) *= pow(asp_ratio_tar,0.5);
176  }
177  }
178 
180  const IntegrationPoint &ip, int comp)
181  {
182  Vector pos(3);
183  T.Transform(ip, pos);
184  K = 0.;
185  if (metric != 14 && metric != 85)
186  {
187  const double xc = pos(0) - 0.5, yc = pos(1) - 0.5;
188  const double r = sqrt(xc*xc + yc*yc);
189  double r1 = 0.15; double r2 = 0.35; double sf=30.0;
190 
191  const double tan1 = std::tanh(sf*(r-r1)),
192  tan2 = std::tanh(sf*(r-r2));
193  double tan1d = 0., tan2d = 0.;
194  if (r > 0.001)
195  {
196  tan1d = (1.-tan1*tan1)*(sf)/r,
197  tan2d = (1.-tan2*tan2)*(sf)/r;
198  }
199 
200  K(0, 1) = 0.0;
201  K(1, 0) = 0.0;
202  K(1, 1) = 1.0;
203  if (comp == 0) { K(0, 0) = tan1d*xc - tan2d*xc; }
204  else if (comp == 1) { K(0, 0) = tan1d*yc - tan2d*yc; }
205  }
206  }
207 };
208 
209 // Additional IntegrationRules that can be used with the --quad-type option.
212 
213 // Defined with respect to the icf mesh.
214 double weight_fun(const Vector &x)
215 {
216  const double r = sqrt(x(0)*x(0) + x(1)*x(1) + 1e-12);
217  const double den = 0.002;
218  double l2 = 0.2 + 0.5*std::tanh((r-0.16)/den) - 0.5*std::tanh((r-0.17)/den)
219  + 0.5*std::tanh((r-0.23)/den) - 0.5*std::tanh((r-0.24)/den);
220  return l2;
221 }
222 
223 // Used for the adaptive limiting examples.
224 double adapt_lim_fun(const Vector &x)
225 {
226  const double xc = x(0) - 0.1, yc = x(1) - 0.2;
227  const double r = sqrt(xc*xc + yc*yc);
228  double r1 = 0.45; double r2 = 0.55; double sf=30.0;
229  double val = 0.5*(1+std::tanh(sf*(r-r1))) - 0.5*(1+std::tanh(sf*(r-r2)));
230 
231  val = std::max(0.,val);
232  val = std::min(1.,val);
233  return val;
234 }
235 
236 void DiffuseField(GridFunction &field, int smooth_steps)
237 {
238  //Setup the Laplacian operator
239  BilinearForm *Lap = new BilinearForm(field.FESpace());
241  Lap->Assemble();
242  Lap->Finalize();
243 
244  //Setup the smoothing operator
245  DSmoother *S = new DSmoother(0,1.0,smooth_steps);
246  S->iterative_mode = true;
247  S->SetOperator(Lap->SpMat());
248 
249  Vector tmp(field.Size());
250  tmp = 0.0;
251  S->Mult(tmp, field);
252 
253  delete S;
254  delete Lap;
255 }
256 
257 #ifdef MFEM_USE_MPI
258 void DiffuseField(ParGridFunction &field, int smooth_steps)
259 {
260  //Setup the Laplacian operator
261  ParBilinearForm *Lap = new ParBilinearForm(field.ParFESpace());
263  Lap->Assemble();
264  Lap->Finalize();
265  HypreParMatrix *A = Lap->ParallelAssemble();
266 
267  HypreSmoother *S = new HypreSmoother(*A,0,smooth_steps);
268  S->iterative_mode = true;
269 
270  Vector tmp(A->Width());
271  field.SetTrueVector();
272  Vector fieldtrue = field.GetTrueVector();
273  tmp = 0.0;
274  S->Mult(tmp, fieldtrue);
275 
276  field.SetFromTrueDofs(fieldtrue);
277 
278  delete S;
279  delete Lap;
280 }
281 #endif
void discrete_aspr_3d(const Vector &x, Vector &v)
virtual void EvalGrad(DenseMatrix &K, ElementTransformation &T, const IntegrationPoint &ip, int comp)
Evaluate the derivative of the matrix coefficient with respect to comp in the element described by T ...
Class for grid function - Vector with associated FE space.
Definition: gridfunc.hpp:30
double discrete_size_2d(const Vector &x)
Data type for scaled Jacobi-type smoother of sparse matrix.
IntegrationRules IntRulesLo(0, Quadrature1D::GaussLobatto)
void Assemble(int skip_zeros=1)
Assembles the form i.e. sums over all domain/bdr integrators.
void SetSize(int s)
Resize the vector to size s.
Definition: vector.hpp:459
void DiffuseField(GridFunction &field, int smooth_steps)
int Width() const
Get the width (size of input) of the Operator. Synonym with NumCols().
Definition: operator.hpp:71
virtual void SetOperator(const Operator &a)
Set/update the solver for the given operator.
aka closed Newton-Cotes
Definition: intrules.hpp:296
Container class for integration rules.
Definition: intrules.hpp:309
Data type dense matrix using column-major storage.
Definition: densemat.hpp:23
int Size() const
Returns the size of the vector.
Definition: vector.hpp:160
bool iterative_mode
If true, use the second argument of Mult() as an initial guess.
Definition: operator.hpp:638
HypreParMatrix * ParallelAssemble()
Returns the matrix assembled on the true dofs, i.e. P^t A P.
virtual void Mult(const HypreParVector &b, HypreParVector &x) const
Relax the linear system Ax=b.
Definition: hypre.cpp:2369
void SetTrueVector()
Shortcut for calling GetTrueDofs() with GetTrueVector() as argument.
Definition: gridfunc.hpp:136
virtual void SetFromTrueDofs(const Vector &tv)
Set the GridFunction from the given true-dof vector.
Definition: pgridfunc.hpp:164
void Assemble(int skip_zeros=1)
Assemble the local matrix.
double discrete_ori_2d(const Vector &x)
Parallel smoothers in hypre.
Definition: hypre.hpp:592
FiniteElementSpace * FESpace()
Definition: gridfunc.hpp:595
const Vector & GetTrueVector() const
Read only access to the (optional) internal true-dof Vector.
Definition: gridfunc.hpp:124
virtual void Eval(DenseMatrix &K, ElementTransformation &T, const IntegrationPoint &ip)
Evaluate the matrix coefficient in the element described by T at the point ip, storing the result in ...
virtual void Mult(const Vector &x, Vector &y) const
Matrix vector multiplication with Jacobi smoother.
HessianCoefficient(int dim, int metric_id)
double weight_fun(const Vector &x)
double adapt_lim_fun(const Vector &x)
double material_indicator_2d(const Vector &x)
A &quot;square matrix&quot; operator for the associated FE space and BLFIntegrators The sum of all the BLFInteg...
double discrete_aspr_2d(const Vector &x)
Class for integration point with weight.
Definition: intrules.hpp:25
virtual void Finalize(int skip_zeros=1)
Finalizes the matrix initialization.
int dim
Definition: ex24.cpp:53
IntegrationRules IntRulesCU(0, Quadrature1D::ClosedUniform)
void AddDomainIntegrator(BilinearFormIntegrator *bfi)
Adds new Domain Integrator. Assumes ownership of bfi.
Class for parallel bilinear form.
Vector data type.
Definition: vector.hpp:51
virtual void Transform(const IntegrationPoint &, Vector &)=0
Transform integration point from reference coordinates to physical coordinates and store them in the ...
const SparseMatrix & SpMat() const
Returns a const reference to the sparse matrix.
Class for parallel grid function.
Definition: pgridfunc.hpp:32
Wrapper for hypre&#39;s ParCSR matrix class.
Definition: hypre.hpp:181
ParFiniteElementSpace * ParFESpace() const
Definition: pgridfunc.hpp:108