3D EM based on Schur decomposition
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.

EMSchur3DBase.h 13KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438
  1. // ===========================================================================
  2. //
  3. // Filename: EMSchur3DBase.h
  4. //
  5. // Created: 09/20/2013 04:35:57 PM
  6. // Compiler: Tested with g++, icpc, and MSVC 2010
  7. //
  8. // Author: Trevor Irons (ti)
  9. //
  10. // Organisation: University of Utah,
  11. // Colorado School of Mines
  12. // US Geological Survey
  13. //
  14. // Email: tirons@egi.utah.edu
  15. //
  16. // ===========================================================================
  17. /**
  18. @file
  19. @author Trevor Irons
  20. @date 09/20/2013
  21. @version $Id$
  22. **/
  23. #ifndef EMSCHUR3DBASE_INC
  24. #define EMSCHUR3DBASE_INC
  25. #include <LemmaCore>
  26. #include <FDEM1D>
  27. //#include "LemmaObject.h"
  28. //#include "rectilineargrid.h"
  29. //#include "RectilinearGridVTKExporter.h"
  30. //#include "ASCIIParser.h"
  31. //#include "AEMSurvey.h"
  32. //#include "receiverpoints.h"
  33. //#include "layeredearthem.h"
  34. //#include "emearth1d.h"
  35. #include "timer.h"
  36. #include <Eigen/Sparse>
  37. #include "bicgstab.h"
  38. // Solvers
  39. #ifdef HAVE_PASTIX
  40. #include <Eigen/PaStiXSupport>
  41. #endif
  42. #ifdef HAVE_METIS
  43. #include <Eigen/MetisSupport>
  44. #endif
  45. #ifdef HAVE_SUPERLU
  46. #include <Eigen/SuperLUSupport>
  47. #endif
  48. #ifdef HAVE_SUPERLUMT
  49. #include <Eigen/SuperLUMTSupport>
  50. #endif
  51. #ifdef HAVE_SPQR
  52. #include <Eigen/SPQRSupport>
  53. #endif
  54. // Cholmod Support won't compile typedef issue
  55. // #ifdef HAVE_CHOLMOD
  56. // #include <Eigen/CholmodSupport>
  57. // #endif
  58. //
  59. // // Cholmod Support won't compile
  60. // #ifdef HAVE_UMFPACK
  61. // #include <Eigen/UmfPackSupport>
  62. // #endif
  63. namespace Lemma {
  64. /**
  65. \defgroup EMSchur3DBase EMSchur3DBase
  66. Provides 3D solution to Maxwell's equations.
  67. */
  68. enum SOLVER{ SPARSELU, SimplicialLLT, SimplicialLDLT, BiCGStab, SparseQR };
  69. /**
  70. @class EMSchur3DBase
  71. \ingroup EMSchur3DBase
  72. \brief Provides a 3D solution to Maxwell's equations.
  73. \details 3D finite difference solution to maxwells equations
  74. using a SCHUR decomposition on a staggered grid.
  75. Performs a Schur decomposition on the vector scalar formulation of
  76. Maxwell's equations.
  77. \f[
  78. -\nabla^2 (\mathbf{A}) - \jmath \omega \mu \sigma \mathbf{A} - \mu \sigma \nabla (\phi) = - \mu \mathbf{J}_s
  79. \f]
  80. Which can be written in the functional form
  81. \f[ \begin{pmatrix}
  82. -\nabla^2 + \jmath \omega \mu \sigma & \mu \sigma \nabla \\
  83. \nabla \cdot & 0
  84. \end{pmatrix}
  85. \begin{pmatrix} \mathbf{A} \\ \phi \end{pmatrix}
  86. = \begin{pmatrix} \mathbf{s}_E \\ 0 \end{pmatrix}
  87. \f]
  88. Using the notation
  89. \f[ \begin{pmatrix}
  90. \mathbf{C} & \mathbf{B} \\
  91. \mathbf{D} & \mathbf{0}
  92. \end{pmatrix} \begin{pmatrix} \mathbf{A} \\ \phi \end{pmatrix} =
  93. \begin{pmatrix} \mathbf{s}_E \\ 0 \end{pmatrix}
  94. \f]
  95. Which is decomposed for seperate solutions to \f$ \mathbf{A}, \phi \f$ using a Schur decomposition
  96. \f[ \begin{matrix}
  97. \mathbf{D}\mathbf{C}^{-1}\mathbf{B} \phi & = \mathbf{D} \mathbf{C}^{-1} \mathbf{s}_E \\
  98. \mathbf{C}\mathbf{A} & = \mathbf{s}_E - \mathbf{G} \phi
  99. \end{matrix}
  100. \f]
  101. Where \f$ \mathbf{B} = \mu \sigma \mathbf{D}^T \f$. Additional algorithmic details may be found at
  102. @verbatim
  103. @inproceedings{doi:10.1190/segam2012-0896.1,
  104. author = {Trevor Irons and Yaoguo Li and Jason R. McKenna},
  105. title = {3D frequency-domain electromagnetics modeling using decoupled scalar and vector potentials},
  106. booktitle = {SEG Technical Program Expanded Abstracts 2012},
  107. chapter = {112},
  108. year = {2012},
  109. pages = {1-6},
  110. doi = {10.1190/segam2012-0896.1},
  111. URL = {http://library.seg.org/doi/abs/10.1190/segam2012-0896.1},
  112. eprint = {http://library.seg.org/doi/pdf/10.1190/segam2012-0896.1}
  113. }
  114. @endverbatim
  115. */
  116. //template< class Solver >
  117. class EMSchur3DBase : public LemmaObject {
  118. friend std::ostream &operator<<(std::ostream &stream,
  119. const EMSchur3DBase &ob);
  120. protected:
  121. struct ctor_key {};
  122. public:
  123. // ==================== LIFECYCLE =======================
  124. /** Default protected constructor, use New */
  125. explicit EMSchur3DBase ( const ctor_key& );
  126. /** Default protected constructor, use New */
  127. explicit EMSchur3DBase ( const YAML::Node& node, const ctor_key& );
  128. /** Default protected destructor, use Delete */
  129. virtual ~EMSchur3DBase ();
  130. /**
  131. * Initialises antenna to contain no points, with no current
  132. * and no frequency. NumberOfTurns set to 1
  133. */
  134. static std::shared_ptr<EMSchur3DBase> NewSP();
  135. /**
  136. * Provides deep copy
  137. */
  138. virtual std::shared_ptr<EMSchur3DBase> Clone() const ;
  139. /**
  140. * Uses YAML to serialize this object.
  141. * @return a YAML::Node
  142. */
  143. YAML::Node Serialize() const;
  144. /**
  145. * Constructs an object from a YAML::Node.
  146. */
  147. static std::shared_ptr<EMSchur3DBase> DeSerialize( const YAML::Node& node );
  148. // ==================== OPERATORS =======================
  149. // ==================== OPERATIONS =======================
  150. /* Performs the solution
  151. * This is thread safe. TODO, but where should the results go? Just to file?
  152. * Who does the parsing? Actually I think this method is the wrong place to talk
  153. * about that. This is just a big red button. The details should be worked out in private
  154. * methods, that this could well call. Still though, where should the damn results go. What if
  155. * someone wants to use them *right now*, and not go through file IO. This is a good cause for
  156. * fixing the model class. So the result will be a final RectGrid (or Model) where the results live.
  157. * THEN users can call a seperate WriteToVTK or whatever based on that.
  158. */
  159. void Solve( );
  160. // ==================== ACCESS =======================
  161. /** Sets the RectilinearGrid to use
  162. * @param[in] Grid is a pointer to the Grid to be used.
  163. */
  164. void SetRectilinearGrid( std::shared_ptr<RectilinearGrid> Grid);
  165. /** Sets the RectilinearGrid to use
  166. * @param[in] Grid is a pointer to the Grid to be used.
  167. */
  168. //void SetAEMSurvey(AEMSurvey* Survey);
  169. /** Sets the prefix for results files (.log and .vtr) the source fiducial is added as well
  170. */
  171. void SetResFileName(const std::string& filename);
  172. /** Sets the solver to use to invert the C matrix. This is done many times. If you are reusing the same matrix
  173. for numerous simulations, it may be benefitial to use the direct (SPARSELU) solver. For one off calculations the BiCGSTAB
  174. is a good choice. Default is SPARSELU.
  175. */
  176. void SetCSolver(const SOLVER& CSOLVER);
  177. /**
  178. * Sets the LayredEarthEM model that will be used for the primary field calculation as well as deterimining the
  179. * bacground conductivity file.
  180. @ @param[in] LayModel is a pointer to the LayeredEarthEM model to use.
  181. */
  182. void SetLayeredEarthEM( std::shared_ptr<LayeredEarthEM> LayModel );
  183. /**
  184. * Loads a MeshTools conductivity model.
  185. * @param[in] fname is the file to load.
  186. */
  187. void LoadMeshToolsConductivityModel( const std::string& fname );
  188. /**
  189. * Writes out results to into a vtkRectilinearGrid file
  190. * @param[in] fname is the file name that is created, the .vtr suffix is added.
  191. */
  192. void WriteVTKResults( const std::string& fname, Eigen::Ref<VectorXcr> A,
  193. Eigen::Ref<VectorXcr> Se, Eigen::Ref<VectorXcr> E0, Eigen::Ref<VectorXcr> E,
  194. Eigen::Ref<VectorXcr> Phi, Eigen::Ref<VectorXcr> ADiv, Eigen::Ref<VectorXcr> ADiv2,
  195. Eigen::Ref<VectorXcr> B
  196. );
  197. // ==================== INQUIRY =======================
  198. /** Returns the name of the underlying class, similiar to Python's type */
  199. virtual std::string GetName() const;
  200. protected:
  201. // ==================== LIFECYCLE =======================
  202. //private:
  203. /** Builds the C matrix */
  204. void BuildC( Real*** sigmax, Real*** sigmay, Real*** sigmaz, const int& iw);
  205. /* Builds the C matrix as a block real system. Benefits of this are the availability of an
  206. * LDL^T decomposition. Also, as complex number in C++ are templates and will necessarily have
  207. * real and imaginary parts, this formulation will have a reduced cost, due to less computations
  208. * with the zero valued imaginary parts (off-diagonals)
  209. * The \f$C \in I^3\f$ matrix is instead written as
  210. * [ C_r C_i ] [ x_r ] = [ b_r ]
  211. * [ C_i -C_r ] [ x_i ] [ b_i ]
  212. */
  213. //void BuildCReal( Real*** sigmax, Real*** sigmay, Real*** sigmaz, const int& iw);
  214. /** Builds the C matrix */
  215. void BuildCPreconditioner( const int& iw);
  216. /** Builds the C matrix */
  217. virtual void BuildCDirectSolver( )=0;
  218. /** Fills the actual points on the grid that 1D primary field calculations need to be made.
  219. @todo This is a little stupid as all threads share the same points. Stupid in that right now
  220. this is done for every calculation. Not a huge amount of time is spent here, I suppose some extra memory
  221. though. We need to extend
  222. EmEARTH to be able to input a grid so that points are not explicitly needed like
  223. this. This requires some care as calcs are made on faces.
  224. Alternatively, the bins function of FieldPoints could be extended quite easily.
  225. This may be the way to do this.
  226. */
  227. void FillPoints( std::shared_ptr<FieldPoints> Points );
  228. /** Builds D/G
  229. */
  230. void BuildD( );
  231. /** Used to manage tradititional C 3D array */
  232. template <typename T>
  233. void Allocate3DScalar(T ***&Array, T init) {
  234. Array = new T**[nx];
  235. for (int ix=0; ix<nx; ix++){
  236. Array[ix] = new T*[ny];
  237. for (int iy=0; iy<ny; iy++){
  238. Array[ix][iy] = new T[nz];
  239. for (int iz=0; iz<nz; iz++) Array[ix][iy][iz] = init;
  240. }
  241. }
  242. return;
  243. }
  244. /** Used to manage tradititional C 3D array */
  245. template <typename T>
  246. void Delete3DScalar(T ***&Array) {
  247. for (int ix=0; ix<nx; ix++){
  248. for (int iy=0; iy<ny; iy++){
  249. delete [] Array[ix][iy];
  250. }
  251. delete [] Array[ix];
  252. }
  253. delete [] Array;
  254. Array = NULL;
  255. return;
  256. }
  257. /**
  258. * This is called just before solve and gets all shared objects ready to go
  259. */
  260. void Setup( );
  261. /** Solves a single source problem. This method is thread safe.
  262. * @param[in] Source is the source term for generating primary fields
  263. * @param[in] isource is the source index
  264. */
  265. virtual void SolveSource( DipoleSource* Source , const int& isource)=0;
  266. /** Computes the primary field
  267. */
  268. void PrimaryField( std::shared_ptr<DipoleSource> Source, std::shared_ptr<FieldPoints> dpoint);
  269. /**
  270. * Fills the vectors that are called in
  271. */
  272. void FillSourceTerms( Eigen::Ref<VectorXcr> ms,
  273. Eigen::Ref<VectorXcr> Se, Eigen::Ref<VectorXcr> E0,
  274. std::shared_ptr<FieldPoints> dpoint, const Real& omega );
  275. /** Computes the curl of A on the staggered grid
  276. */
  277. VectorXcr StaggeredGridCurl(Eigen::Ref<VectorXcr> A);
  278. // ==================== DATA MEMBERS =========================
  279. /** Grid over which operators are active */
  280. std::shared_ptr<RectilinearGrid> Grid;
  281. /* Used to help dump results */
  282. //std::shared_ptr<RectilinearGridVTKExporter> VTKGridExporter;
  283. /** Class containing information about the AEM survey */
  284. //AEMSurvey* Survey;
  285. std::shared_ptr<LayeredEarthEM> LayModel;
  286. /** Matrix objects representing discrete operators
  287. */
  288. Eigen::SparseMatrix<Complex, Eigen::ColMajor>* Cvec;
  289. Eigen::SparseMatrix<Complex, Eigen::ColMajor> D;
  290. /** Squeezed matrices for reduced phi grid
  291. */
  292. Eigen::SparseMatrix<Complex, Eigen::ColMajor>* Cvec_s;
  293. Eigen::SparseMatrix<Complex, Eigen::ColMajor> D_s;
  294. /** number of cells in x, set when RectilinearGrid is attached */
  295. int nx;
  296. /** number of cells in y, set when RectilinearGrid is attached */
  297. int ny;
  298. /** number of cells in z, set when RectilinearGrid is attached */
  299. int nz;
  300. /** number of fields/faces in x, unwrapped x */
  301. int unx;
  302. /** number of fields/faces in y, unwrapped y */
  303. int uny;
  304. /** number of fields/faces in z, unwrapped z */
  305. int unz;
  306. /** number of cell centres, unwrapped scalars */
  307. int uns;
  308. /** name for log files and VTK files */
  309. std::string ResFile;
  310. /** frequency of source */
  311. VectorXr Omegas;
  312. /** Conductivity model */
  313. //Complex ***sigma_jwe;
  314. Real ***sigma;
  315. /** Conductivity model minus reference model */
  316. //Complex ***sigmap;
  317. Real ***sigmap;
  318. /** rectilinear grid spacing in x */
  319. VectorXr hx;
  320. /** rectilinear grid spacing in y */
  321. VectorXr hy;
  322. /** rectilinear grid spacing in z */
  323. VectorXr hz;
  324. /** inverse of hx */
  325. VectorXr ihx;
  326. /** inverse of hx squared */
  327. VectorXr ihx2;
  328. /** inverse of hy */
  329. VectorXr ihy;
  330. /** inverse of hy squared */
  331. VectorXr ihy2;
  332. /** inverse of hz */
  333. VectorXr ihz;
  334. /** inverse of hz squared */
  335. VectorXr ihz2;
  336. /** Marker for air cells */
  337. VectorXi MAC;
  338. /** Marker for air cells */
  339. std::vector<int> idx;
  340. private:
  341. static constexpr auto CName = "EMSchur3DBase";
  342. }; // ----- end of class EMSchur3DBase -----
  343. }
  344. #endif // ----- #ifndef EMSCHUR3BASE_INC -----