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 9.8KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274
  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 <Eigen/IterativeLinearSolvers>
  50. #include <Eigen/SparseLU>
  51. #include "DCSurvey.h"
  52. namespace Lemma {
  53. /** @defgroup FEM4EllipticPDE
  54. @brief Finite element solver for elliptic PDE's
  55. @details Uses Galerkin method.
  56. */
  57. /** \addtogroup FEM4EllipticPDE
  58. @{
  59. */
  60. // ===================================================================
  61. // Class: FEM4EllipticPDE
  62. /**
  63. @class FEM4EllipticPDE
  64. @brief Galerkin FEM solver for Elliptic PDE problems
  65. @details Solves problems of the type
  66. \f[ - \nabla \cdot ( \sigma(\mathbf{r}) \nabla u(\mathbf{r}) ) = g(\mathbf{r})
  67. \f]
  68. */
  69. class FEM4EllipticPDE : public LemmaObject {
  70. friend std::ostream &operator<<(std::ostream &stream,
  71. const FEM4EllipticPDE &ob);
  72. public:
  73. // ==================== LIFECYCLE =======================
  74. /** Returns a pointer to a new object of type FEM4EllipticPDE.
  75. * It allocates all necessary memory.
  76. */
  77. static FEM4EllipticPDE* New();
  78. /** Deletes this object. Delete also disconnects any
  79. * attachments to this object.
  80. */
  81. void Delete();
  82. // ==================== OPERATORS =======================
  83. // ==================== OPERATIONS =======================
  84. /** Sets the vtkImplictFunction that is used as sigma in
  85. \f$ - \nabla \cdot ( \sigma(\mathbf{r}) \nabla u(\mathbf{r}) ) = g(\mathbf{r})
  86. \f$. If this is not set, sigma evaluates as 1.
  87. */
  88. void SetSigmaFunction(vtkImplicitFunction* sigma);
  89. /** Sets the step to be used in applying boundary conditions.
  90. */
  91. void SetBoundaryStep(const Real& h);
  92. /** Sets the vtkImplictFunction that is used
  93. */
  94. void SetGFunction(vtkImplicitFunction* G);
  95. /** Uses a function pointer instead of a vtkImplicitFunction. For functions that can
  96. * be written in closed form, this can be much faster and more accurate then interpolating
  97. * off a grid.
  98. */
  99. void SetGFunction( Real (*pFcn3)(const Real&, const Real&, const Real&) );
  100. /** Sets the vtkDataSet that defines the calculation grid.
  101. */
  102. void SetGrid(vtkDataSet* Grid);
  103. /** Sets up a DC problem with a survey
  104. * @param[in] ij is the injection index
  105. */
  106. void SetupDC(DCSurvey* Survey, const int& ij);
  107. /** Sets up the potential problem from the VTK file
  108. */
  109. void SetupPotential();
  110. /** Performs solution */
  111. void Solve(const std::string& fname);
  112. /** Performs solution */
  113. void SolveOLD(const std::string& fname);
  114. // ==================== ACCESS =======================
  115. // ==================== INQUIRY =======================
  116. protected:
  117. // ==================== LIFECYCLE =======================
  118. /// Default protected constructor.
  119. FEM4EllipticPDE (const std::string& cname);
  120. /** Default protected constructor. */
  121. ~FEM4EllipticPDE ();
  122. /**
  123. * @copybrief LemmaObject::Release()
  124. * @copydetails LemmaObject::Release()
  125. */
  126. void Release();
  127. // ==================== OPERATIONS =======================
  128. /** Used internally to construct stiffness matrix, A
  129. */
  130. void ConstructAMatrix();
  131. /* Used internally to construct the load vector, g
  132. */
  133. //void ConstructLoadVector();
  134. // This function is copyright (C) 2012 Joseph Capriotti
  135. /**
  136. Returns a vtkIdList of points in member vtkGrid connected to point ido.
  137. @param[in] ido is the point of interest
  138. @param[in] grid is a pointer to a vtkDataSet
  139. */
  140. vtkSmartPointer<vtkIdList> GetConnectedPoints(const int& id0);
  141. /** Uses Simpson's rule to numerically integrate
  142. * (Shamelessly adapted from http://en.wikipedia.org/wiki/Simpson's_rule
  143. */
  144. Real CompositeSimpsons(vtkImplicitFunction* f, Real r0[3], Real r1[3], int numInt);
  145. /** Uses Simpson's rule to numerically integrate
  146. * (Shamelessly adapted from http://en.wikipedia.org/wiki/Simpson's_rule
  147. * this method is for a static f
  148. */
  149. Real CompositeSimpsons(const Real& f, Real r0[3], Real r1[3], int numInt);
  150. /**
  151. * Uses Simpson's rule to numerically integrate two functions in 1 dimension over the points
  152. * r0 and r1
  153. */
  154. Real CompositeSimpsons2(vtkImplicitFunction* f, Real r0[3], Real r1[3], int numInt);
  155. /**
  156. * Uses Simpson's rule to numerically integrate two functions in 1 dimension over the points
  157. * r0 and r1, uses the Hat function and the function pointer fPtr
  158. */
  159. Real CompositeSimpsons2( Real (*fPtr)(const Real&, const Real&, const Real&), Real r0[3], Real r1[3], int numInt);
  160. /**
  161. * Uses Simpson's rule to numerically integrate two functions in 1 dimension over the points
  162. * r0 and r1, uses the Hat function. Assumes a constand function value f
  163. */
  164. Real CompositeSimpsons2( const Real& f, Real r0[3], Real r1[3], int numInt);
  165. /**
  166. * Uses Simpson's rule to numerically integrate two functions in 1 dimension over the points
  167. * r0 and r1, uses the Hat function. Assumes a linear function from f0 to f1.
  168. */
  169. Real CompositeSimpsons2( const Real& f0, const Real& f1, Real r0[3], Real r1[3], int numInt);
  170. /**
  171. * Uses Simpson's rule to numerically integrate three functions in 1 dimension over the points
  172. * r0 and r1, uses two Hat functions. Assumes a linear function from f0 to f1.
  173. */
  174. Real CompositeSimpsons3( const Real& f0, const Real& f1, Real r0[3], Real r1[3], int numInt);
  175. /** Hat function for use in Galerkin FEM method, actually a 1/2 Hat.
  176. */
  177. //Real hat(Real r[3], Real ri[3], Real rm1[3] );
  178. Real Hat(const Vector3r &r, const Vector3r& ri, const Vector3r& rm1 );
  179. /** Computes distance between r0 and r1
  180. */
  181. Real dist(Real r0[3], Real r1[3]);
  182. /** Computes distance between r0 and r1
  183. */
  184. Real dist(const Vector3r& r0, const Vector3r& r1);
  185. // ==================== DATA MEMBERS =========================
  186. /** The effective step at the boundary condition. Defaults to 1
  187. */
  188. Real BndryH;
  189. /** The sigma value at the boundary condition. Defaults to 1
  190. */
  191. Real BndrySigma;
  192. /** Implicit function defining the \f$\sigma\f$ term in
  193. \f$ - \nabla \cdot ( \sigma(\mathbf{r}) \nabla u(\mathbf{r}) ) = g(\mathbf{r})
  194. \f$.
  195. If set to 1, the Poisson equation is solved. Any type of vtkImplicitFunction may
  196. be used, including vtkImplicitDataSet.
  197. */
  198. vtkSmartPointer<vtkImplicitFunction> vtkSigma;
  199. /** Implicit function defining the \f$g\f$ term in
  200. \f$ - \nabla \cdot ( \sigma(\mathbf{r}) \nabla u(\mathbf{r}) ) = g(\mathbf{r})
  201. \f$.
  202. Typically this is used to define the source term.
  203. */
  204. vtkSmartPointer<vtkImplicitFunction> vtkG;
  205. /** Any type of vtkDataSet may be used as a calculation Grid */
  206. vtkSmartPointer<vtkDataSet> vtkGrid;
  207. /** Function pointer to function describing g */
  208. Real (*gFcn3)(const Real&, const Real&, const Real&);
  209. Eigen::SparseMatrix<Real> A;
  210. VectorXr g;
  211. private:
  212. }; // ----- end of class FEM4EllipticPDE -----
  213. /*! @} */
  214. } // ----- end of Lemma name -----
  215. #endif // ----- #ifndef FEM4ELLIPTICPDE_INC -----