3D EM based on Schur decomposition
Você não pode selecionar mais de 25 tópicos Os tópicos devem começar com uma letra ou um número, podem incluir traços ('-') e podem ter até 35 caracteres.

EMSchur3DBase.h 14KB

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