3D EM based on Schur decomposition
Du kannst nicht mehr als 25 Themen auswählen Themen müssen mit entweder einem Buchstaben oder einer Ziffer beginnen. Sie können Bindestriche („-“) enthalten und bis zu 35 Zeichen lang sein.

EMSchur3DBase.h 14KB

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