Galerkin FEM for elliptic PDEs
Ви не можете вибрати більше 25 тем Теми мають розпочинатися з літери або цифри, можуть містити дефіси (-) і не повинні перевищувати 35 символів.

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325
  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 FEM4EllipticPDE* New();
  80. /** Deletes this object. Delete also disconnects any
  81. * attachments to this object.
  82. */
  83. void Delete();
  84. #ifdef HAVE_YAMLCPP
  85. /**
  86. * Uses YAML to serialize this object.
  87. * @return a YAML::Node
  88. */
  89. YAML::Node Serialize() const;
  90. /**
  91. * Constructs an object from a YAML::Node.
  92. */
  93. static FEM4EllipticPDE* DeSerialize(const YAML::Node& node);
  94. #endif
  95. // ==================== OPERATORS =======================
  96. // ==================== OPERATIONS =======================
  97. /** Sets the vtkImplictFunction that is used as sigma in
  98. \f$ - \nabla \cdot ( \sigma(\mathbf{r}) \nabla u(\mathbf{r}) ) = g(\mathbf{r})
  99. \f$. If this is not set, sigma evaluates as 1.
  100. */
  101. void SetSigmaFunction(vtkImplicitFunction* sigma);
  102. /** Sets the step to be used in applying boundary conditions.
  103. */
  104. void SetBoundaryStep(const Real& h);
  105. /** Sets the vtkImplictFunction that is used
  106. */
  107. void SetGFunction(vtkImplicitFunction* G);
  108. /** Uses a function pointer instead of a vtkImplicitFunction. For functions that can
  109. * be written in closed form, this can be much faster and more accurate then interpolating
  110. * off a grid.
  111. */
  112. void SetGFunction( Real (*pFcn3)(const Real&, const Real&, const Real&) );
  113. /* Sets the vtkDataSet that defines the calculation grid.
  114. */
  115. //void SetGrid(vtkDataSet* Grid);
  116. /** Sets the vtkDataSet that defines the calculation grid.
  117. */
  118. void SetGrid(vtkUnstructuredGrid* Grid);
  119. /**
  120. * Read grid in from VTK file
  121. */
  122. void SetVTUGridFile( std::string& vtkGridFile );
  123. /**
  124. * Read grid in from VTK file
  125. */
  126. void SetVTUGridFile( char* vtkGridFile );
  127. /** Sets up a DC problem with a survey
  128. * @param[in] ij is the injection index
  129. */
  130. void SetupDC(DCSurvey* Survey, const int& ij);
  131. /** Sets up the potential problem from the VTK file
  132. */
  133. void SetupPotential();
  134. /** Performs solution */
  135. void Solve(const std::string& fname);
  136. /* Performs solution */
  137. //void SolveOLD2(const std::string& fname);
  138. // ==================== ACCESS =======================
  139. // ==================== INQUIRY =======================
  140. protected:
  141. // ==================== LIFECYCLE =======================
  142. /// Default protected constructor.
  143. FEM4EllipticPDE (const std::string& cname);
  144. #ifdef HAVE_YAMLCPP
  145. /** Protected DeDerializing constructor, use factory DeSerialize method*/
  146. FEM4EllipticPDE (const YAML::Node& node);
  147. #endif
  148. /** Default protected constructor. */
  149. ~FEM4EllipticPDE ();
  150. /**
  151. * @copybrief LemmaObject::Release()
  152. * @copydetails LemmaObject::Release()
  153. */
  154. void Release();
  155. // ==================== OPERATIONS =======================
  156. /** Used internally to construct stiffness matrix, A
  157. */
  158. void ConstructAMatrix();
  159. /* Used internally to construct the load vector, g
  160. */
  161. //void ConstructLoadVector();
  162. // This function is copyright (C) 2012 Joseph Capriotti
  163. /**
  164. Returns a vtkIdList of points in member vtkGrid connected to point ido.
  165. @param[in] ido is the point of interest
  166. @param[in] grid is a pointer to a vtkDataSet
  167. */
  168. vtkSmartPointer<vtkIdList> GetConnectedPoints(const int& id0);
  169. /** Uses Simpson's rule to numerically integrate
  170. * (Shamelessly adapted from http://en.wikipedia.org/wiki/Simpson's_rule
  171. */
  172. Real CompositeSimpsons(vtkImplicitFunction* f, Real r0[3], Real r1[3], int numInt);
  173. /** Uses Simpson's rule to numerically integrate
  174. * (Shamelessly adapted from http://en.wikipedia.org/wiki/Simpson's_rule
  175. * this method is for a static f
  176. */
  177. Real CompositeSimpsons(const Real& f, Real r0[3], Real r1[3], int numInt);
  178. /**
  179. * Uses Simpson's rule to numerically integrate two functions in 1 dimension over the points
  180. * r0 and r1
  181. */
  182. Real CompositeSimpsons2(vtkImplicitFunction* f, Real r0[3], Real r1[3], int numInt);
  183. /**
  184. * Uses Simpson's rule to numerically integrate two functions in 1 dimension over the points
  185. * r0 and r1, uses the Hat function and the function pointer fPtr
  186. */
  187. Real CompositeSimpsons2( Real (*fPtr)(const Real&, const Real&, const Real&), Real r0[3], Real r1[3], int numInt);
  188. /**
  189. * Uses Simpson's rule to numerically integrate two functions in 1 dimension over the points
  190. * r0 and r1, uses the Hat function. Assumes a constand function value f
  191. */
  192. Real CompositeSimpsons2( const Real& f, Real r0[3], Real r1[3], int numInt);
  193. /**
  194. * Uses Simpson's rule to numerically integrate two functions in 1 dimension over the points
  195. * r0 and r1, uses the Hat function. Assumes a linear function from f0 to f1.
  196. */
  197. Real CompositeSimpsons2( const Real& f0, const Real& f1, Real r0[3], Real r1[3], int numInt);
  198. /**
  199. * Uses Simpson's rule to numerically integrate three functions in 1 dimension over the points
  200. * r0 and r1, uses two Hat functions. Assumes a linear function from f0 to f1.
  201. */
  202. Real CompositeSimpsons3( const Real& f0, const Real& f1, Real r0[3], Real r1[3], int numInt);
  203. /** Hat function for use in Galerkin FEM method, actually a 1/2 Hat.
  204. */
  205. //Real hat(Real r[3], Real ri[3], Real rm1[3] );
  206. Real Hat(const Vector3r &r, const Vector3r& ri, const Vector3r& rm1 );
  207. /** Computes distance between r0 and r1
  208. */
  209. Real dist(Real r0[3], Real r1[3]);
  210. /** Computes distance between r0 and r1
  211. */
  212. Real dist(const Vector3r& r0, const Vector3r& r1);
  213. /** Sets up the potential problem from the VTK file
  214. */
  215. void SetupLineSourcePotential();
  216. /** Sets up the potential problem from the VTK file
  217. */
  218. void SetupSurfaceSourcePotential();
  219. /** Sets up the potential problem from the VTK file
  220. */
  221. void SetupPointSourcePotential();
  222. /** Sets up the potential problem from the VTK file
  223. */
  224. void SetupVolumeSourcePotential();
  225. // ==================== DATA MEMBERS =========================
  226. /** The effective step at the boundary condition. Defaults to 1
  227. */
  228. Real BndryH;
  229. /** The sigma value at the boundary condition. Defaults to 1
  230. */
  231. Real BndrySigma;
  232. /** Implicit function defining the \f$\sigma\f$ term in
  233. \f$ - \nabla \cdot ( \sigma(\mathbf{r}) \nabla u(\mathbf{r}) ) = g(\mathbf{r})
  234. \f$.
  235. If set to 1, the Poisson equation is solved. Any type of vtkImplicitFunction may
  236. be used, including vtkImplicitDataSet.
  237. */
  238. vtkSmartPointer<vtkImplicitFunction> vtkSigma;
  239. /** Implicit function defining the \f$g\f$ term in
  240. \f$ - \nabla \cdot ( \sigma(\mathbf{r}) \nabla u(\mathbf{r}) ) = g(\mathbf{r})
  241. \f$.
  242. Typically this is used to define the source term.
  243. */
  244. vtkSmartPointer<vtkImplicitFunction> vtkG;
  245. /** Any type of vtkDataSet may be used as a calculation Grid */
  246. vtkSmartPointer<vtkUnstructuredGrid> vtkGrid;
  247. /** Function pointer to function describing g */
  248. Real (*gFcn3)(const Real&, const Real&, const Real&);
  249. Eigen::SparseMatrix<Real> A;
  250. VectorXr g;
  251. private:
  252. }; // ----- end of class FEM4EllipticPDE -----
  253. /*! @} */
  254. } // ----- end of Lemma name -----
  255. #endif // ----- #ifndef FEM4ELLIPTICPDE_INC -----