Galerkin FEM for elliptic PDEs
You can not select more than 25 topics Topics must start with a letter or number, can include dashes ('-') and can be up to 35 characters long.

FEM4EllipticPDE.h 11KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318
  1. // ===========================================================================
  2. //
  3. // Filename: FEM4EllipticPDE.h
  4. //
  5. // Created: 08/16/12 18:19:41
  6. // Compiler: Tested with g++, icpc, and MSVC 2010
  7. //
  8. // Author: Trevor Irons (ti)
  9. //
  10. // Organisation: Colorado School of Mines (CSM)
  11. // United States Geological Survey (USGS)
  12. //
  13. // Email: tirons@mines.edu, tirons@usgs.gov
  14. //
  15. // This program is free software: you can redistribute it and/or modify
  16. // it under the terms of the GNU General Public License as published by
  17. // the Free Software Foundation, either version 3 of the License, or
  18. // (at your option) any later version.
  19. //
  20. // This program is distributed in the hope that it will be useful,
  21. // but WITHOUT ANY WARRANTY; without even the implied warranty of
  22. // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
  23. // GNU General Public License for more details.
  24. //
  25. // You should have received a copy of the GNU General Public License
  26. // along with this program. If not, see <http://www.gnu.org/licenses/>.
  27. //
  28. // ===========================================================================
  29. /**
  30. @file
  31. @author Trevor Irons
  32. @date 08/16/12
  33. @version $Rev$
  34. **/
  35. #ifndef FEM4ELLIPTICPDE_INC
  36. #define FEM4ELLIPTICPDE_INC
  37. #include "LemmaObject.h"
  38. #include "vtkImplicitFunction.h"
  39. #include "vtkDataSet.h"
  40. #include "vtkXMLDataSetWriter.h"
  41. #include "vtkSmartPointer.h"
  42. #include "vtkIdList.h"
  43. #include "vtkCell.h"
  44. #include "vtkDoubleArray.h"
  45. #include "vtkPointData.h"
  46. #include "vtkDataSetSurfaceFilter.h"
  47. #include "vtkCellArray.h"
  48. #include "vtkCellData.h"
  49. #include <vtkXMLUnstructuredGridReader.h>
  50. #include <vtkUnstructuredGrid.h>
  51. #include <Eigen/IterativeLinearSolvers>
  52. #include <Eigen/SparseLU>
  53. #include "DCSurvey.h"
  54. namespace Lemma {
  55. /** @defgroup FEM4EllipticPDE
  56. @brief Finite element solver for elliptic PDE's
  57. @details Uses Galerkin method.
  58. */
  59. /** \addtogroup FEM4EllipticPDE
  60. @{
  61. */
  62. // ===================================================================
  63. // Class: FEM4EllipticPDE
  64. /**
  65. @class FEM4EllipticPDE
  66. @brief Galerkin FEM solver for Elliptic PDE problems
  67. @details Solves problems of the type
  68. \f[ - \nabla \cdot ( \sigma(\mathbf{r}) \nabla u(\mathbf{r}) ) = g(\mathbf{r})
  69. \f]
  70. */
  71. class FEM4EllipticPDE : public LemmaObject {
  72. friend std::ostream &operator<<(std::ostream &stream,
  73. const FEM4EllipticPDE &ob);
  74. public:
  75. // ==================== LIFECYCLE =======================
  76. /** Returns a pointer to a new object of type FEM4EllipticPDE.
  77. * It allocates all necessary memory.
  78. */
  79. static std::shared_ptr<FEM4EllipticPDE> NewSP();
  80. /**
  81. * Default locked constructor.
  82. */
  83. explicit FEM4EllipticPDE (const ctor_key& key);
  84. /** Protected DeDerializing constructor, use factory DeSerialize method*/
  85. FEM4EllipticPDE (const YAML::Node& node, const ctor_key& key);
  86. /** Default protected constructor. */
  87. ~FEM4EllipticPDE ();
  88. /**
  89. * Uses YAML to serialize this object.
  90. * @return a YAML::Node
  91. */
  92. YAML::Node Serialize() const;
  93. /**
  94. * Constructs an object from a YAML::Node.
  95. */
  96. static std::shared_ptr<FEM4EllipticPDE> DeSerialize(const YAML::Node& node);
  97. // ==================== OPERATORS =======================
  98. // ==================== OPERATIONS =======================
  99. /** Sets the vtkImplictFunction that is used as sigma in
  100. \f$ - \nabla \cdot ( \sigma(\mathbf{r}) \nabla u(\mathbf{r}) ) = g(\mathbf{r})
  101. \f$. If this is not set, sigma evaluates as 1.
  102. */
  103. void SetSigmaFunction(vtkImplicitFunction* sigma);
  104. /** Sets the step to be used in applying boundary conditions.
  105. */
  106. void SetBoundaryStep(const Real& h);
  107. /** Sets the vtkImplictFunction that is used
  108. */
  109. void SetGFunction(vtkImplicitFunction* G);
  110. /** Uses a function pointer instead of a vtkImplicitFunction. For functions that can
  111. * be written in closed form, this can be much faster and more accurate then interpolating
  112. * off a grid.
  113. */
  114. void SetGFunction( Real (*pFcn3)(const Real&, const Real&, const Real&) );
  115. /* Sets the vtkDataSet that defines the calculation grid.
  116. */
  117. //void SetGrid(vtkDataSet* Grid);
  118. /** Sets the vtkDataSet that defines the calculation grid.
  119. */
  120. void SetGrid(vtkUnstructuredGrid* Grid);
  121. /**
  122. * Read grid in from VTK file
  123. */
  124. void SetVTUGridFile( const std::string& vtkGridFile );
  125. /**
  126. * Read grid in from VTK file
  127. */
  128. void SetVTUGridFile( const char* vtkGridFile );
  129. /** Sets up a DC problem with a survey
  130. * @param[in] ij is the injection index
  131. */
  132. void SetupDC(DCSurvey* Survey, const int& ij);
  133. /** Sets up the potential problem from the VTK file
  134. */
  135. void SetupPotential();
  136. /** Performs solution */
  137. void Solve(const std::string& fname);
  138. /* Performs solution */
  139. //void SolveOLD2(const std::string& fname);
  140. // ==================== ACCESS =======================
  141. // ==================== INQUIRY =======================
  142. /** Returns the name of the underlying class, similiar to Python's type */
  143. virtual std::string GetName() const {
  144. return this->CName;
  145. }
  146. protected:
  147. // ==================== LIFECYCLE =======================
  148. // ==================== OPERATIONS =======================
  149. /** Used internally to construct stiffness matrix, A
  150. */
  151. void ConstructAMatrix();
  152. /* Used internally to construct the load vector, g
  153. */
  154. //void ConstructLoadVector();
  155. // This function is copyright (C) 2012 Joseph Capriotti
  156. /**
  157. Returns a vtkIdList of points in member vtkGrid connected to point ido.
  158. @param[in] ido is the point of interest
  159. @param[in] grid is a pointer to a vtkDataSet
  160. */
  161. vtkSmartPointer<vtkIdList> GetConnectedPoints(const int& id0);
  162. /** Uses Simpson's rule to numerically integrate
  163. * (Shamelessly adapted from http://en.wikipedia.org/wiki/Simpson's_rule
  164. */
  165. Real CompositeSimpsons(vtkImplicitFunction* f, Real r0[3], Real r1[3], int numInt);
  166. /** Uses Simpson's rule to numerically integrate
  167. * (Shamelessly adapted from http://en.wikipedia.org/wiki/Simpson's_rule
  168. * this method is for a static f
  169. */
  170. Real CompositeSimpsons(const Real& f, Real r0[3], Real r1[3], int numInt);
  171. /**
  172. * Uses Simpson's rule to numerically integrate two functions in 1 dimension over the points
  173. * r0 and r1
  174. */
  175. Real CompositeSimpsons2(vtkImplicitFunction* f, Real r0[3], Real r1[3], int numInt);
  176. /**
  177. * Uses Simpson's rule to numerically integrate two functions in 1 dimension over the points
  178. * r0 and r1, uses the Hat function and the function pointer fPtr
  179. */
  180. Real CompositeSimpsons2( Real (*fPtr)(const Real&, const Real&, const Real&), Real r0[3], Real r1[3], int numInt);
  181. /**
  182. * Uses Simpson's rule to numerically integrate two functions in 1 dimension over the points
  183. * r0 and r1, uses the Hat function. Assumes a constand function value f
  184. */
  185. Real CompositeSimpsons2( const Real& f, Real r0[3], Real r1[3], int numInt);
  186. /**
  187. * Uses Simpson's rule to numerically integrate two functions in 1 dimension over the points
  188. * r0 and r1, uses the Hat function. Assumes a linear function from f0 to f1.
  189. */
  190. Real CompositeSimpsons2( const Real& f0, const Real& f1, Real r0[3], Real r1[3], int numInt);
  191. /**
  192. * Uses Simpson's rule to numerically integrate three functions in 1 dimension over the points
  193. * r0 and r1, uses two Hat functions. Assumes a linear function from f0 to f1.
  194. */
  195. Real CompositeSimpsons3( const Real& f0, const Real& f1, Real r0[3], Real r1[3], int numInt);
  196. /** Hat function for use in Galerkin FEM method, actually a 1/2 Hat.
  197. */
  198. //Real hat(Real r[3], Real ri[3], Real rm1[3] );
  199. Real Hat(const Vector3r &r, const Vector3r& ri, const Vector3r& rm1 );
  200. /** Computes distance between r0 and r1
  201. */
  202. Real dist(Real r0[3], Real r1[3]);
  203. /** Computes distance between r0 and r1
  204. */
  205. Real dist(const Vector3r& r0, const Vector3r& r1);
  206. /** Sets up the potential problem from the VTK file
  207. */
  208. void SetupLineSourcePotential();
  209. /** Sets up the potential problem from the VTK file
  210. */
  211. void SetupSurfaceSourcePotential();
  212. /** Sets up the potential problem from the VTK file
  213. */
  214. void SetupPointSourcePotential();
  215. /** Sets up the potential problem from the VTK file
  216. */
  217. void SetupVolumeSourcePotential();
  218. // ==================== DATA MEMBERS =========================
  219. /** The effective step at the boundary condition. Defaults to 1
  220. */
  221. Real BndryH;
  222. /** The sigma value at the boundary condition. Defaults to 1
  223. */
  224. Real BndrySigma;
  225. /** Implicit function defining the \f$\sigma\f$ term in
  226. \f$ - \nabla \cdot ( \sigma(\mathbf{r}) \nabla u(\mathbf{r}) ) = g(\mathbf{r})
  227. \f$.
  228. If set to 1, the Poisson equation is solved. Any type of vtkImplicitFunction may
  229. be used, including vtkImplicitDataSet.
  230. */
  231. vtkSmartPointer<vtkImplicitFunction> vtkSigma;
  232. /** Implicit function defining the \f$g\f$ term in
  233. \f$ - \nabla \cdot ( \sigma(\mathbf{r}) \nabla u(\mathbf{r}) ) = g(\mathbf{r})
  234. \f$.
  235. Typically this is used to define the source term.
  236. */
  237. vtkSmartPointer<vtkImplicitFunction> vtkG;
  238. /** Any type of vtkDataSet may be used as a calculation Grid */
  239. vtkSmartPointer<vtkUnstructuredGrid> vtkGrid;
  240. /** Function pointer to function describing g */
  241. Real (*gFcn3)(const Real&, const Real&, const Real&);
  242. Eigen::SparseMatrix<Real> A;
  243. VectorXr g;
  244. private:
  245. static constexpr auto CName = "FEM4EllipticPDE";
  246. }; // ----- end of class FEM4EllipticPDE -----
  247. /*! @} */
  248. } // ----- end of Lemma name -----
  249. #endif // ----- #ifndef FEM4ELLIPTICPDE_INC -----