123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444 |
- // ===========================================================================
- //
- // Filename: EMSchur3DBase.h
- //
- // Created: 09/20/2013 04:35:57 PM
- // Compiler: Tested with g++, icpc, and MSVC 2010
- //
- // Author: Trevor Irons (ti)
- //
- // Organisation: University of Utah,
- // Colorado School of Mines
- // US Geological Survey
- //
- // Email: tirons@egi.utah.edu
- //
- // ===========================================================================
-
- /**
- @file
- @author Trevor Irons
- @date 09/20/2013
- @version $Id$
- **/
- #pragma once
-
- #ifndef EMSCHUR3DBASE_INC
- #define EMSCHUR3DBASE_INC
-
- #include <LemmaCore>
- #include "EMSchur3DConfig.h"
- #include <FDEM1D>
-
-
- //#include "LemmaObject.h"
- //#include "rectilineargrid.h"
- //#include "RectilinearGridVTKExporter.h"
- //#include "ASCIIParser.h"
- //#include "AEMSurvey.h"
- //#include "receiverpoints.h"
- //#include "layeredearthem.h"
- //#include "emearth1d.h"
-
- #include "timer.h"
- #include <Eigen/Sparse>
- //#include "bicgstab.h"
-
- // Solvers
- #ifdef HAVE_PASTIX
- #include <Eigen/PaStiXSupport>
- #endif
-
- #ifdef HAVE_METIS
- #include <Eigen/MetisSupport>
- #endif
-
- #ifdef HAVE_SUPERLU
- #include <Eigen/SuperLUSupport>
- #endif
-
- #ifdef HAVE_SUPERLUMT
- #include <Eigen/SuperLUMTSupport>
- #endif
-
- #ifdef HAVE_SPQR
- #include <Eigen/SPQRSupport>
- #endif
-
- // Cholmod Support won't compile typedef issue
- // #ifdef HAVE_CHOLMOD
- // #include <Eigen/CholmodSupport>
- // #endif
- //
- // // Cholmod Support won't compile
- // #ifdef HAVE_UMFPACK
- // #include <Eigen/UmfPackSupport>
- // #endif
-
- namespace Lemma {
-
- /**
- \defgroup EMSchur3DBase EMSchur3DBase
- Provides 3D solution to Maxwell's equations.
- */
-
- enum SOLVER{ SPARSELU, SimplicialLLT, SimplicialLDLT, BiCGStab, SparseQR };
-
-
- /**
- @class EMSchur3DBase
- \ingroup EMSchur3DBase
- \brief Provides a 3D solution to Maxwell's equations.
- \details 3D finite difference solution to maxwells equations
- using a SCHUR decomposition on a staggered grid.
- Performs a Schur decomposition on the vector scalar formulation of
- Maxwell's equations.
- \f[
- -\nabla^2 (\mathbf{A}) - \jmath \omega \mu \sigma \mathbf{A} - \mu \sigma \nabla (\phi) = - \mu \mathbf{J}_s
- \f]
-
- Which can be written in the functional form
- \f[ \begin{pmatrix}
- -\nabla^2 + \jmath \omega \mu \sigma & \mu \sigma \nabla \\
- \nabla \cdot & 0
- \end{pmatrix}
- \begin{pmatrix} \mathbf{A} \\ \phi \end{pmatrix}
- = \begin{pmatrix} \mathbf{s}_E \\ 0 \end{pmatrix}
- \f]
- Using the notation
- \f[ \begin{pmatrix}
- \mathbf{C} & \mathbf{B} \\
- \mathbf{D} & \mathbf{0}
- \end{pmatrix} \begin{pmatrix} \mathbf{A} \\ \phi \end{pmatrix} =
- \begin{pmatrix} \mathbf{s}_E \\ 0 \end{pmatrix}
- \f]
- Which is decomposed for seperate solutions to \f$ \mathbf{A}, \phi \f$ using a Schur decomposition
- \f[ \begin{matrix}
- \mathbf{D}\mathbf{C}^{-1}\mathbf{B} \phi & = \mathbf{D} \mathbf{C}^{-1} \mathbf{s}_E \\
- \mathbf{C}\mathbf{A} & = \mathbf{s}_E - \mathbf{G} \phi
- \end{matrix}
- \f]
-
- Where \f$ \mathbf{B} = \mu \sigma \mathbf{D}^T \f$. Additional algorithmic details may be found at
- @verbatim
- @inproceedings{doi:10.1190/segam2012-0896.1,
- author = {Trevor Irons and Yaoguo Li and Jason R. McKenna},
- title = {3D frequency-domain electromagnetics modeling using decoupled scalar and vector potentials},
- booktitle = {SEG Technical Program Expanded Abstracts 2012},
- chapter = {112},
- year = {2012},
- pages = {1-6},
- doi = {10.1190/segam2012-0896.1},
- URL = {http://library.seg.org/doi/abs/10.1190/segam2012-0896.1},
- eprint = {http://library.seg.org/doi/pdf/10.1190/segam2012-0896.1}
- }
- @endverbatim
- */
-
- //template< class Solver >
- class EMSchur3DBase : public LemmaObject {
-
- friend std::ostream &operator<<(std::ostream &stream,
- const EMSchur3DBase &ob);
-
- protected:
-
- //template<typename U>
- //friend class EMSchur3D;
-
- public:
-
- // ==================== LIFECYCLE =======================
-
- /** Default protected constructor, use New */
- explicit EMSchur3DBase ( const ctor_key& );
-
- /** Default protected constructor, use New */
- explicit EMSchur3DBase ( const YAML::Node& node, const ctor_key& );
-
- /** Default protected destructor, use Delete */
- virtual ~EMSchur3DBase ( );
-
- /**
- * Initialises antenna to contain no points, with no current
- * and no frequency. NumberOfTurns set to 1
- */
- static std::shared_ptr<EMSchur3DBase> NewSP();
-
- /*
- * Provides deep copy
- */
- //virtual std::shared_ptr<EMSchur3DBase> Clone() const ;
-
- /**
- * Uses YAML to serialize this object.
- * @return a YAML::Node
- */
- YAML::Node Serialize() const;
-
- /**
- * Constructs an object from a YAML::Node.
- */
- static std::shared_ptr<EMSchur3DBase> DeSerialize( const YAML::Node& node );
-
- // ==================== OPERATORS =======================
-
- // ==================== OPERATIONS =======================
-
- /* Performs the solution
- * This is thread safe. TODO, but where should the results go? Just to file?
- * Who does the parsing? Actually I think this method is the wrong place to talk
- * about that. This is just a big red button. The details should be worked out in private
- * methods, that this could well call. Still though, where should the damn results go. What if
- * someone wants to use them *right now*, and not go through file IO. This is a good cause for
- * fixing the model class. So the result will be a final RectGrid (or Model) where the results live.
- * THEN users can call a seperate WriteToVTK or whatever based on that.
- */
- void Solve( );
-
- // ==================== ACCESS =======================
-
- /** Sets the RectilinearGrid to use
- * @param[in] Grid is a pointer to the Grid to be used.
- */
- void SetRectilinearGrid( std::shared_ptr<RectilinearGrid> Grid);
-
- /** Sets the RectilinearGrid to use
- * @param[in] Grid is a pointer to the Grid to be used.
- */
- //void SetAEMSurvey(AEMSurvey* Survey);
-
- /** Sets the prefix for results files (.log and .vtr) the source fiducial is added as well
- */
- void SetResFileName(const std::string& filename);
-
- /** Sets the solver to use to invert the C matrix. This is done many times. If you are reusing the same matrix
- for numerous simulations, it may be benefitial to use the direct (SPARSELU) solver. For one off calculations the BiCGSTAB
- is a good choice. Default is SPARSELU.
- */
- void SetCSolver(const SOLVER& CSOLVER);
-
- /**
- * Sets the LayredEarthEM model that will be used for the primary field calculation as well as deterimining the
- * bacground conductivity file.
- @ @param[in] LayModel is a pointer to the LayeredEarthEM model to use.
- */
- void SetLayeredEarthEM( std::shared_ptr<LayeredEarthEM> LayModel );
-
- /**
- * Loads a MeshTools conductivity model.
- * @param[in] fname is the file to load.
- */
- void LoadMeshToolsConductivityModel( const std::string& fname );
-
- /**
- * Writes out results to into a vtkRectilinearGrid file
- * @param[in] fname is the file name that is created, the .vtr suffix is added.
- */
- void WriteVTKResults( const std::string& fname, Eigen::Ref<VectorXcr> A,
- Eigen::Ref<VectorXcr> Se, Eigen::Ref<VectorXcr> E0, Eigen::Ref<VectorXcr> E,
- Eigen::Ref<VectorXcr> Phi, Eigen::Ref<VectorXcr> ADiv, Eigen::Ref<VectorXcr> ADiv2,
- Eigen::Ref<VectorXcr> B
- );
-
- // ==================== INQUIRY =======================
-
- /** Returns the name of the underlying class, similiar to Python's type */
- virtual std::string GetName() const {
- return this->CName;
- }
-
- protected:
-
- // ==================== LIFECYCLE =======================
-
- //private:
-
- /** Builds the C matrix */
- void BuildC( Real*** sigmax, Real*** sigmay, Real*** sigmaz, const int& iw);
-
- /* Builds the C matrix as a block real system. Benefits of this are the availability of an
- * LDL^T decomposition. Also, as complex number in C++ are templates and will necessarily have
- * real and imaginary parts, this formulation will have a reduced cost, due to less computations
- * with the zero valued imaginary parts (off-diagonals)
- * The \f$C \in I^3\f$ matrix is instead written as
- * [ C_r C_i ] [ x_r ] = [ b_r ]
- * [ C_i -C_r ] [ x_i ] [ b_i ]
- */
- //void BuildCReal( Real*** sigmax, Real*** sigmay, Real*** sigmaz, const int& iw);
-
- /** Builds the C matrix */
- void BuildCPreconditioner( const int& iw);
-
- /** Builds the C matrix */
- virtual void BuildCDirectSolver( )=0;
-
- /** Fills the actual points on the grid that 1D primary field calculations need to be made.
- @todo This is a little stupid as all threads share the same points. Stupid in that right now
- this is done for every calculation. Not a huge amount of time is spent here, I suppose some extra memory
- though. We need to extend
- EmEARTH to be able to input a grid so that points are not explicitly needed like
- this. This requires some care as calcs are made on faces.
- Alternatively, the bins function of FieldPoints could be extended quite easily.
- This may be the way to do this.
- */
- void FillPoints( std::shared_ptr<FieldPoints> Points );
-
- /** Builds D/G
- */
- void BuildD( );
-
- /** Used to manage tradititional C 3D array */
- template <typename T>
- void Allocate3DScalar(T ***&Array, T init) {
- Array = new T**[nx];
- for (int ix=0; ix<nx; ix++){
- Array[ix] = new T*[ny];
- for (int iy=0; iy<ny; iy++){
- Array[ix][iy] = new T[nz];
- for (int iz=0; iz<nz; iz++) Array[ix][iy][iz] = init;
- }
- }
- return;
- }
-
- /** Used to manage tradititional C 3D array */
- template <typename T>
- void Delete3DScalar(T ***&Array) {
- for (int ix=0; ix<nx; ix++){
- for (int iy=0; iy<ny; iy++){
- delete [] Array[ix][iy];
- }
- delete [] Array[ix];
- }
- delete [] Array;
- Array = NULL;
- return;
- }
-
- /**
- * This is called just before solve and gets all shared objects ready to go
- */
- void Setup( );
-
- /** Solves a single source problem. This method is thread safe.
- * @param[in] Source is the source term for generating primary fields
- * @param[in] isource is the source index
- */
- virtual void SolveSource( std::shared_ptr< DipoleSource > Source , const int& isource)=0;
-
- /** Computes the primary field
- */
- void PrimaryField( std::shared_ptr<DipoleSource> Source, std::shared_ptr<FieldPoints> dpoint);
-
- /**
- * Fills the vectors that are called in
- */
- void FillSourceTerms( Eigen::Ref<VectorXcr> ms,
- Eigen::Ref<VectorXcr> Se, Eigen::Ref<VectorXcr> E0,
- std::shared_ptr<FieldPoints> dpoint, const Real& omega );
-
- /** Computes the curl of A on the staggered grid
- */
- VectorXcr StaggeredGridCurl(Eigen::Ref<VectorXcr> A);
-
- // ==================== DATA MEMBERS =========================
-
- /** Grid over which operators are active */
- std::shared_ptr<RectilinearGrid> Grid;
-
- /* Used to help dump results */
- //std::shared_ptr<RectilinearGridVTKExporter> VTKGridExporter;
-
- /** Class containing information about the AEM survey */
- //AEMSurvey* Survey;
-
- std::shared_ptr<LayeredEarthEM> LayModel;
-
- /** Matrix objects representing discrete operators
- */
- Eigen::SparseMatrix<Complex, Eigen::ColMajor>* Cvec;
- Eigen::SparseMatrix<Complex, Eigen::ColMajor> D;
-
- /** Squeezed matrices for reduced phi grid
- */
- Eigen::SparseMatrix<Complex, Eigen::ColMajor>* Cvec_s;
- Eigen::SparseMatrix<Complex, Eigen::ColMajor> D_s;
-
- /** number of cells in x, set when RectilinearGrid is attached */
- int nx;
-
- /** number of cells in y, set when RectilinearGrid is attached */
- int ny;
-
- /** number of cells in z, set when RectilinearGrid is attached */
- int nz;
-
- /** number of fields/faces in x, unwrapped x */
- int unx;
-
- /** number of fields/faces in y, unwrapped y */
- int uny;
-
- /** number of fields/faces in z, unwrapped z */
- int unz;
-
- /** number of cell centres, unwrapped scalars */
- int uns;
-
- /** name for log files and VTK files */
- std::string ResFile;
-
- /** frequency of source */
- VectorXr Omegas;
-
- /** Conductivity model */
- //Complex ***sigma_jwe;
- Real ***sigma;
-
- /** Conductivity model minus reference model */
- //Complex ***sigmap;
- Real ***sigmap;
-
- /** rectilinear grid spacing in x */
- VectorXr hx;
-
- /** rectilinear grid spacing in y */
- VectorXr hy;
-
- /** rectilinear grid spacing in z */
- VectorXr hz;
-
- /** inverse of hx */
- VectorXr ihx;
-
- /** inverse of hx squared */
- VectorXr ihx2;
-
- /** inverse of hy */
- VectorXr ihy;
-
- /** inverse of hy squared */
- VectorXr ihy2;
-
- /** inverse of hz */
- VectorXr ihz;
-
- /** inverse of hz squared */
- VectorXr ihz2;
-
- /** Marker for air cells */
- VectorXi MAC;
-
- /** Marker for air cells */
- std::vector<int> idx;
-
- private:
-
- static constexpr auto CName = "EMSchur3DBase";
-
- }; // ----- end of class EMSchur3DBase -----
-
- }
-
- #endif // ----- #ifndef EMSCHUR3BASE_INC -----
|