Galerkin FEM for elliptic PDEs
Вы не можете выбрать более 25 тем Темы должны начинаться с буквы или цифры, могут содержать дефисы(-) и должны содержать не более 35 символов.

LinearMag.cpp 10.0KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285
  1. /* This file is part of Lemma, a geophysical modelling and inversion API.
  2. * More information is available at http://lemmasoftware.org
  3. */
  4. /* This Source Code Form is subject to the terms of the Mozilla Public
  5. * License, v. 2.0. If a copy of the MPL was not distributed with this
  6. * file, You can obtain one at http://mozilla.org/MPL/2.0/.
  7. */
  8. /**
  9. * @file
  10. * @date 03/21/2016 02:10:08 PM
  11. * @version $Id$
  12. * @author Trevor Irons (ti)
  13. * @email tirons@egi.utah.edu
  14. * @copyright Copyright (c) 2016, University of Utah
  15. * @copyright Copyright (c) 2016, Lemma Software, LLC
  16. */
  17. #include "LinearMag.h"
  18. namespace Lemma {
  19. // ==================== FRIEND METHODS =====================
  20. #ifdef HAVE_YAMLCPP
  21. std::ostream &operator << (std::ostream &stream, const LinearMag &ob) {
  22. stream << ob.Serialize() << "\n---\n"; // End of doc --- as a direct stream should encapsulate thingy
  23. return stream;
  24. }
  25. #else
  26. std::ostream &operator<<(std::ostream &stream, const LinearMag& ob) {
  27. stream << *(FEM4EllipticPDE*)(&ob);
  28. return stream;
  29. }
  30. #endif
  31. // ==================== LIFECYCLE =======================
  32. //--------------------------------------------------------------------------------------
  33. // Class: LinearMag
  34. // Method: LinearMag
  35. // Description: constructor (protected)
  36. //--------------------------------------------------------------------------------------
  37. LinearMag::LinearMag (const std::string& name) : FEM4EllipticPDE(name) {
  38. } // ----- end of method LinearMag::LinearMag (constructor) -----
  39. #ifdef HAVE_YAMLCPP
  40. //--------------------------------------------------------------------------------------
  41. // Class: LinearMag
  42. // Method: LinearMag
  43. // Description: DeSerializing constructor (protected)
  44. //--------------------------------------------------------------------------------------
  45. LinearMag::LinearMag (const YAML::Node& node) : FEM4EllipticPDE(node) {
  46. } // ----- end of method LinearMag::LinearMag (constructor) -----
  47. #endif
  48. //--------------------------------------------------------------------------------------
  49. // Class: LinearMag
  50. // Method: New()
  51. // Description: public constructor
  52. //--------------------------------------------------------------------------------------
  53. LinearMag* LinearMag::New() {
  54. LinearMag* Obj = new LinearMag("LinearMag");
  55. Obj->AttachTo(Obj);
  56. return Obj;
  57. }
  58. //--------------------------------------------------------------------------------------
  59. // Class: LinearMag
  60. // Method: ~LinearMag
  61. // Description: destructor (protected)
  62. //--------------------------------------------------------------------------------------
  63. LinearMag::~LinearMag () {
  64. } // ----- end of method LinearMag::~LinearMag (destructor) -----
  65. //--------------------------------------------------------------------------------------
  66. // Class: LinearMag
  67. // Method: Delete
  68. // Description: public destructor
  69. //--------------------------------------------------------------------------------------
  70. void LinearMag::Delete() {
  71. this->DetachFrom(this);
  72. }
  73. //--------------------------------------------------------------------------------------
  74. // Class: LinearMag
  75. // Method: Release
  76. // Description: destructor (protected)
  77. //--------------------------------------------------------------------------------------
  78. void LinearMag::Release() {
  79. delete this;
  80. }
  81. #ifdef HAVE_YAMLCPP
  82. //--------------------------------------------------------------------------------------
  83. // Class: LinearMag
  84. // Method: Serialize
  85. //--------------------------------------------------------------------------------------
  86. YAML::Node LinearMag::Serialize ( ) const {
  87. YAML::Node node = FEM4EllipticPDE::Serialize();;
  88. node.SetTag( this->Name );
  89. // FILL IN CLASS SPECIFICS HERE
  90. node["B0"] = B0;
  91. return node;
  92. } // ----- end of method LinearMag::Serialize -----
  93. //--------------------------------------------------------------------------------------
  94. // Class: LinearMag
  95. // Method: DeSerialize
  96. //--------------------------------------------------------------------------------------
  97. LinearMag* LinearMag::DeSerialize ( const YAML::Node& node ) {
  98. LinearMag* Object = new LinearMag(node);
  99. Object->AttachTo(Object);
  100. DESERIALIZECHECK( node, Object )
  101. return Object ;
  102. } // ----- end of method LinearMag::DeSerialize -----
  103. #endif
  104. //--------------------------------------------------------------------------------------
  105. // Class: LinearMag
  106. // Method: SetInducingMagField
  107. //--------------------------------------------------------------------------------------
  108. void LinearMag::SetInducingMagField ( const Real& intensity, const Real& inc,
  109. const Real& dec, const MAGUNITS& U ) {
  110. B0(0) = intensity * std::cos(inc) * std::cos(dec); // northing
  111. B0(1) = intensity * std::cos(inc) * std::sin(dec); // easting
  112. B0(2) = intensity * std::sin(inc); // z
  113. ScaleB0(U);
  114. return ;
  115. } // ----- end of method LinearMag::SetInducingMagField -----
  116. //--------------------------------------------------------------------------------------
  117. // Class: LinearMag
  118. // Method: SetInducingMagFieldVector
  119. //--------------------------------------------------------------------------------------
  120. void LinearMag::SetInducingMagFieldVector ( const Vector3r& BB0, const MAGUNITS& U ) {
  121. B0 = BB0;
  122. ScaleB0(U);
  123. return ;
  124. } // ----- end of method LinearMag::SetInducingMagFieldVector -----
  125. //--------------------------------------------------------------------------------------
  126. // Class: LinearMag
  127. // Method: ScaleB0
  128. //--------------------------------------------------------------------------------------
  129. void LinearMag::ScaleB0 ( const MAGUNITS& U ) {
  130. switch ( U ) {
  131. case TESLA:
  132. break;
  133. case NANOTESLA:
  134. B0 *= 1e-9;
  135. break;
  136. case GAUSS:
  137. B0 *= 1e-4;
  138. break;
  139. }
  140. return ;
  141. } // ----- end of method LinearMag::ScaleB0 -----
  142. //--------------------------------------------------------------------------------------
  143. // Class: LinearMag
  144. // Method: CalculateRHS
  145. //--------------------------------------------------------------------------------------
  146. void LinearMag::CalculateRHS ( const std::string& susName ) {
  147. std::cout << "Calculating RHS...";
  148. std::cout.flush();
  149. if ( !vtkGrid->GetCellData()->GetScalars(susName.c_str()) ) {
  150. std::string err("No cell data by name ");
  151. err.append(susName);
  152. throw std::runtime_error(err.c_str());
  153. }
  154. if (!vtkGrid->GetNumberOfPoints()) {
  155. throw std::runtime_error("Number of points zero in input grid!");
  156. }
  157. vtkDoubleArray* G = vtkDoubleArray::New();
  158. G->SetNumberOfComponents(1);
  159. G->SetNumberOfTuples( vtkGrid->GetNumberOfPoints() );
  160. G->SetName("G");
  161. //g.resize(vtkGrid->GetNumberOfPoints());
  162. VectorXr GG = VectorXr::Zero( vtkGrid->GetNumberOfPoints() );
  163. // Iterate over all the points or all of the cells?
  164. for (int ic=0; ic < vtkGrid->GetNumberOfCells(); ++ic) {
  165. if ( vtkGrid->GetCell(ic)->GetNumberOfPoints() != 4 ) {
  166. throw std::runtime_error("Non-tetrahedral mesh encountered!");
  167. }
  168. Real cellSus = vtkGrid->GetCellData()->GetScalars(susName.c_str())->GetTuple(ic)[0];
  169. Eigen::Matrix<Real, 4, 4> C = Eigen::Matrix<Real, 4, 4>::Zero() ;
  170. for (int ii=0; ii<4; ++ii) {
  171. double* pts = vtkGrid->GetCell(ic)->GetPoints()->GetPoint(ii);
  172. C(ii, 0) = 1;
  173. C(ii, 1) = pts[0];
  174. C(ii, 2) = pts[1];
  175. C(ii, 3) = pts[2];
  176. }
  177. /* The indices */
  178. vtkIdList* Ids = vtkGrid->GetCell(ic)->GetPointIds();
  179. int ID[4];
  180. ID[0] = Ids->GetId(0);
  181. ID[1] = Ids->GetId(1);
  182. ID[2] = Ids->GetId(2);
  183. ID[3] = Ids->GetId(3);
  184. /* the 4 faces of the tetrahedra
  185. ID[0] ID[1] ID[2]
  186. ID[0] ID[1] ID[3]
  187. ID[0] ID[2] ID[3]
  188. ID[1] ID[2] ID[3]
  189. */
  190. // Face 0, ID 0,1,2
  191. Eigen::Matrix<Real, 3, 2> CC = Eigen::Matrix<Real, 3, 2>::Ones() ;
  192. CC.col(1) = C.row(0).tail<3>() - C.row(1).tail<3>();
  193. CC.col(2) = C.row(0).tail<3>() - C.row(2).tail<3>();
  194. Vector3r nhat = CC.col(1).cross(CC.col(2));
  195. nhat.array() /= nhat.norm();
  196. Real flux = cellSus*nhat.dot(B0);
  197. GG(ID[0]) += flux;
  198. GG(ID[1]) += flux;
  199. GG(ID[2]) += flux;
  200. // Face 1, ID 0,1,3
  201. {
  202. CC.col(1) = C.row(0).tail<3>() - C.row(1).tail<3>();
  203. CC.col(2) = C.row(0).tail<3>() - C.row(3).tail<3>();
  204. Vector3r nhat = CC.col(1).cross(CC.col(2));
  205. nhat.array() /= nhat.norm();
  206. Real flux = cellSus*nhat.dot(B0);
  207. GG(ID[0]) += flux;
  208. GG(ID[1]) += flux;
  209. GG(ID[3]) += flux;
  210. }
  211. // Face 2, ID 0,2,3
  212. {
  213. CC.col(1) = C.row(0).tail<3>() - C.row(2).tail<3>();
  214. CC.col(2) = C.row(0).tail<3>() - C.row(3).tail<3>();
  215. Vector3r nhat = CC.col(1).cross(CC.col(2));
  216. nhat.array() /= nhat.norm();
  217. Real flux = cellSus*nhat.dot(B0);
  218. GG(ID[0]) += flux;
  219. GG(ID[2]) += flux;
  220. GG(ID[3]) += flux;
  221. }
  222. // Face 3, ID 1,2,3
  223. {
  224. CC.col(1) = C.row(1).tail<3>() - C.row(2).tail<3>();
  225. CC.col(2) = C.row(1).tail<3>() - C.row(3).tail<3>();
  226. Vector3r nhat = CC.col(1).cross(CC.col(2));
  227. nhat.array() /= nhat.norm();
  228. Real flux = cellSus*nhat.dot(B0);
  229. GG(ID[1]) += flux;
  230. GG(ID[2]) += flux;
  231. GG(ID[3]) += flux;
  232. }
  233. }
  234. for (int ip=0; ip<vtkGrid->GetNumberOfPoints(); ++ip) {
  235. G->InsertTuple1( ip, GG(ip) );
  236. }
  237. vtkGrid->GetPointData()->AddArray( G );
  238. vtkGrid->GetPointData()->SetScalars( G );
  239. std::cout << "finished" << std::endl;
  240. return ;
  241. } // ----- end of method LinearMag::CalculateRHS -----
  242. } // ----- end of Lemma name -----