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

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