MFEM  v4.6.0
Finite element discretization library
pml.hpp
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 #include "mfem.hpp"
13 
14 namespace mfem
15 {
16 
17 /// Class for setting up a simple Cartesian PML region
19 {
20 private:
21  Mesh *mesh;
22 
23  /** 2D array of size (dim,2) representing the length of the PML
24  region in each direction */
25  Array2D<double> length;
26 
27  /// 2D array of size (dim,2) representing the Computational Domain Boundary
28  Array2D<double> comp_dom_bdr;
29 
30  /// 2D array of size (dim,2) representing the Domain Boundary
31  Array2D<double> dom_bdr;
32 
33  /** Integer Array identifying elements in the pml
34  0: in the pml, 1: not in the pml */
35  Array<int> elems;
36 
37  /// Compute Domain and Computational Domain Boundaries
38  void SetBoundaries();
39 
40 public:
41  /** Constructor of the PML region using the mesh @a mesh_ and
42  the 2D array of size (dim,2) @a length_ which represents the
43  length of the PML in each direction. */
44  CartesianPML(Mesh *mesh_, const Array2D<double> &length_);
45 
46  int dim;
47  double omega;
48  // Default values for Maxwell
49  double epsilon = 1.0;
50  double mu = 1.0;
51  /// Return Computational Domain Boundary
52  const Array2D<double> & GetCompDomainBdr() {return comp_dom_bdr;}
53 
54  /// Return Domain Boundary
55  const Array2D<double> & GetDomainBdr() {return dom_bdr;}
56 
57  /// Return Marker list for elements
58  const Array<int> & GetMarkedPMLElements() {return elems;}
59 
60  /// Mark element in the PML region
61  void SetAttributes(Mesh *mesh_, Array<int> * attrNonPML = nullptr,
62  Array<int> * attrPML = nullptr);
63 
64  void SetOmega(double omega_) {omega = omega_;}
65  void SetEpsilonAndMu(double epsilon_, double mu_)
66  {
67  epsilon = epsilon_;
68  mu = mu_;
69  }
70  /// PML complex stretching function
71  void StretchFunction(const Vector &x, std::vector<std::complex<double>> &dxs);
72 };
73 
74 
76 {
77 private:
78  CartesianPML * pml = nullptr;
79  double (*Function)(const Vector &, CartesianPML * );
80 public:
81  PmlCoefficient(double (*F)(const Vector &, CartesianPML *), CartesianPML * pml_)
82  : pml(pml_), Function(F)
83  {}
84  virtual double Eval(ElementTransformation &T, const IntegrationPoint &ip)
85  {
86  double x[3];
87  Vector transip(x, 3);
88  T.Transform(ip, transip);
89  return ((*Function)(transip, pml));
90  }
91 };
92 
94 {
95 private:
96  CartesianPML * pml = nullptr;
97  void (*Function)(const Vector &, CartesianPML *, DenseMatrix &);
98 public:
99  PmlMatrixCoefficient(int dim, void(*F)(const Vector &, CartesianPML *,
100  DenseMatrix &),
101  CartesianPML * pml_)
102  : MatrixCoefficient(dim), pml(pml_), Function(F)
103  {}
105  const IntegrationPoint &ip)
106  {
107  double x[3];
108  Vector transip(x, 3);
109  T.Transform(ip, transip);
110  K.SetSize(height, width);
111  (*Function)(transip, pml, K);
112  }
113 };
114 
115 /// PML stretching functions: See https://doi.org/10.1006/jcph.1994.1159
116 // Helmholtz
117 double detJ_r_function(const Vector & x, CartesianPML * pml);
118 double detJ_i_function(const Vector & x, CartesianPML * pml);
119 double abs_detJ_2_function(const Vector & x, CartesianPML * pml);
120 
121 void Jt_J_detJinv_r_function(const Vector & x, CartesianPML * pml,
122  DenseMatrix & M);
123 void Jt_J_detJinv_i_function(const Vector & x, CartesianPML * pml,
124  DenseMatrix & M);
125 void abs_Jt_J_detJinv_2_function(const Vector & x, CartesianPML * pml,
126  DenseMatrix & M);
127 
128 // Maxwell
129 // |J| J^-1 J^-T
130 void detJ_Jt_J_inv_r_function(const Vector &x, CartesianPML * pml,
131  DenseMatrix &M);
132 void detJ_Jt_J_inv_i_function(const Vector &x, CartesianPML * pml,
133  DenseMatrix &M);
134 void abs_detJ_Jt_J_inv_2_function(const Vector &x, CartesianPML * pml,
135  DenseMatrix &M);
136 
137 } // namespace mfem
const Array2D< double > & GetDomainBdr()
Return Domain Boundary.
Definition: pml.hpp:55
const Array2D< double > & GetCompDomainBdr()
Return Computational Domain Boundary.
Definition: pml.hpp:52
Data type dense matrix using column-major storage.
Definition: densemat.hpp:23
void Jt_J_detJinv_i_function(const Vector &x, CartesianPML *pml, DenseMatrix &M)
Definition: pml.cpp:203
void SetOmega(double omega_)
Definition: pml.hpp:64
void abs_detJ_Jt_J_inv_2_function(const Vector &x, CartesianPML *pml, DenseMatrix &M)
Definition: pml.cpp:272
void detJ_Jt_J_inv_r_function(const Vector &x, CartesianPML *pml, DenseMatrix &M)
Definition: pml.cpp:238
double detJ_i_function(const Vector &x, CartesianPML *pml)
Definition: pml.cpp:166
CartesianPML(Mesh *mesh_, const Array2D< double > &length_)
Definition: pml.cpp:17
virtual double Eval(ElementTransformation &T, const IntegrationPoint &ip)
Evaluate the coefficient in the element described by T at the point ip.
Definition: pml.hpp:84
PmlCoefficient(double(*F)(const Vector &, CartesianPML *), CartesianPML *pml_)
Definition: pml.hpp:81
void abs_Jt_J_detJinv_2_function(const Vector &x, CartesianPML *pml, DenseMatrix &M)
Definition: pml.cpp:219
void SetAttributes(Mesh *mesh_, Array< int > *attrNonPML=nullptr, Array< int > *attrPML=nullptr)
Mark element in the PML region.
Definition: pml.cpp:68
double abs_detJ_2_function(const Vector &x, CartesianPML *pml)
Definition: pml.cpp:176
void SetEpsilonAndMu(double epsilon_, double mu_)
Definition: pml.hpp:65
void detJ_Jt_J_inv_i_function(const Vector &x, CartesianPML *pml, DenseMatrix &M)
Definition: pml.cpp:255
Base class Coefficients that optionally depend on space and time. These are used by the BilinearFormI...
Definition: coefficient.hpp:41
Base class for Matrix Coefficients that optionally depend on time and space.
double detJ_r_function(const Vector &x, CartesianPML *pml)
PML stretching functions: See https://doi.org/10.1006/jcph.1994.1159.
Definition: pml.cpp:156
Class for setting up a simple Cartesian PML region.
Definition: pml.hpp:18
double omega
Definition: pml.hpp:47
Class for integration point with weight.
Definition: intrules.hpp:31
int dim
Definition: ex24.cpp:53
const Array< int > & GetMarkedPMLElements()
Return Marker list for elements.
Definition: pml.hpp:58
void Jt_J_detJinv_r_function(const Vector &x, CartesianPML *pml, DenseMatrix &M)
Definition: pml.cpp:187
double epsilon
Definition: pml.hpp:49
Vector data type.
Definition: vector.hpp:58
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 ...
Definition: pml.hpp:104
void SetSize(int s)
Change the size of the DenseMatrix to s x s.
Definition: densemat.hpp:105
PmlMatrixCoefficient(int dim, void(*F)(const Vector &, CartesianPML *, DenseMatrix &), CartesianPML *pml_)
Definition: pml.hpp:99
void StretchFunction(const Vector &x, std::vector< std::complex< double >> &dxs)
PML complex stretching function.
Definition: pml.cpp:128